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-2024 OpenCFD Ltd.
10  Copyright (C) 2024 Haakan Nilsson
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 Application
29  transformPoints
30 
31 Group
32  grpMeshManipulationUtilities
33 
34 Description
35  Transforms the mesh points in the polyMesh directory according to the
36  translate, rotate and scale options.
37 
38 Usage
39  Options are:
40 
41  -time value
42  Specify the time to search from and apply the transformation
43  (default is latest)
44 
45  -recentre
46  Recentre using the bounding box centre before other operations
47 
48  -translate vector
49  Translate the points by the given vector before rotations
50 
51  -rotate (vector vector)
52  Rotate the points from the first vector to the second
53 
54  -rotate-angle (vector angle)
55  Rotate angle degrees about vector axis.
56 
57  -rotate-x angle
58  Rotate (degrees) about x-axis.
59 
60  -rotate-y angle
61  Rotate (degrees) about y-axis.
62 
63  -rotate-z angle
64  Rotate (degrees) about z-axis.
65 
66  or -yawPitchRoll : (yaw pitch roll) degrees
67  or -rollPitchYaw : (roll pitch yaw) degrees
68 
69  -scale scalar|vector
70  Scale the points by the given scalar or vector on output.
71 
72  The any or all of the three options may be specified and are processed
73  in the above order.
74 
75  -cylToCart (originVector axisVector directionVector)
76  Tranform cylindrical coordinates to cartesian coordinates
77 
78  With -rotateFields (in combination with -rotate/yawPitchRoll/rollPitchYaw)
79  it will also read & transform vector & tensor fields.
80 
81 Note
82  roll (rotation about x)
83  pitch (rotation about y)
84  yaw (rotation about z)
85 
86 \*---------------------------------------------------------------------------*/
87 
88 #include "argList.H"
89 #include "Time.H"
90 #include "fvMesh.H"
91 #include "volFields.H"
92 #include "surfaceFields.H"
93 #include "pointFields.H"
94 #include "ReadFields.H"
95 #include "regionProperties.H"
96 #include "transformField.H"
98 #include "axisAngleRotation.H"
100 #include "cylindricalCS.H"
101 
102 using namespace Foam;
103 using namespace Foam::coordinateRotations;
104 
105 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
106 
107 template<class GeoField>
108 void readAndRotateFields
109 (
110  PtrList<GeoField>& flds,
111  const fvMesh& mesh,
112  const dimensionedTensor& rotT,
113  const IOobjectList& objects
114 )
115 {
116  ReadFields(mesh, objects, flds);
117  for (GeoField& fld : flds)
118  {
119  Info<< "Transforming " << fld.name() << endl;
120  transform(fld, rotT, fld);
121  }
122 }
123 
124 
125 void rotateFields
126 (
127  const word& regionName,
128  const Time& runTime,
129  const tensor& rotT
130 )
131 {
132  // Need dimensionedTensor for geometric fields
133  const dimensionedTensor T(rotT);
134 
135  #include "createRegionMesh.H"
136 
137  // Read objects in time directory
138  IOobjectList objects(mesh, runTime.timeName());
139 
140  // Read vol fields.
141 
143  readAndRotateFields(vsFlds, mesh, T, objects);
144 
146  readAndRotateFields(vvFlds, mesh, T, objects);
147 
149  readAndRotateFields(vstFlds, mesh, T, objects);
150 
151  PtrList<volSymmTensorField> vsymtFlds;
152  readAndRotateFields(vsymtFlds, mesh, T, objects);
153 
155  readAndRotateFields(vtFlds, mesh, T, objects);
156 
157  // Read surface fields.
158 
160  readAndRotateFields(ssFlds, mesh, T, objects);
161 
163  readAndRotateFields(svFlds, mesh, T, objects);
164 
166  readAndRotateFields(sstFlds, mesh, T, objects);
167 
169  readAndRotateFields(ssymtFlds, mesh, T, objects);
170 
172  readAndRotateFields(stFlds, mesh, T, objects);
173 
174  mesh.write();
175 }
176 
177 
178 // Retrieve scaling option
179 // - size 0 : no scaling
180 // - size 1 : uniform scaling
181 // - size 3 : non-uniform scaling
182 List<scalar> getScalingOpt(const word& optName, const argList& args)
183 {
184  // readListIfPresent handles single or multiple values
185  // - accept 1 or 3 values
186 
187  List<scalar> scaling;
188  args.readListIfPresent(optName, scaling);
189 
190  if (scaling.size() == 1)
191  {
192  // Uniform scaling
193  }
194  else if (scaling.size() == 3)
195  {
196  // Non-uniform, but may actually be uniform
197  if
198  (
199  equal(scaling[0], scaling[1])
200  && equal(scaling[0], scaling[2])
201  )
202  {
203  scaling.resize(1);
204  }
205  }
206  else if (!scaling.empty())
207  {
208  FatalError
209  << "Incorrect number of components, must be 1 or 3." << nl
210  << " -" << optName << ' ' << args[optName].c_str() << endl
211  << exit(FatalError);
212  }
213 
214  if (scaling.size() == 1 && equal(scaling[0], 1))
215  {
216  // Scale factor 1 == no scaling
217  scaling.clear();
218  }
219 
220  // Zero and negative scaling are permitted
221 
222  return scaling;
223 }
224 
225 
226 void applyScaling(pointField& points, const List<scalar>& scaling)
227 {
228  if (scaling.size() == 1)
229  {
230  Info<< "Scaling points uniformly by " << scaling[0] << nl;
231  points *= scaling[0];
232  }
233  else if (scaling.size() == 3)
234  {
235  const vector factor(scaling[0], scaling[1], scaling[2]);
236  Info<< "Scaling points by " << factor << nl;
237  cmptMultiply(points, points, factor);
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  );
335  (
336  "cylToCart",
337  "(originVec axisVec directionVec)",
338  "Tranform cylindrical coordinates to cartesian coordinates"
339  );
340 
341 
342  // Compatibility with surfaceTransformPoints
343  argList::addOptionCompat("scale", {"write-scale", 0});
344 
345  #include "addAllRegionOptions.H"
346  #include "setRootCase.H"
347 
348  const bool doRotateFields = args.found("rotateFields");
349 
350  // Verify that an operation has been specified
351  {
352  const List<word> operationNames
353  ({
354  "recentre",
355  "translate",
356  "rotate",
357  "rotate-angle",
358  "rotate-x",
359  "rotate-y",
360  "rotate-z",
361  "rollPitchYaw",
362  "yawPitchRoll",
363  "scale",
364  "cylToCart"
365  });
366 
367  if (!args.count(operationNames))
368  {
369  FatalError
370  << "No operation supplied, "
371  << "use at least one of the following:" << nl
372  << " ";
373 
374  for (const auto& opName : operationNames)
375  {
376  FatalError
377  << " -" << opName;
378  }
379 
380  FatalError
381  << nl << exit(FatalError);
382  }
383  }
384 
385  // ------------------------------------------------------------------------
386 
387  #include "createTime.H"
388 
389  if (args.found("time"))
390  {
391  if (args["time"] == "constant")
392  {
393  runTime.setTime(instant(0, "constant"), 0);
394  }
395  else
396  {
397  const scalar timeValue = args.get<scalar>("time");
398  runTime.setTime(instant(timeValue), 0);
399  }
400  }
401 
402  // Handle -allRegions, -regions, -region
403  #include "getAllRegionOptions.H"
404 
405  // ------------------------------------------------------------------------
406 
407  forAll(regionNames, regioni)
408  {
409  const word& regionName = regionNames[regioni];
410  const fileName meshDir
411  (
413  );
414 
415  if (regionNames.size() > 1)
416  {
417  Info<< "region=" << regionName << nl;
418  }
419 
421  (
422  IOobject
423  (
424  "points",
425  runTime.findInstance(meshDir, "points"),
426  meshDir,
427  runTime,
431  )
432  );
433 
434 
435  // Begin operations
436 
437  vector v;
438  if (args.found("recentre"))
439  {
440  v = boundBox(points).centre();
441  Info<< "Adjust centre " << v << " -> (0 0 0)" << endl;
442  points -= v;
443  }
444 
445  if (args.readIfPresent("translate", v))
446  {
447  Info<< "Translating points by " << v << endl;
448  points += v;
449  }
450 
451  vector rotationCentre;
452  bool useRotationCentre = args.readIfPresent("centre", rotationCentre);
453  if (args.found("auto-centre") && !useRotationCentre)
454  {
455  useRotationCentre = true;
456  rotationCentre = boundBox(points).centre();
457  }
458 
459  if (useRotationCentre)
460  {
461  Info<< "Set centre of rotation to " << rotationCentre << endl;
462  points -= rotationCentre;
463  }
464 
465 
466  // Get a rotation specification
467 
468  tensor rot(Zero);
469  bool useRotation(false);
470 
471  if (args.found("rotate"))
472  {
473  Pair<vector> n1n2
474  (
475  args.lookup("rotate")()
476  );
477  n1n2[0].normalise();
478  n1n2[1].normalise();
479 
480  rot = rotationTensor(n1n2[0], n1n2[1]);
481  useRotation = true;
482  }
483  else if (args.found("rotate-angle"))
484  {
485  const Tuple2<vector, scalar> rotAxisAngle
486  (
487  args.lookup("rotate-angle")()
488  );
489 
490  const vector& axis = rotAxisAngle.first();
491  const scalar angle = rotAxisAngle.second();
492 
493  Info<< "Rotating points " << nl
494  << " about " << axis << nl
495  << " angle " << angle << nl;
496 
497  rot = axisAngle::rotation(axis, angle, true);
498  useRotation = true;
499  }
500  else if (args.found("rotate-x"))
501  {
502  const scalar angle = args.get<scalar>("rotate-x");
503 
504  Info<< "Rotating points about x-axis: " << angle << nl;
505 
506  rot = axisAngle::rotation(vector::X, angle, true);
507  useRotation = true;
508  }
509  else if (args.found("rotate-y"))
510  {
511  const scalar angle = args.get<scalar>("rotate-y");
512 
513  Info<< "Rotating points about y-axis: " << angle << nl;
514 
515  rot = axisAngle::rotation(vector::Y, angle, true);
516  useRotation = true;
517  }
518  else if (args.found("rotate-z"))
519  {
520  const scalar angle = args.get<scalar>("rotate-z");
521 
522  Info<< "Rotating points about z-axis: " << angle << nl;
523 
524  rot = axisAngle::rotation(vector::Z, angle, true);
525  useRotation = true;
526  }
527  else if (args.readIfPresent("rollPitchYaw", v))
528  {
529  Info<< "Rotating points by" << nl
530  << " roll " << v.x() << nl
531  << " pitch " << v.y() << nl
532  << " yaw " << v.z() << nl;
533 
534  rot = euler::rotation(euler::eulerOrder::ROLL_PITCH_YAW, v, true);
535  useRotation = true;
536  }
537  else if (args.readIfPresent("yawPitchRoll", v))
538  {
539  Info<< "Rotating points by" << nl
540  << " yaw " << v.x() << nl
541  << " pitch " << v.y() << nl
542  << " roll " << v.z() << nl;
543 
544  rot = euler::rotation(euler::eulerOrder::YAW_PITCH_ROLL, v, true);
545  useRotation = true;
546  }
547 
548  if (useRotation)
549  {
550  Info<< "Rotating points by " << rot << endl;
551  transform(points, rot, points);
552 
553  if (doRotateFields)
554  {
555  rotateFields(regionName, runTime, rot);
556  }
557  }
558 
559  if (useRotationCentre)
560  {
561  Info<< "Unset centre of rotation from " << rotationCentre << endl;
562  points += rotationCentre;
563  }
564 
565  // Output scaling
566  applyScaling(points, getScalingOpt("scale", args));
567 
568  if (args.found("cylToCart"))
569  {
570  vectorField n1n2(args.lookup("cylToCart")());
571  n1n2[1].normalise();
572  n1n2[2].normalise();
573 
574  cylindricalCS ccs
575  (
576  "ccs",
577  n1n2[0],
578  n1n2[1],
579  n1n2[2]
580  );
581 
582  points = ccs.globalPosition(points);
583  }
584 
585  // More precision (for points data)
587 
588  Info<< "Writing points into directory "
589  << runTime.relativePath(points.path()) << nl
590  << endl;
591  points.write();
592  }
593 
594  Info<< nl << "End" << nl << endl;
595 
596  return 0;
597 }
598 
599 
600 // ************************************************************************* //
Foam::surfaceFields.
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
A class for handling file names.
Definition: fileName.H:72
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
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians...
Definition: cylindricalCS.H:67
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:153
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.
bool equal(const T &a, const T &b)
Compare two values for equality.
Definition: label.H:164
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:675
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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:103
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:374
static void addOptionCompat(const word &optName, std::pair< const char *, int > compat)
Specify an alias for the option name.
Definition: argList.C:418
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
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:421
Spatial transformation functions for primitive fields.
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
word findInstance(const fileName &directory, const word &name=word::null, IOobjectOption::readOption rOpt=IOobjectOption::MUST_READ, const word &stopInstance=word::null, const bool constant_fallback=true) const
Return time instance (location) of directory containing the file name (eg, used in reading mesh data)...
Definition: Time.C:725
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:130
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:385
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:902
static unsigned int minPrecision(unsigned int prec) noexcept
Set the minimum default precision.
Definition: IOstream.H:440
ITstream lookup(const word &optName) const
Return an input stream from the named option.
Definition: argListI.H:177
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1113
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
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 a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
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.
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.
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
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:78
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:841
Required Classes.
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:186
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
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:180
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:2185
Do not request registration (bool: false)
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:127