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-2022 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 "regionProperties.H"
72 #include "volFields.H"
73 #include "surfaceFields.H"
74 #include "ReadFields.H"
75 #include "pointFields.H"
76 #include "pointSet.H"
77 #include "transformField.H"
79 #include "unitConversion.H"
80 
81 using namespace Foam;
82 
83 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85 template<class GeoField>
86 void readAndRotateFields
87 (
88  PtrList<GeoField>& flds,
89  const fvMesh& mesh,
90  const tensor& T,
91  const IOobjectList& objects
92 )
93 {
94  ReadFields(mesh, objects, flds);
95  forAll(flds, i)
96  {
97  Info<< "Transforming " << flds[i].name() << endl;
98  dimensionedTensor dimT("t", flds[i].dimensions(), T);
99  transform(flds[i], dimT, flds[i]);
100  }
101 }
102 
103 
104 void rotateFields(const argList& args, const Time& runTime, const tensor& T)
105 {
106  #include "createNamedMesh.H"
107 
108  // Read objects in time directory
109  IOobjectList objects(mesh, runTime.timeName());
110 
111  // Read vol fields.
113  readAndRotateFields(vsFlds, mesh, T, objects);
115  readAndRotateFields(vvFlds, mesh, T, objects);
117  readAndRotateFields(vstFlds, mesh, T, objects);
119  readAndRotateFields(vsymtFlds, mesh, T, objects);
121  readAndRotateFields(vtFlds, mesh, T, objects);
122 
123  // Read surface fields.
125  readAndRotateFields(ssFlds, mesh, T, objects);
127  readAndRotateFields(svFlds, mesh, T, objects);
129  readAndRotateFields(sstFlds, mesh, T, objects);
131  readAndRotateFields(ssymtFlds, mesh, T, objects);
133  readAndRotateFields(stFlds, mesh, T, objects);
134 
135  mesh.write();
136 }
137 
138 
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 
141 int main(int argc, char *argv[])
142 {
144  (
145  "Transforms a mesh by translation, rotation and/or scaling.\n"
146  "The <transformations> are listed comma-separated in a string "
147  "and executed in sequence.\n\n"
148  "transformations:\n"
149  " translate=<vector> "
150  "translation by vector, e.g. (1 2 3)\n"
151  " rotate=(<n1> <n2>) "
152  "rotation from unit vector n1 to n2\n"
153  " Rx=<angle> "
154  "rotation by given angle [deg], e.g. 90, about x-axis\n"
155  " Ry=<angle> "
156  "rotation by given angle [deg] about y-axis\n"
157  " Rz=<angle> "
158  "rotation by given angle [deg] about z-axis\n"
159  " Ra=<axis vector> <angle> "
160  "rotation by given angle [deg] about specified axis\n"
161  " scale=<vector> "
162  "scale by factors from vector in x, y, z directions,\n"
163  " "
164  "e.g. (0.001 0.001 0.001) to scale from mm to m\n\n"
165  "example:\n"
166  " transformPoints "
167  "\"translate=(1.2 0 0), Rx=90, translate=(-1.2 0 0)\""
168  );
169 
170  argList::validArgs.append("transformations");
171 
173  (
174  "rotateFields",
175  "transform vector and tensor fields"
176  );
177 
179  (
180  "pointSet",
181  "pointSet",
182  "Point set to limit the transformation to"
183  );
184 
185  #include "addRegionOption.H"
186  #include "addAllRegionsOption.H"
187  #include "setRootCase.H"
188  #include "createTime.H"
189 
190  const string transformationString(args[1]);
191 
192  #include "createTransforms.H"
193 
194  const wordList regionNames(selectRegionNames(args, runTime));
195 
196  const bool doRotateFields = args.optionFound("rotateFields");
197 
198  word pointSetName = word::null;
199  const bool doPointSet = args.optionReadIfPresent("pointSet", pointSetName);
200 
201  if (doRotateFields && doPointSet)
202  {
204  << "Rotation of fields across the entire mesh, and limiting the "
205  << "transformation of points to a set, cannot be done "
206  << "simultaneously" << exit(FatalError);
207  }
208 
209  forAll(regionNames, regioni)
210  {
211  const word& regionName = regionNames[regioni];
212  const word& regionDir = Foam::regionDir(regionName);
213 
214  fileName meshDir(regionDir/polyMesh::meshSubDir);
215 
217  (
218  IOobject
219  (
220  "points",
221  runTime.findInstance(meshDir, "points"),
222  meshDir,
223  runTime,
226  false
227  )
228  );
229 
230  if (doPointSet)
231  {
232  const labelList setPointIDs
233  (
234  pointSet
235  (
236  IOobject
237  (
238  pointSetName,
239  runTime.findInstance(meshDir/"sets", word::null),
240  polyMesh::meshSubDir/"sets",
241  runTime,
244  false
245  )
246  ).toc()
247  );
248 
249  pointField setPoints(UIndirectList<point>(points, setPointIDs));
250 
251  transforms.transformPosition(setPoints, setPoints);
252 
253  UIndirectList<point>(points, setPointIDs) = setPoints;
254  }
255  else
256  {
258  }
259 
260  if (doRotateFields)
261  {
262  rotateFields(args, runTime, transforms.T());
263  }
264 
265  // Set the precision of the points data to 10
267 
268  Info<< "Writing points into directory " << points.path() << nl << endl;
269  points.write();
270  }
271 
272  Info<< "End\n" << endl;
273 
274  return 0;
275 }
276 
277 
278 // ************************************************************************* //
Foam::surfaceFields.
PtrList< surfaceSphericalTensorField > sstFlds
transformer transforms
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
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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
A set of point labels.
Definition: pointSet.H:48
error FatalError
Spatial transformation functions for FieldFields.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Unit conversion functions.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:328
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:458
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
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:689
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
Field reading functions for post-processing utilities.
Generic dimensioned Type class.
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1658
PtrList< volVectorField > vvFlds
Definition: readVolFields.H:6
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
PtrList< volScalarField > vsFlds
Definition: readVolFields.H:3
Spatial transformation functions for primitive fields.
const word & regionDir(const word &regionName)
const pointField & points
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
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:128
static const word null
An empty word.
Definition: word.H:77
PtrList< volSymmTensorField > vsymtFlds
Definition: readVolFields.H:12
PtrList< surfaceTensorField > stFlds
const tensor & T() const
Return the transformation tensor.
Definition: transformerI.H:94
static const char nl
Definition: Ostream.H:260
objects
PtrList< surfaceScalarField > ssFlds
PtrList< surfaceVectorField > svFlds
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:95
A List with indirect addressing.
Definition: fvMatrix.H:106
messageStream Info
vector transformPosition(const vector &v) const
Transform the given position.
Definition: transformerI.H:153
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:202
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
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:483