Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
12License
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
28Application
29 transformPoints
30
31Group
32 grpMeshManipulationUtilities
33
34Description
35 Transforms the mesh points in the polyMesh directory according to the
36 translate, rotate and scale options.
37
38Usage
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
81Note
82 roll (rotation about x)
83 pitch (rotation about y)
84 yaw (rotation about z)
85
86 - with -rotate and two exactly opposing vectors it will actually mirror
87 the geometry. Use any of the other rotation options instead or use
88 two steps, each 90 degrees.
89
90\*---------------------------------------------------------------------------*/
91
92#include "argList.H"
93#include "Time.H"
94#include "fvMesh.H"
95#include "volFields.H"
96#include "surfaceFields.H"
97#include "pointFields.H"
98#include "ReadFields.H"
99#include "regionProperties.H"
100#include "transformField.H"
102#include "axisAngleRotation.H"
104#include "cylindricalCS.H"
105
106using namespace Foam;
107using namespace Foam::coordinateRotations;
108
109// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
110
111template<class GeoField>
112void readAndRotateFields
113(
114 PtrList<GeoField>& flds,
115 const fvMesh& mesh,
116 const dimensionedTensor& rotT,
117 const IOobjectList& objects
118)
119{
120 ReadFields(mesh, objects, flds);
121 for (GeoField& fld : flds)
122 {
123 Info<< "Transforming " << fld.name() << endl;
124 transform(fld, rotT, fld);
125 }
126}
127
128
129void rotateFields
130(
131 const word& regionName,
132 const Time& runTime,
133 const tensor& rotT
134)
135{
136 // Need dimensionedTensor for geometric fields
137 const dimensionedTensor T(rotT);
138
139 #include "createRegionMesh.H"
140
141 // Read objects in time directory
142 IOobjectList objects(mesh, runTime.timeName());
143
144 // Read vol fields.
145
147 readAndRotateFields(vsFlds, mesh, T, objects);
148
150 readAndRotateFields(vvFlds, mesh, T, objects);
151
153 readAndRotateFields(vstFlds, mesh, T, objects);
154
156 readAndRotateFields(vsymtFlds, mesh, T, objects);
157
159 readAndRotateFields(vtFlds, mesh, T, objects);
160
161 // Read surface fields.
162
164 readAndRotateFields(ssFlds, mesh, T, objects);
165
167 readAndRotateFields(svFlds, mesh, T, objects);
168
170 readAndRotateFields(sstFlds, mesh, T, objects);
171
173 readAndRotateFields(ssymtFlds, mesh, T, objects);
174
176 readAndRotateFields(stFlds, mesh, T, objects);
177
178 mesh.write();
179}
180
181
182// Retrieve scaling option
183// - size 0 : no scaling
184// - size 1 : uniform scaling
185// - size 3 : non-uniform scaling
186List<scalar> getScalingOpt(const word& optName, const argList& args)
187{
188 // readListIfPresent handles single or multiple values
189 // - accept 1 or 3 values
190
191 List<scalar> scaling;
192 args.readListIfPresent(optName, scaling);
193
194 if (scaling.size() == 1)
195 {
196 // Uniform scaling
197 }
198 else if (scaling.size() == 3)
199 {
200 // Non-uniform, but may actually be uniform
201 if
202 (
203 equal(scaling[0], scaling[1])
204 && equal(scaling[0], scaling[2])
205 )
206 {
207 scaling.resize(1);
208 }
209 }
210 else if (!scaling.empty())
211 {
213 << "Incorrect number of components, must be 1 or 3." << nl
214 << " -" << optName << ' ' << args[optName].c_str() << endl
215 << exit(FatalError);
216 }
217
218 if (scaling.size() == 1 && equal(scaling[0], 1))
219 {
220 // Scale factor 1 == no scaling
221 scaling.clear();
222 }
223
224 // Zero and negative scaling are permitted
225
226 return scaling;
227}
228
229
230void applyScaling(pointField& points, const List<scalar>& scaling)
231{
232 if (scaling.size() == 1)
233 {
234 Info<< "Scaling points uniformly by " << scaling[0] << nl;
235 points *= scaling[0];
236 }
237 else if (scaling.size() == 3)
238 {
239 const vector factor(scaling[0], scaling[1], scaling[2]);
240 Info<< "Scaling points by " << factor << nl;
241 cmptMultiply(points, points, factor);
242 }
243}
244
245
246// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247
248int main(int argc, char *argv[])
249{
251 (
252 "Transform (translate / rotate / scale) mesh points.\n"
253 "Note: roll=rotate about x, pitch=rotate about y, yaw=rotate about z"
254 );
256 (
257 "time",
258 "time",
259 "Specify the time to search from and apply the transformation"
260 " (default is latest)"
261 );
263 (
264 "recentre",
265 "Recentre the bounding box before other operations"
266 );
268 (
269 "translate",
270 "vector",
271 "Translate by specified <vector> before rotations"
272 );
274 (
275 "auto-centre",
276 "Use bounding box centre as centre for rotations"
277 );
279 (
280 "centre",
281 "point",
282 "Use specified <point> as centre for rotations"
283 );
284 argList::addOptionCompat("auto-centre", {"auto-origin", 2206});
285 argList::addOptionCompat("centre", {"origin", 2206});
286
288 (
289 "rotate",
290 "(vectorA vectorB)",
291 "Rotate from <vectorA> to <vectorB> - eg, '((1 0 0) (0 0 1))'"
292 );
294 (
295 "rotate-angle",
296 "(vector angle)",
297 "Rotate <angle> degrees about <vector> - eg, '((1 0 0) 45)'"
298 );
300 (
301 "rotate-x", "deg",
302 "Rotate (degrees) about x-axis"
303 );
305 (
306 "rotate-y", "deg",
307 "Rotate (degrees) about y-axis"
308 );
310 (
311 "rotate-z", "deg",
312 "Rotate (degrees) about z-axis"
313 );
315 (
316 "rollPitchYaw",
317 "vector",
318 "Rotate by '(roll pitch yaw)' degrees"
319 );
321 (
322 "yawPitchRoll",
323 "vector",
324 "Rotate by '(yaw pitch roll)' degrees"
325 );
327 (
328 "rotateFields",
329 "Read and transform vector and tensor fields too"
330 );
332 (
333 "scale",
334 "scalar | vector",
335 "Scale by the specified amount - Eg, for uniform [mm] to [m] scaling "
336 "use either '(0.001 0.001 0.001)' or simply '0.001'"
337 );
339 (
340 "cylToCart",
341 "(originVec axisVec directionVec)",
342 "Tranform cylindrical coordinates to cartesian coordinates"
343 );
344
345
346 // Compatibility with surfaceTransformPoints
347 argList::addOptionCompat("scale", {"write-scale", 0});
348
349 #include "addAllRegionOptions.H"
350 #include "setRootCase.H"
351
352 const bool doRotateFields = args.found("rotateFields");
353
354 // Verify that an operation has been specified
355 {
356 const List<word> operationNames
357 ({
358 "recentre",
359 "translate",
360 "rotate",
361 "rotate-angle",
362 "rotate-x",
363 "rotate-y",
364 "rotate-z",
365 "rollPitchYaw",
366 "yawPitchRoll",
367 "scale",
368 "cylToCart"
369 });
370
371 if (!args.count(operationNames))
372 {
374 << "No operation supplied, "
375 << "use at least one of the following:" << nl
376 << " ";
377
378 for (const auto& opName : operationNames)
379 {
381 << " -" << opName;
382 }
383
385 << nl << exit(FatalError);
386 }
387 }
388
389 // ------------------------------------------------------------------------
390
391 #include "createTime.H"
392
393 if (args.found("time"))
394 {
395 if (args["time"] == "constant")
396 {
397 runTime.setTime(instant(0, "constant"), 0);
398 }
399 else
400 {
401 const scalar timeValue = args.get<scalar>("time");
402 runTime.setTime(instant(timeValue), 0);
403 }
404 }
405
406 // Handle -allRegions, -regions, -region
407 #include "getAllRegionOptions.H"
408
409 // ------------------------------------------------------------------------
410
411 forAll(regionNames, regioni)
412 {
413 const word& regionName = regionNames[regioni];
414 const fileName meshDir
415 (
417 );
418
419 if (regionNames.size() > 1)
420 {
421 Info<< "region=" << regionName << nl;
422 }
423
425 (
427 (
428 "points",
429 runTime.findInstance(meshDir, "points"),
430 meshDir,
431 runTime,
435 )
436 );
437
438
439 // Begin operations
440
441 vector v;
442 if (args.found("recentre"))
443 {
444 v = boundBox(points).centre();
445 Info<< "Adjust centre " << v << " -> (0 0 0)" << endl;
446 points -= v;
447 }
448
449 if (args.readIfPresent("translate", v))
450 {
451 Info<< "Translating points by " << v << endl;
452 points += v;
453 }
454
455 vector rotationCentre;
456 bool useRotationCentre = args.readIfPresent("centre", rotationCentre);
457 if (args.found("auto-centre") && !useRotationCentre)
458 {
459 useRotationCentre = true;
460 rotationCentre = boundBox(points).centre();
461 }
462
463 if (useRotationCentre)
464 {
465 Info<< "Set centre of rotation to " << rotationCentre << endl;
466 points -= rotationCentre;
467 }
468
469
470 // Get a rotation specification
471
472 tensor rot(Zero);
473 bool useRotation(false);
474
475 if (args.found("rotate"))
476 {
477 Pair<vector> n1n2
478 (
479 args.lookup("rotate")()
480 );
481 n1n2[0].normalise();
482 n1n2[1].normalise();
483
484 rot = rotationTensor(n1n2[0], n1n2[1]);
485 useRotation = true;
486 }
487 else if (args.found("rotate-angle"))
488 {
489 const Tuple2<vector, scalar> rotAxisAngle
490 (
491 args.lookup("rotate-angle")()
492 );
493
494 const vector& axis = rotAxisAngle.first();
495 const scalar angle = rotAxisAngle.second();
496
497 Info<< "Rotating points " << nl
498 << " about " << axis << nl
499 << " angle " << angle << nl;
500
501 rot = axisAngle::rotation(axis, angle, true);
502 useRotation = true;
503 }
504 else if (args.found("rotate-x"))
505 {
506 const scalar angle = args.get<scalar>("rotate-x");
507
508 Info<< "Rotating points about x-axis: " << angle << nl;
509
510 rot = axisAngle::rotation(vector::X, angle, true);
511 useRotation = true;
512 }
513 else if (args.found("rotate-y"))
514 {
515 const scalar angle = args.get<scalar>("rotate-y");
516
517 Info<< "Rotating points about y-axis: " << angle << nl;
518
519 rot = axisAngle::rotation(vector::Y, angle, true);
520 useRotation = true;
521 }
522 else if (args.found("rotate-z"))
523 {
524 const scalar angle = args.get<scalar>("rotate-z");
525
526 Info<< "Rotating points about z-axis: " << angle << nl;
527
528 rot = axisAngle::rotation(vector::Z, angle, true);
529 useRotation = true;
530 }
531 else if (args.readIfPresent("rollPitchYaw", v))
532 {
533 Info<< "Rotating points by" << nl
534 << " roll " << v.x() << nl
535 << " pitch " << v.y() << nl
536 << " yaw " << v.z() << nl;
537
538 rot = euler::rotation(euler::eulerOrder::ROLL_PITCH_YAW, v, true);
539 useRotation = true;
540 }
541 else if (args.readIfPresent("yawPitchRoll", v))
542 {
543 Info<< "Rotating points by" << nl
544 << " yaw " << v.x() << nl
545 << " pitch " << v.y() << nl
546 << " roll " << v.z() << nl;
547
548 rot = euler::rotation(euler::eulerOrder::YAW_PITCH_ROLL, v, true);
549 useRotation = true;
550 }
551
552 if (useRotation)
553 {
554 Info<< "Rotating points by " << rot << endl;
555 transform(points, rot, points);
556
557 if (doRotateFields)
558 {
559 rotateFields(regionName, runTime, rot);
560 }
561 }
562
563 if (useRotationCentre)
564 {
565 Info<< "Unset centre of rotation from " << rotationCentre << endl;
566 points += rotationCentre;
567 }
568
569 // Output scaling
570 applyScaling(points, getScalingOpt("scale", args));
571
572 if (args.found("cylToCart"))
573 {
574 vectorField n1n2(args.lookup("cylToCart")());
575 n1n2[1].normalise();
576 n1n2[2].normalise();
577
578 cylindricalCS ccs
579 (
580 "ccs",
581 n1n2[0],
582 n1n2[1],
583 n1n2[2]
584 );
585
586 points = ccs.globalPosition(points);
587 }
588
589 // More precision (for points data)
591
592 Info<< "Writing points into directory "
593 << runTime.relativePath(points.path()) << nl
594 << endl;
595 points.write();
596 }
597
598 Info<< nl << "End" << nl << endl;
599
600 return 0;
601}
602
603
604// ************************************************************************* //
Field reading functions for post-processing utilities.
Required Classes.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static unsigned int minPrecision(unsigned int prec) noexcept
Set the minimum default precision.
Definition IOstream.H:459
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition Tuple2.H:51
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Extract command arguments and options from the supplied argc and argv parameters.
Definition argList.H:119
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition argList.C:389
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:400
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
static void addOptionCompat(const word &optName, std::pair< const char *, int > compat)
Specify an alias for the option name.
Definition argList.C:433
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
point centre() const
The centre (midpoint) of the bounding box.
Definition boundBoxI.H:186
static tensor rotation(const vector &axis, const scalar angle, bool degrees=false)
The rotation tensor for given axis/angle.
static tensor rotation(const vector &angles, bool degrees=false)
Rotation tensor calculated for the intrinsic Euler angles in z-x-z order.
A class for handling file names.
Definition fileName.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name.
Definition instant.H:56
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir).
Definition polyMesh.C:826
Tensor of scalars, i.e. Tensor<scalar>.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
const volScalarField & T
dynamicFvMesh & mesh
engineTime & runTime
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
Required Variables.
wordList regionNames
const pointField & points
Namespace for coordinate system rotations.
Namespace for OpenFOAM.
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.
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition transform.H:47
bool equal(const T &a, const T &b)
Compare two values for equality.
Definition label.H:180
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
messageStream Info
Information stream (stdout output on master, null elsewhere).
vectorIOField pointIOField
pointIOField is a vectorIOField.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
coordSystem::cylindrical cylindricalCS
Compatibility typedef 1806.
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.
Spatial transformation functions for primitive fields.
Spatial transformation functions for GeometricField.