rotateMesh.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) 2019-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Application
28 rotateMesh
29
30Group
31 grpMeshManipulationUtilities
32
33Description
34 Rotates the mesh and fields from the direction n1 to direction n2.
35
36\*---------------------------------------------------------------------------*/
37
38#include "argList.H"
39#include "timeSelector.H"
40#include "Time.H"
41#include "fvMesh.H"
42#include "volFields.H"
43#include "surfaceFields.H"
44#include "regionProperties.H"
46#include "IOobjectList.H"
47
48using namespace Foam;
49
50template<class GeoField>
51void ReadAndRotateFields
52(
53 const fvMesh& mesh,
54 const IOobjectList& objects,
55 const dimensionedTensor& rotT
56)
57{
58 // Objects of field type
59 IOobjectList fields(objects.lookupClass<GeoField>());
60
61 forAllConstIters(fields, fieldIter)
62 {
63 GeoField fld(*fieldIter(), mesh);
64 Info<< " Rotating " << fld.name() << endl;
65 transform(fld, rotT, fld);
66 fld.write();
67 }
68}
69
70
71void rotateFields
72(
73 const fvMesh& mesh,
74 const Time& runTime,
75 const tensor& rotationT
76)
77{
78 // Need dimensionedTensor for geometric fields
79 const dimensionedTensor rotT(rotationT);
80
81 // Search for list of objects for this time
82 IOobjectList objects(mesh, runTime.timeName());
83
84 ReadAndRotateFields<volVectorField>(mesh, objects, rotT);
85 ReadAndRotateFields<volSphericalTensorField>(mesh, objects, rotT);
86 ReadAndRotateFields<volSymmTensorField>(mesh, objects, rotT);
87 ReadAndRotateFields<volTensorField>(mesh, objects, rotT);
88
89 ReadAndRotateFields<surfaceVectorField>(mesh, objects, rotT);
90 ReadAndRotateFields<surfaceSphericalTensorField>(mesh, objects, rotT);
91 ReadAndRotateFields<surfaceSymmTensorField>(mesh, objects, rotT);
92 ReadAndRotateFields<surfaceTensorField>(mesh, objects, rotT);
93}
94
95
96// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97
98int main(int argc, char *argv[])
99{
100 argList::addNote
101 (
102 "Rotate mesh points and vector/tensor fields\n"
103 "Rotation from the <from> vector to the <to> vector"
104 );
105
106 timeSelector::addOptions();
107
108 argList::addArgument("from", "The vector to rotate from");
109 argList::addArgument("to", "The vector to rotate to");
110
111 #include "addAllRegionOptions.H"
112 #include "setRootCase.H"
113
114 const vector n1(args.get<vector>(1).normalise());
115 const vector n2(args.get<vector>(2).normalise());
116
117 const tensor rotT(rotationTensor(n1, n2));
118
119 // ------------------------------------------------------------------------
120
121 #include "createTime.H"
122
123 // Handle -allRegions, -regions, -region
124 #include "getAllRegionOptions.H"
125
126 // ------------------------------------------------------------------------
127
128 #include "createNamedMeshes.H"
129
130 forAll(regionNames, regioni)
131 {
132 const word& regionName = regionNames[regioni];
133 const fileName meshDir
134 (
135 polyMesh::regionName(regionName)/polyMesh::meshSubDir
136 );
137
138 if (regionNames.size() > 1)
139 {
140 Info<< "region=" << regionName << nl;
141 }
142
144 (
146 (
147 "points",
148 runTime.findInstance(meshDir, "points"),
149 meshDir,
150 runTime,
151 IOobject::MUST_READ,
152 IOobject::NO_WRITE,
153 false
154 )
155 );
156
157 points = transform(rotT, points);
158
159 // Set the precision of the points data to 10
160 IOstream::defaultPrecision(max(10u, IOstream::defaultPrecision()));
161
162 Info<< "Writing points into directory "
163 << runTime.relativePath(points.path()) << nl
164 << endl;
165 points.write();
166 }
167
168 instantList timeDirs = timeSelector::select0(runTime, args);
169
170 forAll(timeDirs, timei)
171 {
172 runTime.setTime(timeDirs[timei], timei);
173
174 Info<< "Time = " << runTime.timeName() << endl;
175
176 forAll(regionNames, regioni)
177 {
178 const word& regionName = regionNames[regioni];
179 if (regionNames.size() > 1)
180 {
181 Info<< "region=" << regionName << nl;
182 }
183
184 rotateFields(meshes[regioni], runTime, rotT);
185 }
186 }
187
188 Info<< "\nEnd\n" << endl;
189
190 return 0;
191}
192
193
194// ************************************************************************* //
Y[inertIndex] max(0.0)
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
A primitive field of type <T> with automated input and output.
Definition: IOField.H:58
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:59
IOobjectList lookupClass(const char *clsName) const
The list of IOobjects with the given headerClassName.
Definition: IOobjectList.C:313
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:278
Generic dimensioned Type class.
A class for handling file names.
Definition: fileName.H:76
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
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:68
dynamicFvMesh & mesh
engineTime & runTime
Foam::word regionName(Foam::polyMesh::defaultRegion)
Required Variables.
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
wordList regionNames
const pointField & points
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
Foam::argList args(argc, argv)
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Foam::surfaceFields.
Spatial transformation functions for GeometricField.