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 | Copyright (C) 2011-2015 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  rotateMesh
26 
27 Description
28  Rotates the mesh and fields from the direction n1 to direction n2.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "argList.H"
33 #include "timeSelector.H"
34 #include "Time.H"
35 #include "fvMesh.H"
36 #include "volFields.H"
37 #include "surfaceFields.H"
39 #include "IOobjectList.H"
40 
41 using namespace Foam;
42 
43 template<class GeometricField>
44 void RotateFields
45 (
46  const fvMesh& mesh,
47  const IOobjectList& objects,
48  const tensor& T
49 )
50 {
51  // Search list of objects for volScalarFields
52  IOobjectList fields(objects.lookupClass(GeometricField::typeName));
53 
54  forAllIter(IOobjectList, fields, fieldIter)
55  {
56  Info<< " Rotating " << fieldIter()->name() << endl;
57 
58  GeometricField theta(*fieldIter(), mesh);
59  transform(theta, dimensionedTensor(T), theta);
60  theta.write();
61  }
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 int main(int argc, char *argv[])
68 {
70 
71  argList::validArgs.append("n1");
72  argList::validArgs.append("n2");
73 
74  #include "setRootCase.H"
75  #include "createTime.H"
76 
77  vector n1(args.argRead<vector>(1));
78  n1 /= mag(n1);
79 
80  vector n2(args.argRead<vector>(2));
81  n2 /= mag(n2);
82 
83  tensor T(rotationTensor(n1, n2));
84 
85  {
87  (
88  IOobject
89  (
90  "points",
91  runTime.findInstance(polyMesh::meshSubDir, "points"),
93  runTime,
96  false
97  )
98  );
99 
100  points = transform(T, points);
101 
102  // Set the precision of the points data to 10
104 
105  Info<< "Writing points into directory " << points.path() << nl << endl;
106  points.write();
107  }
108 
109 
111 
112  #include "createMesh.H"
113 
114  forAll(timeDirs, timeI)
115  {
116  runTime.setTime(timeDirs[timeI], timeI);
117 
118  Info<< "Time = " << runTime.timeName() << endl;
119 
120  // Search for list of objects for this time
121  IOobjectList objects(mesh, runTime.timeName());
122 
123  RotateFields<volVectorField>(mesh, objects, T);
124  RotateFields<volSphericalTensorField>(mesh, objects, T);
125  RotateFields<volSymmTensorField>(mesh, objects, T);
126  RotateFields<volTensorField>(mesh, objects, T);
127 
128  RotateFields<surfaceVectorField>(mesh, objects, T);
129  RotateFields<surfaceSphericalTensorField>(mesh, objects, T);
130  RotateFields<surfaceSymmTensorField>(mesh, objects, T);
131  RotateFields<surfaceTensorField>(mesh, objects, T);
132  }
133 
134  Info<< "\nEnd\n" << endl;
135 
136  return 0;
137 }
138 
139 
140 // ************************************************************************* //
Foam::argList::validArgs
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:143
volFields.H
Foam::Tensor
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
Definition: complexI.H:224
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::IOField< vector >
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
fields
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p - rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:127
transformGeometricField.H
Spatial transformation functions for FieldFields.
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
IOobjectList.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::transform
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
argList.H
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::IOobjectList::lookupClass
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:199
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
setRootCase.H
timeDirs
static instantList timeDirs
Definition: globalFoam.H:44
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:461
Foam::timeSelector::addOptions
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
createMesh.H
points
const pointField & points
Definition: gmvOutputHeader.H:1
timeSelector.H
createTime.H
Foam::argList::argRead
T argRead(const label index) const
Read a value from the argument at index.
Definition: argListI.H:177
Foam::rotationTensor
tensor rotationTensor(const vector &n1, const vector &n2)
Definition: transform.H:46
Foam::timeSelector::select0
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:253
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
args
Foam::argList args(argc, argv)
Foam::dimensionedTensor
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
Definition: dimensionedTensor.H:49