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-2021 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 word& regionDir =
134  (
136  );
137  const fileName meshDir = regionDir/polyMesh::meshSubDir;
138 
139  if (regionNames.size() > 1)
140  {
141  Info<< "region=" << regionName << nl;
142  }
143 
145  (
146  IOobject
147  (
148  "points",
149  runTime.findInstance(meshDir, "points"),
150  meshDir,
151  runTime,
154  false
155  )
156  );
157 
158  points = transform(rotT, points);
159 
160  // Set the precision of the points data to 10
162 
163  Info<< "Writing points into directory "
164  << runTime.relativePath(points.path()) << nl
165  << endl;
166  points.write();
167  }
168 
170 
171  forAll(timeDirs, timei)
172  {
173  runTime.setTime(timeDirs[timei], timei);
174 
175  Info<< "Time = " << runTime.timeName() << endl;
176 
177  forAll(regionNames, regioni)
178  {
179  const word& regionName = regionNames[regioni];
180  if (regionNames.size() > 1)
181  {
182  Info<< "region=" << regionName << nl;
183  }
184 
185  rotateFields(meshes[regioni], runTime, rotT);
186  }
187  }
188 
189  Info<< "\nEnd\n" << endl;
190 
191  return 0;
192 }
193 
194 
195 // ************************************************************************* //
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
volFields.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::Tensor
A templated (3 x 3) tensor of objects of <T> derived from MatrixSpace.
Definition: complexI.H:268
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
regionProperties.H
Foam::fileName
A class for handling file names.
Definition: fileName.H:71
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:57
Foam::polyMesh::defaultRegion
static word defaultRegion
Definition: polyMesh.H:314
transformGeometricField.H
Spatial transformation functions for GeometricField.
Foam::polyMesh::meshSubDir
static word meshSubDir
Definition: polyMesh.H:317
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
IOobjectList.H
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
surfaceFields.H
Foam::surfaceFields.
Foam::Vector::normalise
Vector< Cmpt > & normalise()
Definition: VectorI.H:116
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Definition: dimensionSet.C:514
Foam::argList::get
T get(const label index) const
Definition: argListI.H:271
regionNames
wordList regionNames
Definition: getAllRegionOptions.H:31
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::argList::addArgument
static void addArgument(const string &argName, const string &usage="")
Definition: argList.C:294
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::TimePaths::relativePath
fileName relativePath(const fileName &input, const bool caseTag=false) const
Definition: TimePathsI.H:80
Foam::Info
messageStream Info
argList.H
meshes
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
fld
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))
Definition: gmvOutputLagrangian.H:23
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:36
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:81
fvMesh.H
Foam
Definition: atmBoundaryLayer.C:26
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:51
Foam::Time::findInstance
word findInstance(const fileName &dir, const word &name=word::null, const IOobject::readOption rOpt=IOobject::MUST_READ, const word &stopInstance=word::null) const
Definition: Time.C:790
getAllRegionOptions.H
Priority.
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision() noexcept
Definition: IOstream.H:338
Time.H
setRootCase.H
addAllRegionOptions.H
Foam::nl
constexpr char nl
Definition: Ostream.H:424
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::timeSelector::addOptions
static void addOptions(const bool constant=true, const bool withZero=false)
Definition: timeSelector.C:95
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::Time::setTime
virtual void setTime(const Time &t)
Definition: Time.C:996
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::word::null
static const word null
Definition: word.H:78
timeSelector.H
createTime.H
Foam::rotationTensor
tensor rotationTensor(const vector &n1, const vector &n2)
Definition: transform.H:47
regionDir
const word & regionDir
Definition: findMeshDefinitionDict.H:28
Foam::timeSelector::select0
static instantList select0(Time &runTime, const argList &args)
Definition: timeSelector.C:228
Foam::IOobjectList::lookupClass
IOobjectList lookupClass(const char *clsName) const
Definition: IOobjectList.C:283
args
Foam::argList args(argc, argv)
createNamedMeshes.H
Required Variables.
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:181