transformPoints.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  transformPoints
26 
27 Description
28  Transform (translate, rotate, scale) the mesh points, and optionally also
29  any vector and tensor fields.
30 
31 Usage
32  \b transformPoints "<transformations>" [OPTION]
33 
34  Supported transformations:
35  - \par translate=<translation vector>
36  Translational transformation by given vector
37  - \par rotate=(<n1 vector> <n2 vector>)
38  Rotational transformation from unit vector n1 to n2
39  - \par Rx=<angle [deg] about x-axis>
40  Rotational transformation by given angle about x-axis
41  - \par Ry=<angle [deg] about y-axis>
42  Rotational transformation by given angle about y-axis
43  - \par Rz=<angle [deg] about z-axis>
44  Rotational transformation by given angle about z-axis
45  - \par Ra=<axis vector> <angle [deg] about axis>
46  Rotational transformation by given angle about given axis
47  - \par scale=<x-y-z scaling vector>
48  Anisotropic scaling by the given vector in the x, y, z
49  coordinate directions
50 
51  Options:
52  - \par -rotateFields \n
53  Additionally transform vector and tensor fields.
54 
55  Example usage:
56  transformPoints \
57  "translate=(-0.05 -0.05 0), \
58  Rz=45, \
59  translate=(0.05 0.05 0)"
60 
61 See also
62  Foam::transformer
63  surfaceTransformPoints
64 
65 \*---------------------------------------------------------------------------*/
66 
67 #include "argList.H"
68 #include "fvMesh.H"
69 #include "regionProperties.H"
70 #include "volFields.H"
71 #include "surfaceFields.H"
72 #include "ReadFields.H"
73 #include "pointFields.H"
74 #include "transformField.H"
76 #include "unitConversion.H"
77 
78 using namespace Foam;
79 
80 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
81 
82 template<class GeoField>
83 void readAndRotateFields
84 (
85  PtrList<GeoField>& flds,
86  const fvMesh& mesh,
87  const tensor& T,
88  const IOobjectList& objects
89 )
90 {
91  ReadFields(mesh, objects, flds);
92  forAll(flds, i)
93  {
94  Info<< "Transforming " << flds[i].name() << endl;
95  dimensionedTensor dimT("t", flds[i].dimensions(), T);
96  transform(flds[i], dimT, flds[i]);
97  }
98 }
99 
100 
101 void rotateFields(const argList& args, const Time& runTime, const tensor& T)
102 {
103  #include "createNamedMesh.H"
104 
105  // Read objects in time directory
106  IOobjectList objects(mesh, runTime.timeName());
107 
108  // Read vol fields.
110  readAndRotateFields(vsFlds, mesh, T, objects);
112  readAndRotateFields(vvFlds, mesh, T, objects);
114  readAndRotateFields(vstFlds, mesh, T, objects);
116  readAndRotateFields(vsymtFlds, mesh, T, objects);
118  readAndRotateFields(vtFlds, mesh, T, objects);
119 
120  // Read surface fields.
122  readAndRotateFields(ssFlds, mesh, T, objects);
124  readAndRotateFields(svFlds, mesh, T, objects);
126  readAndRotateFields(sstFlds, mesh, T, objects);
128  readAndRotateFields(ssymtFlds, mesh, T, objects);
130  readAndRotateFields(stFlds, mesh, T, objects);
131 
132  mesh.write();
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 
138 int main(int argc, char *argv[])
139 {
140  const wordList supportedTransformations
141  (
142  {"translate", "rotate", "Rx", "Ry", "Rz", "Ra", "scale"}
143  );
144 
145  {
146  OStringStream supportedTransformationsStr;
147  supportedTransformationsStr << supportedTransformations << endl;
148 
150  (
151  "Transforms a mesh e.g.\n"
152  "transformPoints "
153  "\"translate=(-0.586 0 -0.156), "
154  "Ry=3.485, "
155  "translate=(0.586 0 0.156)\"\n\n"
156  "Supported transformations " + supportedTransformationsStr.str()
157  );
158  }
159 
160  argList::validArgs.append("transformations");
161 
163  (
164  "rotateFields",
165  "transform vector and tensor fields"
166  );
167 
168  #include "addRegionOption.H"
169  #include "addAllRegionsOption.H"
170  #include "setRootCase.H"
171  #include "createTime.H"
172 
173  const string transformationString(args[1]);
174 
175  #include "createTransforms.H"
176 
177  const wordList regionNames(selectRegionNames(args, runTime));
178 
179  const bool doRotateFields = args.optionFound("rotateFields");
180 
181  forAll(regionNames, regioni)
182  {
183  const word& regionName = regionNames[regioni];
184  const word& regionDir = Foam::regionDir(regionName);
185 
186  fileName meshDir(regionDir/polyMesh::meshSubDir);
187 
189  (
190  IOobject
191  (
192  "points",
193  runTime.findInstance(meshDir, "points"),
194  meshDir,
195  runTime,
198  false
199  )
200  );
201 
203 
204  if (doRotateFields)
205  {
206  rotateFields(args, runTime, transforms.T());
207  }
208 
209  // Set the precision of the points data to 10
211 
212  Info<< "Writing points into directory " << points.path() << nl << endl;
213  points.write();
214  }
215 
216  Info<< "End\n" << endl;
217 
218  return 0;
219 }
220 
221 
222 // ************************************************************************* //
Foam::surfaceFields.
PtrList< surfaceSphericalTensorField > sstFlds
transformer transforms
PtrList< volTensorField > vtFlds
Definition: readVolFields.H:15
PtrList< volSphericalTensorField > vstFlds
Definition: readVolFields.H:9
wordList ReadFields(const Mesh &mesh, const IOobjectList &objects, PtrList< GeoField > &fields, const bool syncPar=true)
Read all fields of the specified type.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A class for handling file names.
Definition: fileName.H:79
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Spatial transformation functions for FieldFields.
Unit conversion functions.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:312
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:458
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
word findInstance(const fileName &dir, const word &name=word::null, const IOobject::readOption rOpt=IOobject::MUST_READ, const word &stopInstance=word::null) const
Return the location of "dir" containing the file "name".
Definition: Time.C:659
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
Field reading functions for post-processing utilities.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
Generic dimensioned Type class.
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1093
PtrList< volVectorField > vvFlds
Definition: readVolFields.H:6
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
PtrList< volScalarField > vsFlds
Definition: readVolFields.H:3
Spatial transformation functions for primitive fields.
const word & regionDir(const word &regionName)
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:102
PtrList< surfaceSymmTensorField > ssymtFlds
PtrList< volSymmTensorField > vsymtFlds
Definition: readVolFields.H:12
PtrList< surfaceTensorField > stFlds
const tensor & T() const
Return the transformation tensor.
Definition: transformerI.H:94
static const char nl
Definition: Ostream.H:260
objects
PtrList< surfaceScalarField > ssFlds
PtrList< surfaceVectorField > svFlds
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
string str() const
Return the string.
messageStream Info
vector transformPosition(const vector &v) const
Transform the given position.
Definition: transformerI.H:153
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Output to memory buffer stream.
Definition: OStringStream.H:49
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
wordList selectRegionNames(const argList &args, const Time &runTime)
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477