transformPoints.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) 2017-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  transformPoints
29 
30 Group
31  grpMeshManipulationUtilities
32 
33 Description
34  Transforms the mesh points in the polyMesh directory according to the
35  translate, rotate and scale options.
36 
37 Usage
38  Options are:
39 
40  -time value
41  Specify the time to search from and apply the transformation
42  (default is latest)
43 
44  -recentre
45  Recentre using the bounding box centre before other operations
46 
47  -translate vector
48  Translate the points by the given vector before rotations
49 
50  -rotate (vector vector)
51  Rotate the points from the first vector to the second
52 
53  -rotate-angle (vector angle)
54  Rotate angle degrees about vector axis.
55 
56  -rotate-x angle
57  Rotate (degrees) about x-axis.
58 
59  -rotate-y angle
60  Rotate (degrees) about y-axis.
61 
62  -rotate-z angle
63  Rotate (degrees) about z-axis.
64 
65  or -yawPitchRoll : (yaw pitch roll) degrees
66  or -rollPitchYaw : (roll pitch yaw) degrees
67 
68  -scale scalar|vector
69  Scale the points by the given scalar or vector on output.
70 
71  The any or all of the three options may be specified and are processed
72  in the above order.
73 
74  With -rotateFields (in combination with -rotate/yawPitchRoll/rollPitchYaw)
75  it will also read & transform vector & tensor fields.
76 
77 Note
78  roll (rotation about x)
79  pitch (rotation about y)
80  yaw (rotation about z)
81 
82 \*---------------------------------------------------------------------------*/
83 
84 #include "argList.H"
85 #include "Time.H"
86 #include "fvMesh.H"
87 #include "volFields.H"
88 #include "surfaceFields.H"
89 #include "pointFields.H"
90 #include "ReadFields.H"
91 #include "regionProperties.H"
92 #include "transformField.H"
94 #include "axisAngleRotation.H"
96 
97 using namespace Foam;
98 using namespace Foam::coordinateRotations;
99 
100 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
101 
102 template<class GeoField>
103 void readAndRotateFields
104 (
105  PtrList<GeoField>& flds,
106  const fvMesh& mesh,
107  const dimensionedTensor& rotT,
108  const IOobjectList& objects
109 )
110 {
111  ReadFields(mesh, objects, flds);
112  for (GeoField& fld : flds)
113  {
114  Info<< "Transforming " << fld.name() << endl;
115  transform(fld, rotT, fld);
116  }
117 }
118 
119 
120 void rotateFields
121 (
122  const word& regionName,
123  const Time& runTime,
124  const tensor& rotT
125 )
126 {
127  // Need dimensionedTensor for geometric fields
128  const dimensionedTensor T(rotT);
129 
130  #include "createRegionMesh.H"
131 
132  // Read objects in time directory
133  IOobjectList objects(mesh, runTime.timeName());
134 
135  // Read vol fields.
136 
138  readAndRotateFields(vsFlds, mesh, T, objects);
139 
141  readAndRotateFields(vvFlds, mesh, T, objects);
142 
144  readAndRotateFields(vstFlds, mesh, T, objects);
145 
146  PtrList<volSymmTensorField> vsymtFlds;
147  readAndRotateFields(vsymtFlds, mesh, T, objects);
148 
150  readAndRotateFields(vtFlds, mesh, T, objects);
151 
152  // Read surface fields.
153 
155  readAndRotateFields(ssFlds, mesh, T, objects);
156 
158  readAndRotateFields(svFlds, mesh, T, objects);
159 
161  readAndRotateFields(sstFlds, mesh, T, objects);
162 
164  readAndRotateFields(ssymtFlds, mesh, T, objects);
165 
167  readAndRotateFields(stFlds, mesh, T, objects);
168 
169  mesh.write();
170 }
171 
172 
173 // Retrieve scaling option
174 // - size 0 : no scaling
175 // - size 1 : uniform scaling
176 // - size 3 : non-uniform scaling
177 List<scalar> getScalingOpt(const word& optName, const argList& args)
178 {
179  // readListIfPresent handles single or multiple values
180  // - accept 1 or 3 values
181 
182  List<scalar> scaling;
183  args.readListIfPresent(optName, scaling);
184 
185  if (scaling.size() == 1)
186  {
187  // Uniform scaling
188  }
189  else if (scaling.size() == 3)
190  {
191  // Non-uniform, but may actually be uniform
192  if
193  (
194  equal(scaling[0], scaling[1])
195  && equal(scaling[0], scaling[2])
196  )
197  {
198  scaling.resize(1);
199  }
200  }
201  else if (!scaling.empty())
202  {
203  FatalError
204  << "Incorrect number of components, must be 1 or 3." << nl
205  << " -" << optName << ' ' << args[optName].c_str() << endl
206  << exit(FatalError);
207  }
208 
209  if (scaling.size() == 1 && equal(scaling[0], 1))
210  {
211  // Scale factor 1 == no scaling
212  scaling.clear();
213  }
214 
215  // Zero and negative scaling are permitted
216 
217  return scaling;
218 }
219 
220 
221 void applyScaling(pointField& points, const List<scalar>& scaling)
222 {
223  if (scaling.size() == 1)
224  {
225  Info<< "Scaling points uniformly by " << scaling[0] << nl;
226  points *= scaling[0];
227  }
228  else if (scaling.size() == 3)
229  {
230  Info<< "Scaling points by ("
231  << scaling[0] << ' '
232  << scaling[1] << ' '
233  << scaling[2] << ')' << nl;
234 
235  points.replace(vector::X, scaling[0]*points.component(vector::X));
237  points.replace(vector::Z, scaling[2]*points.component(vector::Z));
238  }
239 }
240 
241 
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 
244 int main(int argc, char *argv[])
245 {
247  (
248  "Transform (translate / rotate / scale) mesh points.\n"
249  "Note: roll=rotate about x, pitch=rotate about y, yaw=rotate about z"
250  );
252  (
253  "time",
254  "time",
255  "Specify the time to search from and apply the transformation"
256  " (default is latest)"
257  );
259  (
260  "recentre",
261  "Recentre the bounding box before other operations"
262  );
264  (
265  "translate",
266  "vector",
267  "Translate by specified <vector> before rotations"
268  );
270  (
271  "auto-centre",
272  "Use bounding box centre as centre for rotations"
273  );
275  (
276  "centre",
277  "point",
278  "Use specified <point> as centre for rotations"
279  );
280  argList::addOptionCompat("auto-centre", {"auto-origin", 2206});
281  argList::addOptionCompat("centre", {"origin", 2206});
282 
284  (
285  "rotate",
286  "(vectorA vectorB)",
287  "Rotate from <vectorA> to <vectorB> - eg, '((1 0 0) (0 0 1))'"
288  );
290  (
291  "rotate-angle",
292  "(vector angle)",
293  "Rotate <angle> degrees about <vector> - eg, '((1 0 0) 45)'"
294  );
296  (
297  "rotate-x", "deg",
298  "Rotate (degrees) about x-axis"
299  );
301  (
302  "rotate-y", "deg",
303  "Rotate (degrees) about y-axis"
304  );
306  (
307  "rotate-z", "deg",
308  "Rotate (degrees) about z-axis"
309  );
311  (
312  "rollPitchYaw",
313  "vector",
314  "Rotate by '(roll pitch yaw)' degrees"
315  );
317  (
318  "yawPitchRoll",
319  "vector",
320  "Rotate by '(yaw pitch roll)' degrees"
321  );
323  (
324  "rotateFields",
325  "Read and transform vector and tensor fields too"
326  );
328  (
329  "scale",
330  "scalar | vector",
331  "Scale by the specified amount - Eg, for uniform [mm] to [m] scaling "
332  "use either '(0.001 0.001 0.001)' or simply '0.001'"
333  );
334 
335  // Compatibility with surfaceTransformPoints
336  argList::addOptionCompat("scale", {"write-scale", 0});
337 
338  #include "addAllRegionOptions.H"
339  #include "setRootCase.H"
340 
341  const bool doRotateFields = args.found("rotateFields");
342 
343  // Verify that an operation has been specified
344  {
345  const List<word> operationNames
346  ({
347  "recentre",
348  "translate",
349  "rotate",
350  "rotate-angle",
351  "rotate-x",
352  "rotate-y",
353  "rotate-z",
354  "rollPitchYaw",
355  "yawPitchRoll",
356  "scale"
357  });
358 
359  if (!args.count(operationNames))
360  {
361  FatalError
362  << "No operation supplied, "
363  << "use at least one of the following:" << nl
364  << " ";
365 
366  for (const auto& opName : operationNames)
367  {
368  FatalError
369  << " -" << opName;
370  }
371 
372  FatalError
373  << nl << exit(FatalError);
374  }
375  }
376 
377  // ------------------------------------------------------------------------
378 
379  #include "createTime.H"
380 
381  if (args.found("time"))
382  {
383  if (args["time"] == "constant")
384  {
385  runTime.setTime(instant(0, "constant"), 0);
386  }
387  else
388  {
389  const scalar timeValue = args.get<scalar>("time");
390  runTime.setTime(instant(timeValue), 0);
391  }
392  }
393 
394  // Handle -allRegions, -regions, -region
395  #include "getAllRegionOptions.H"
396 
397  // ------------------------------------------------------------------------
398 
399  forAll(regionNames, regioni)
400  {
401  const word& regionName = regionNames[regioni];
402  const fileName meshDir
403  (
405  );
406 
407  if (regionNames.size() > 1)
408  {
409  Info<< "region=" << regionName << nl;
410  }
411 
413  (
414  IOobject
415  (
416  "points",
417  runTime.findInstance(meshDir, "points"),
418  meshDir,
419  runTime,
422  false
423  )
424  );
425 
426 
427  // Begin operations
428 
429  vector v;
430  if (args.found("recentre"))
431  {
432  v = boundBox(points).centre();
433  Info<< "Adjust centre " << v << " -> (0 0 0)" << endl;
434  points -= v;
435  }
436 
437  if (args.readIfPresent("translate", v))
438  {
439  Info<< "Translating points by " << v << endl;
440  points += v;
441  }
442 
443  vector rotationCentre;
444  bool useRotationCentre = args.readIfPresent("centre", rotationCentre);
445  if (args.found("auto-centre") && !useRotationCentre)
446  {
447  useRotationCentre = true;
448  rotationCentre = boundBox(points).centre();
449  }
450 
451  if (useRotationCentre)
452  {
453  Info<< "Set centre of rotation to " << rotationCentre << endl;
454  points -= rotationCentre;
455  }
456 
457 
458  // Get a rotation specification
459 
460  tensor rot(Zero);
461  bool useRotation(false);
462 
463  if (args.found("rotate"))
464  {
465  Pair<vector> n1n2
466  (
467  args.lookup("rotate")()
468  );
469  n1n2[0].normalise();
470  n1n2[1].normalise();
471 
472  rot = rotationTensor(n1n2[0], n1n2[1]);
473  useRotation = true;
474  }
475  else if (args.found("rotate-angle"))
476  {
477  const Tuple2<vector, scalar> rotAxisAngle
478  (
479  args.lookup("rotate-angle")()
480  );
481 
482  const vector& axis = rotAxisAngle.first();
483  const scalar angle = rotAxisAngle.second();
484 
485  Info<< "Rotating points " << nl
486  << " about " << axis << nl
487  << " angle " << angle << nl;
488 
489  rot = axisAngle::rotation(axis, angle, true);
490  useRotation = true;
491  }
492  else if (args.found("rotate-x"))
493  {
494  const scalar angle = args.get<scalar>("rotate-x");
495 
496  Info<< "Rotating points about x-axis: " << angle << nl;
497 
498  rot = axisAngle::rotation(vector::X, angle, true);
499  useRotation = true;
500  }
501  else if (args.found("rotate-y"))
502  {
503  const scalar angle = args.get<scalar>("rotate-y");
504 
505  Info<< "Rotating points about y-axis: " << angle << nl;
506 
507  rot = axisAngle::rotation(vector::Y, angle, true);
508  useRotation = true;
509  }
510  else if (args.found("rotate-z"))
511  {
512  const scalar angle = args.get<scalar>("rotate-z");
513 
514  Info<< "Rotating points about z-axis: " << angle << nl;
515 
516  rot = axisAngle::rotation(vector::Z, angle, true);
517  useRotation = true;
518  }
519  else if (args.readIfPresent("rollPitchYaw", v))
520  {
521  Info<< "Rotating points by" << nl
522  << " roll " << v.x() << nl
523  << " pitch " << v.y() << nl
524  << " yaw " << v.z() << nl;
525 
526  rot = euler::rotation(euler::eulerOrder::ROLL_PITCH_YAW, v, true);
527  useRotation = true;
528  }
529  else if (args.readIfPresent("yawPitchRoll", v))
530  {
531  Info<< "Rotating points by" << nl
532  << " yaw " << v.x() << nl
533  << " pitch " << v.y() << nl
534  << " roll " << v.z() << nl;
535 
536  rot = euler::rotation(euler::eulerOrder::YAW_PITCH_ROLL, v, true);
537  useRotation = true;
538  }
539 
540  if (useRotation)
541  {
542  Info<< "Rotating points by " << rot << endl;
543  transform(points, rot, points);
544 
545  if (doRotateFields)
546  {
547  rotateFields(regionName, runTime, rot);
548  }
549  }
550 
551  if (useRotationCentre)
552  {
553  Info<< "Unset centre of rotation from " << rotationCentre << endl;
554  points += rotationCentre;
555  }
556 
557  // Output scaling
558  applyScaling(points, getScalingOpt("scale", args));
559 
560 
561  // Set the precision of the points data to 10
563 
564  Info<< "Writing points into directory "
565  << runTime.relativePath(points.path()) << nl
566  << endl;
567  points.write();
568  }
569 
570  Info<< nl << "End" << nl << endl;
571 
572  return 0;
573 }
574 
575 
576 // ************************************************************************* //
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
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh >> &fields, const bool syncPar=true, const bool readOldTime=false)
Read Geometric fields of templated type.
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
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
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:132
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
Spatial transformation functions for GeometricField.
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:841
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1072
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
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:420
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
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:365
static void addOptionCompat(const word &optName, std::pair< const char *, int > compat)
Specify an alias for the option name.
Definition: argList.C:409
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Field reading functions for post-processing utilities.
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
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:587
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:413
Spatial transformation functions for primitive fields.
Foam::word regionName(Foam::polyMesh::defaultRegion)
Namespace for coordinate system rotations.
Definition: axesRotation.C:30
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: instant.H:46
dynamicFvMesh & mesh
const pointField & points
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:109
A class for handling words, derived from Foam::string.
Definition: word.H:63
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:118
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:376
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:977
ITstream lookup(const word &optName) const
Return an input stream from the named option.
Definition: argListI.H:177
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:575
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
bool readListIfPresent(const word &optName, List< T > &list) const
If named option is present, get a List of values treating a single entry like a list of size 1...
Definition: argListI.H:387
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
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))
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name...
Definition: instant.H:53
Required Variables.
bool equal(const T &s1, const T &s2)
Compare two values for equality.
Definition: doubleFloat.H:41
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:271
static tensor rotation(const vector &angles, bool degrees=false)
Rotation tensor calculated for the intrinsic Euler angles in z-x-z order.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
PtrList< volScalarField > & Y
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
static tensor rotation(const vector &axis, const scalar angle, bool degrees=false)
The rotation tensor for given axis/angle.
messageStream Info
Information stream (stdout output on master, null elsewhere)
bool readIfPresent(const word &optName, T &val) const
Read a value from the named option if present.
Definition: argListI.H:316
point centre() const
The centre (midpoint) of the bounding box.
Definition: boundBoxI.H:179
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.
label count(const UList< word > &optionNames) const
Return how many of the specified options were used.
Definition: argList.C:1843
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157