rotateMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  rotateMesh
29 
30 Group
31  grpMeshManipulationUtilities
32 
33 Description
34  Rotates the mesh and fields from the direction n1 to direction n2.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "argList.H"
39 #include "timeSelector.H"
40 #include "Time.H"
41 #include "fvMesh.H"
42 #include "volFields.H"
43 #include "surfaceFields.H"
44 #include "regionProperties.H"
46 #include "IOobjectList.H"
47 
48 using namespace Foam;
49 
50 template<class GeoField>
51 void ReadAndRotateFields
52 (
53  const fvMesh& mesh,
54  const IOobjectList& objects,
55  const dimensionedTensor& rotT
56 )
57 {
58  for (const IOobject& io : objects.csorted<GeoField>())
59  {
60  GeoField fld(io, mesh);
61  Info<< " Rotating " << fld.name() << endl;
62  transform(fld, rotT, fld);
63  fld.write();
64  }
65 }
66 
67 
68 void rotateFields
69 (
70  const fvMesh& mesh,
71  const Time& runTime,
72  const tensor& rotationT
73 )
74 {
75  // Need dimensionedTensor for geometric fields
76  const dimensionedTensor rotT(rotationT);
77 
78  // Search for list of objects for this time
79  IOobjectList objects(mesh, runTime.timeName());
80 
81  ReadAndRotateFields<volVectorField>(mesh, objects, rotT);
82  ReadAndRotateFields<volSphericalTensorField>(mesh, objects, rotT);
83  ReadAndRotateFields<volSymmTensorField>(mesh, objects, rotT);
84  ReadAndRotateFields<volTensorField>(mesh, objects, rotT);
85 
86  ReadAndRotateFields<surfaceVectorField>(mesh, objects, rotT);
87  ReadAndRotateFields<surfaceSphericalTensorField>(mesh, objects, rotT);
88  ReadAndRotateFields<surfaceSymmTensorField>(mesh, objects, rotT);
89  ReadAndRotateFields<surfaceTensorField>(mesh, objects, rotT);
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
94 
95 int main(int argc, char *argv[])
96 {
98  (
99  "Rotate mesh points and vector/tensor fields\n"
100  "Rotation from the <from> vector to the <to> vector"
101  );
102 
104 
105  argList::addArgument("from", "The vector to rotate from");
106  argList::addArgument("to", "The vector to rotate to");
107 
108  #include "addAllRegionOptions.H"
109  #include "setRootCase.H"
110 
111  const vector n1(args.get<vector>(1).normalise());
112  const vector n2(args.get<vector>(2).normalise());
113 
114  const tensor rotT(rotationTensor(n1, n2));
115 
116  // ------------------------------------------------------------------------
117 
118  #include "createTime.H"
119 
120  // Handle -allRegions, -regions, -region
121  #include "getAllRegionOptions.H"
122 
123  // ------------------------------------------------------------------------
124 
125  #include "createNamedMeshes.H"
126 
127  forAll(regionNames, regioni)
128  {
129  const word& regionName = regionNames[regioni];
130  const fileName meshDir
131  (
133  );
134 
135  if (regionNames.size() > 1)
136  {
137  Info<< "region=" << regionName << nl;
138  }
139 
141  (
142  IOobject
143  (
144  "points",
145  runTime.findInstance(meshDir, "points"),
146  meshDir,
147  runTime,
151  )
152  );
153 
154  points = transform(rotT, points);
155 
156  // More precision (for points data)
158 
159  Info<< "Writing points into directory "
160  << runTime.relativePath(points.path()) << nl
161  << endl;
162  points.write();
163  }
164 
166 
167  forAll(timeDirs, timei)
168  {
169  runTime.setTime(timeDirs[timei], timei);
170 
171  Info<< "Time = " << runTime.timeName() << endl;
172 
173  forAll(regionNames, regioni)
174  {
175  const word& regionName = regionNames[regioni];
176  if (regionNames.size() > 1)
177  {
178  Info<< "region=" << regionName << nl;
179  }
180 
181  rotateFields(meshes[regioni], runTime, rotT);
182  }
183  }
184 
185  Info<< "\nEnd\n" << endl;
186 
187  return 0;
188 }
189 
190 
191 // ************************************************************************* //
Foam::surfaceFields.
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
A class for handling file names.
Definition: fileName.H:72
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
Spatial transformation functions for GeometricField.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
fileName relativePath(const fileName &input, const bool caseTag=false) const
Return the input relative to the globalPath by stripping off a leading value of the globalPath...
Definition: TimePathsI.H:103
wordList regionNames
Generic dimensioned Type class.
Ignore writing from objectRegistry::writeObject()
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:47
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dynamicFvMesh & mesh
const pointField & points
word findInstance(const fileName &directory, const word &name=word::null, IOobjectOption::readOption rOpt=IOobjectOption::MUST_READ, const word &stopInstance=word::null, const bool constant_fallback=true) const
Return time instance (location) of directory containing the file name (eg, used in reading mesh data)...
Definition: Time.C:725
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:902
static unsigned int minPrecision(unsigned int prec) noexcept
Set the minimum default precision.
Definition: IOstream.H:440
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options and also set the runTime to the first i...
Definition: timeSelector.C:260
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:271
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:841
Required Classes.
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition: argList.C:351
messageStream Info
Information stream (stdout output on master, null elsewhere)
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Foam::argList args(argc, argv)
Tensor of scalars, i.e. Tensor<scalar>.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
UPtrList< const IOobject > csorted() const
The sorted list of IOobjects with headerClassName == Type::typeName.
A primitive field of type <T> with automated input and output.
Do not request registration (bool: false)
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:101
Namespace for OpenFOAM.
Required Variables.