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-2022 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  // Objects of field type
59  IOobjectList fields(objects.lookupClass<GeoField>());
60 
61  forAllConstIters(fields, fieldIter)
62  {
63  GeoField fld(*fieldIter(), mesh);
64  Info<< " Rotating " << fld.name() << endl;
65  transform(fld, rotT, fld);
66  fld.write();
67  }
68 }
69 
70 
71 void rotateFields
72 (
73  const fvMesh& mesh,
74  const Time& runTime,
75  const tensor& rotationT
76 )
77 {
78  // Need dimensionedTensor for geometric fields
79  const dimensionedTensor rotT(rotationT);
80 
81  // Search for list of objects for this time
82  IOobjectList objects(mesh, runTime.timeName());
83 
84  ReadAndRotateFields<volVectorField>(mesh, objects, rotT);
85  ReadAndRotateFields<volSphericalTensorField>(mesh, objects, rotT);
86  ReadAndRotateFields<volSymmTensorField>(mesh, objects, rotT);
87  ReadAndRotateFields<volTensorField>(mesh, objects, rotT);
88 
89  ReadAndRotateFields<surfaceVectorField>(mesh, objects, rotT);
90  ReadAndRotateFields<surfaceSphericalTensorField>(mesh, objects, rotT);
91  ReadAndRotateFields<surfaceSymmTensorField>(mesh, objects, rotT);
92  ReadAndRotateFields<surfaceTensorField>(mesh, objects, rotT);
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97 
98 int main(int argc, char *argv[])
99 {
101  (
102  "Rotate mesh points and vector/tensor fields\n"
103  "Rotation from the <from> vector to the <to> vector"
104  );
105 
107 
108  argList::addArgument("from", "The vector to rotate from");
109  argList::addArgument("to", "The vector to rotate to");
110 
111  #include "addAllRegionOptions.H"
112  #include "setRootCase.H"
113 
114  const vector n1(args.get<vector>(1).normalise());
115  const vector n2(args.get<vector>(2).normalise());
116 
117  const tensor rotT(rotationTensor(n1, n2));
118 
119  // ------------------------------------------------------------------------
120 
121  #include "createTime.H"
122 
123  // Handle -allRegions, -regions, -region
124  #include "getAllRegionOptions.H"
125 
126  // ------------------------------------------------------------------------
127 
128  #include "createNamedMeshes.H"
129 
130  forAll(regionNames, regioni)
131  {
132  const word& regionName = regionNames[regioni];
133  const fileName meshDir
134  (
136  );
137 
138  if (regionNames.size() > 1)
139  {
140  Info<< "region=" << regionName << nl;
141  }
142 
144  (
145  IOobject
146  (
147  "points",
148  runTime.findInstance(meshDir, "points"),
149  meshDir,
150  runTime,
153  false
154  )
155  );
156 
157  points = transform(rotT, points);
158 
159  // Set the precision of the points data to 10
161 
162  Info<< "Writing points into directory "
163  << runTime.relativePath(points.path()) << nl
164  << endl;
165  points.write();
166  }
167 
169 
170  forAll(timeDirs, timei)
171  {
172  runTime.setTime(timeDirs[timei], timei);
173 
174  Info<< "Time = " << runTime.timeName() << endl;
175 
176  forAll(regionNames, regioni)
177  {
178  const word& regionName = regionNames[regioni];
179  if (regionNames.size() > 1)
180  {
181  Info<< "region=" << regionName << nl;
182  }
183 
184  rotateFields(meshes[regioni], runTime, rotT);
185  }
186  }
187 
188  Info<< "\nEnd\n" << endl;
189 
190  return 0;
191 }
192 
193 
194 // ************************************************************************* //
Foam::surfaceFields.
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:453
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
A class for handling file names.
Definition: fileName.H:71
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
IOobjectList lookupClass(const char *clsName) const
The list of IOobjects with the given headerClassName.
Definition: IOobjectList.C:231
word findInstance(const fileName &dir, const word &name=word::null, IOobjectOption::readOption rOpt=IOobjectOption::MUST_READ, const word &stopInstance=word::null) const
Return time instance (location) of dir that contains the file name (eg, used in reading mesh data)...
Definition: Time.C:777
Spatial transformation functions for GeometricField.
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:841
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:402
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
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:80
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:416
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
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
Foam::word regionName(Foam::polyMesh::defaultRegion)
dynamicFvMesh & mesh
const pointField & points
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:977
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 time name of given scalar time formatted with the given precision.
Definition: Time.C:760
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:234
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition: argList.C:342
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:529
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:166
A primitive field of type <T> with automated input and output.
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:101
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Required Variables.