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-2024 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  - \par -pointSet <name> \n
55  Only transform points in the given point set.
56 
57  Example usage:
58  transformPoints \
59  "translate=(-0.05 -0.05 0), \
60  Rz=45, \
61  translate=(0.05 0.05 0)"
62 
63 See also
64  Foam::transformer
65  surfaceTransformPoints
66 
67 \*---------------------------------------------------------------------------*/
68 
69 #include "argList.H"
70 #include "fvMesh.H"
71 #include "volFields.H"
72 #include "surfaceFields.H"
73 #include "ReadFields.H"
74 #include "pointFields.H"
75 #include "pointSet.H"
76 #include "transformField.H"
78 
79 using namespace Foam;
80 
81 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
83 template<class GeoField>
84 void readAndRotateFields
85 (
86  PtrList<GeoField>& flds,
87  const fvMesh& mesh,
88  const tensor& T,
89  const IOobjectList& objects
90 )
91 {
92  ReadFields(mesh, objects, flds);
93  forAll(flds, i)
94  {
95  Info<< "Transforming " << flds[i].name() << endl;
96  dimensionedTensor dimT("t", flds[i].dimensions(), T);
97  transform(flds[i], dimT, flds[i]);
98  }
99 }
100 
101 
102 void rotateFields(const argList& args, const Time& runTime, const tensor& T)
103 {
105 
106  // Read objects in time directory
107  IOobjectList objects(mesh, runTime.name());
108 
109  // Read vol fields.
111  readAndRotateFields(vsFlds, mesh, T, objects);
113  readAndRotateFields(vvFlds, mesh, T, objects);
115  readAndRotateFields(vstFlds, mesh, T, objects);
117  readAndRotateFields(vsymtFlds, mesh, T, objects);
119  readAndRotateFields(vtFlds, mesh, T, objects);
120 
121  // Read surface fields.
123  readAndRotateFields(ssFlds, mesh, T, objects);
125  readAndRotateFields(svFlds, mesh, T, objects);
127  readAndRotateFields(sstFlds, mesh, T, objects);
129  readAndRotateFields(ssymtFlds, mesh, T, objects);
131  readAndRotateFields(stFlds, mesh, T, objects);
132 
133  mesh.write();
134 }
135 
136 
137 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
138 
139 int main(int argc, char *argv[])
140 {
142  (
143  "Transforms a mesh by translation, rotation and/or scaling.\n"
144  "The <transformations> are listed comma-separated in a string "
145  "and executed in sequence.\n\n"
146  "transformations:\n"
147  " translate=<vector> "
148  "translation by vector, e.g. (1 2 3)\n"
149  " rotate=(<n1> <n2>) "
150  "rotation from unit vector n1 to n2\n"
151  " Rx=<angle> "
152  "rotation by given angle [deg], e.g. 90, about x-axis\n"
153  " Ry=<angle> "
154  "rotation by given angle [deg] about y-axis\n"
155  " Rz=<angle> "
156  "rotation by given angle [deg] about z-axis\n"
157  " Ra=<axis vector> <angle> "
158  "rotation by given angle [deg] about specified axis\n"
159  " scale=<vector> "
160  "scale by factors from vector in x, y, z directions,\n"
161  " "
162  "e.g. (0.001 0.001 0.001) to scale from mm to m\n\n"
163  "example:\n"
164  " transformPoints "
165  "\"translate=(1.2 0 0), Rx=90, translate=(-1.2 0 0)\""
166  );
167 
168  argList::validArgs.append("transformations");
169 
171  (
172  "rotateFields",
173  "transform vector and tensor fields"
174  );
175 
177  (
178  "pointSet",
179  "pointSet",
180  "Point set to limit the transformation to"
181  );
182 
183  #include "addRegionOption.H"
184  #include "addAllRegionsOption.H"
185  #include "setRootCase.H"
186  #include "createTime.H"
187 
188  const string transformationString(args[1]);
189 
190  #include "createTransforms.H"
191 
192  #include "setRegionNames.H"
193 
194  const bool doRotateFields = args.optionFound("rotateFields");
195 
196  word pointSetName = word::null;
197  const bool doPointSet = args.optionReadIfPresent("pointSet", pointSetName);
198 
199  if (doRotateFields && doPointSet)
200  {
202  << "Rotation of fields across the entire mesh, and limiting the "
203  << "transformation of points to a set, cannot be done "
204  << "simultaneously" << exit(FatalError);
205  }
206 
207  forAll(regionNames, regioni)
208  {
209  const word& regionName = regionNames[regioni];
210 
211  const word& regionDir =
213  ? word::null
214  : regionName;
215 
216  const fileName meshDir(regionDir/polyMesh::meshSubDir);
217 
219  (
220  IOobject
221  (
222  "points",
223  runTime.findInstance(meshDir, "points"),
224  meshDir,
225  runTime,
228  false
229  )
230  );
231 
232  if (doPointSet)
233  {
234  const labelList setPointIDs
235  (
236  pointSet
237  (
238  IOobject
239  (
240  pointSetName,
241  runTime.findInstance(meshDir/"sets", word::null),
242  polyMesh::meshSubDir/"sets",
243  runTime,
246  false
247  )
248  ).toc()
249  );
250 
251  pointField setPoints(UIndirectList<point>(points, setPointIDs));
252 
253  transforms.transformPosition(setPoints, setPoints);
254 
255  UIndirectList<point>(points, setPointIDs) = setPoints;
256  }
257  else
258  {
260  }
261 
262  if (doRotateFields)
263  {
264  rotateFields(args, runTime, transforms.T());
265  }
266 
267  // Set the precision of the points data to 10
269 
270  Info<< "Writing points into directory " << points.path() << nl << endl;
271  points.write();
272  }
273 
274  Info<< "End\n" << endl;
275 
276  return 0;
277 }
278 
279 
280 // ************************************************************************* //
Field reading functions for post-processing utilities.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:227
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:458
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
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:658
A List with indirect addressing.
Definition: UIndirectList.H:60
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:103
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
Generic dimensioned Type class.
const word & name() const
Return const reference to name.
A class for handling file names.
Definition: fileName.H:82
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1737
A set of point labels.
Definition: pointSet.H:51
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:267
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:270
vector transformPosition(const vector &v) const
Transform the given position.
Definition: transformerI.H:153
const tensor & T() const
Return the transformation tensor.
Definition: transformerI.H:94
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
transformer transforms
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const pointField & points
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
wordList ReadFields(const Mesh &mesh, const IOobjectList &objects, PtrList< GeoField > &fields, const bool syncPar=true)
Read all fields of the specified type.
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
const word & regionName(const solver &region)
Definition: solver.H:209
messageStream Info
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:504
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
PtrList< surfaceScalarField > ssFlds
PtrList< surfaceTensorField > stFlds
PtrList< surfaceSymmTensorField > ssymtFlds
PtrList< surfaceVectorField > svFlds
objects
PtrList< surfaceSphericalTensorField > sstFlds
PtrList< volSphericalTensorField > vstFlds
Definition: readVolFields.H:9
PtrList< volScalarField > vsFlds
Definition: readVolFields.H:3
PtrList< volTensorField > vtFlds
Definition: readVolFields.H:15
PtrList< volSymmTensorField > vsymtFlds
Definition: readVolFields.H:12
PtrList< volVectorField > vvFlds
Definition: readVolFields.H:6
const Foam::wordList regionNames(args.optionFound("allRegions") ? runTime .controlDict().subDict("regionSolvers").toc() :wordList(1, args.optionFound("region") ? args.optionRead< word >("region") :polyMesh::defaultRegion))
Foam::argList args(argc, argv)
Foam::surfaceFields.
Spatial transformation functions for primitive fields.
Spatial transformation functions for FieldFields.