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-2018 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  Transforms the mesh points in the polyMesh directory according to the
29  translate, rotate and scale options.
30 
31 Usage
32  \b transformPoints [OPTION]
33 
34  Options:
35  - \par -translate <vector> \n
36  Translates the points by the given vector.
37 
38  - \par -rotate (<vector> <vector>) \n
39  Rotates the points from the first vector to the second.
40 
41  - \par -yawPitchRoll (<yawdegrees> <pitchdegrees> <rolldegrees>) \n
42  Alternative rotation specification:
43  yaw (rotation about z)
44  pitch (rotation about y)
45  roll (rotation about x)
46 
47  - \par -rollPitchYaw (<rolldegrees> <pitchdegrees> <yawdegrees>) \n
48  Alternative rotation specification:
49  roll (rotation about x)
50  pitch (rotation about y)
51  yaw (rotation about z)
52 
53  - \par -rotateFields \n
54  In combination with \a -rotate, \a -yawPitchRoll or \a -rollPitchYaw
55  additionally transform vector and tensor fields.
56 
57  - \par -scale <vector> \n
58  Scales the points by the given vector.
59 
60  Any or all of the three transformation option types may be specified and are
61  processed in the above order.
62 
63 \*---------------------------------------------------------------------------*/
64 
65 #include "argList.H"
66 #include "fvMesh.H"
67 #include "regionProperties.H"
68 #include "volFields.H"
69 #include "surfaceFields.H"
70 #include "ReadFields.H"
71 #include "pointFields.H"
72 #include "transformField.H"
74 #include "mathematicalConstants.H"
75 
76 using namespace Foam;
77 using namespace Foam::constant::mathematical;
78 
79 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
80 
81 template<class GeoField>
82 void readAndRotateFields
83 (
84  PtrList<GeoField>& flds,
85  const fvMesh& mesh,
86  const tensor& T,
87  const IOobjectList& objects
88 )
89 {
90  ReadFields(mesh, objects, flds);
91  forAll(flds, i)
92  {
93  Info<< "Transforming " << flds[i].name() << endl;
94  dimensionedTensor dimT("t", flds[i].dimensions(), T);
95  transform(flds[i], dimT, flds[i]);
96  }
97 }
98 
99 
100 void rotateFields(const argList& args, const Time& runTime, const tensor& T)
101 {
102  #include "createNamedMesh.H"
103 
104  // Read objects in time directory
105  IOobjectList objects(mesh, runTime.timeName());
106 
107  // Read vol fields.
108 
110  readAndRotateFields(vsFlds, mesh, T, objects);
111 
113  readAndRotateFields(vvFlds, mesh, T, objects);
114 
116  readAndRotateFields(vstFlds, mesh, T, objects);
117 
119  readAndRotateFields(vsymtFlds, mesh, T, objects);
120 
122  readAndRotateFields(vtFlds, mesh, T, objects);
123 
124  // Read surface fields.
125 
127  readAndRotateFields(ssFlds, mesh, T, objects);
128 
130  readAndRotateFields(svFlds, mesh, T, objects);
131 
133  readAndRotateFields(sstFlds, mesh, T, objects);
134 
136  readAndRotateFields(ssymtFlds, mesh, T, objects);
137 
139  readAndRotateFields(stFlds, mesh, T, objects);
140 
141  mesh.write();
142 }
143 
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
147 int main(int argc, char *argv[])
148 {
150  (
151  "translate",
152  "vector",
153  "translate by the specified <vector> - eg, '(1 0 0)'"
154  );
156  (
157  "rotate",
158  "(vectorA vectorB)",
159  "transform in terms of a rotation between <vectorA> and <vectorB> "
160  "- eg, '( (1 0 0) (0 0 1) )'"
161  );
163  (
164  "rollPitchYaw",
165  "vector",
166  "transform in terms of '(roll pitch yaw)' in degrees"
167  );
169  (
170  "yawPitchRoll",
171  "vector",
172  "transform in terms of '(yaw pitch roll)' in degrees"
173  );
175  (
176  "rotateFields",
177  "read and transform vector and tensor fields too"
178  );
180  (
181  "scale",
182  "vector",
183  "scale by the specified amount - eg, '(0.001 0.001 0.001)' for a "
184  "uniform [mm] to [m] scaling"
185  );
186 
187  #include "addRegionOption.H"
188  #include "addAllRegionsOption.H"
189  #include "setRootCase.H"
190  #include "createTime.H"
191 
192  const wordList regionNames(selectRegionNames(args, runTime));
193 
194  forAll(regionNames, regioni)
195  {
196  const word& regionName = regionNames[regioni];
197  const word& regionDir = Foam::regionDir(regionName);
198 
199  fileName meshDir(regionDir/polyMesh::meshSubDir);
200 
202  (
203  IOobject
204  (
205  "points",
206  runTime.findInstance(meshDir, "points"),
207  meshDir,
208  runTime,
211  false
212  )
213  );
214 
215  const bool doRotateFields = args.optionFound("rotateFields");
216 
217  // this is not actually stringent enough:
218  if (args.options().empty())
219  {
221  << "No options supplied, please use one or more of "
222  "-translate, -rotate or -scale options."
223  << exit(FatalError);
224  }
225 
226  vector v;
227  if (args.optionReadIfPresent("translate", v))
228  {
229  Info<< "Translating points by " << v << endl;
230 
231  points += v;
232  }
233 
234  if (args.optionFound("rotate"))
235  {
236  Pair<vector> n1n2
237  (
238  args.optionLookup("rotate")()
239  );
240  n1n2[0] /= mag(n1n2[0]);
241  n1n2[1] /= mag(n1n2[1]);
242  tensor T = rotationTensor(n1n2[0], n1n2[1]);
243 
244  Info<< "Rotating points by " << T << endl;
245 
246  points = transform(T, points);
247 
248  if (doRotateFields)
249  {
250  rotateFields(args, runTime, T);
251  }
252  }
253  else if (args.optionReadIfPresent("rollPitchYaw", v))
254  {
255  Info<< "Rotating points by" << nl
256  << " roll " << v.x() << nl
257  << " pitch " << v.y() << nl
258  << " yaw " << v.z() << nl;
259 
260  // Convert to radians
261  v *= pi/180.0;
262 
263  quaternion R(quaternion::rotationSequence::XYZ, v);
264 
265  Info<< "Rotating points by quaternion " << R << endl;
266  points = transform(R, points);
267 
268  if (doRotateFields)
269  {
270  rotateFields(args, runTime, R.R());
271  }
272  }
273  else if (args.optionReadIfPresent("yawPitchRoll", v))
274  {
275  Info<< "Rotating points by" << nl
276  << " yaw " << v.x() << nl
277  << " pitch " << v.y() << nl
278  << " roll " << v.z() << nl;
279 
280  // Convert to radians
281  v *= pi/180.0;
282 
283  scalar yaw = v.x();
284  scalar pitch = v.y();
285  scalar roll = v.z();
286 
287  quaternion R = quaternion(vector(0, 0, 1), yaw);
288  R *= quaternion(vector(0, 1, 0), pitch);
289  R *= quaternion(vector(1, 0, 0), roll);
290 
291  Info<< "Rotating points by quaternion " << R << endl;
292  points = transform(R, points);
293 
294  if (doRotateFields)
295  {
296  rotateFields(args, runTime, R.R());
297  }
298  }
299 
300  if (args.optionReadIfPresent("scale", v))
301  {
302  Info<< "Scaling points by " << v << endl;
303 
307  }
308 
309  // Set the precision of the points data to 10
311 
312  Info<< "Writing points into directory " << points.path() << nl << endl;
313  points.write();
314  }
315 
316  Info<< "End\n" << endl;
317 
318  return 0;
319 }
320 
321 
322 // ************************************************************************* //
Foam::surfaceFields.
PtrList< surfaceSphericalTensorField > sstFlds
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
tensor R() const
The rotation tensor corresponding the quaternion.
Definition: quaternionI.H:323
A class for handling file names.
Definition: fileName.H:79
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Spatial transformation functions for FieldFields.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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:461
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
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:649
const Cmpt & z() const
Definition: VectorI.H:87
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:626
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Generic dimensioned Type class.
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:47
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:472
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1035
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
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
const Cmpt & y() const
Definition: VectorI.H:81
PtrList< volScalarField > vsFlds
Definition: readVolFields.H:3
Spatial transformation functions for primitive fields.
const word & regionDir(const word &regionName)
const Foam::HashTable< string > & options() const
Return options.
Definition: argListI.H:90
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
const pointField & points
mathematical constants.
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
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:127
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:60
PtrList< volSymmTensorField > vsymtFlds
Definition: readVolFields.H:12
PtrList< surfaceTensorField > stFlds
const Cmpt & x() const
Definition: VectorI.H:75
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:460
static const char nl
Definition: Ostream.H:265
objects
PtrList< surfaceScalarField > ssFlds
PtrList< surfaceVectorField > svFlds
#define R(A, B, C, D, E, F, K, M)
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
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:117
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
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
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.
Definition: argListI.H:114