initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / mesh / manipulation / transformPoints / transformPoints.C
blobe77fbaf022409d6c7ae4b76a7d3eaff1e98a7937
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     transformPoints
28 Description
29     Transforms the mesh points in the polyMesh directory according to the
30     translate, rotate and scale options.
32 Usage
33     Options are:
35     -translate vector
36         Translates the points by the given vector,
38     -rotate (vector vector)
39         Rotates the points from the first vector to the second,
41      or -yawPitchRoll (yawdegrees pitchdegrees rolldegrees)
42      or -rollPitchYaw (rolldegrees pitchdegrees yawdegrees)
44     -scale vector
45         Scales the points by the given vector.
47     The any or all of the three options may be specified and are processed
48     in the above order.
50     With -rotateFields (in combination with -rotate/yawPitchRoll/rollPitchYaw)
51     it will also read & transform vector & tensor fields.
53     Note:
54     yaw (rotation about z)
55     pitch (rotation about y)
56     roll (rotation about x)
58 \*---------------------------------------------------------------------------*/
60 #include "argList.H"
61 #include "Time.H"
62 #include "fvMesh.H"
63 #include "volFields.H"
64 #include "surfaceFields.H"
65 #include "ReadFields.H"
66 #include "pointFields.H"
67 #include "transformField.H"
68 #include "transformGeometricField.H"
69 #include "IStringStream.H"
71 using namespace Foam;
72 using namespace Foam::mathematicalConstant;
74 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76 template<class GeoField>
77 void readAndRotateFields
79     PtrList<GeoField>& flds,
80     const fvMesh& mesh,
81     const tensor& T,
82     const IOobjectList& objects
85     ReadFields(mesh, objects, flds);
86     forAll(flds, i)
87     {
88         Info<< "Transforming " << flds[i].name() << endl;
89         dimensionedTensor dimT("t", flds[i].dimensions(), T);
90         transform(flds[i], dimT, flds[i]);
91     }
95 void rotateFields(const Time& runTime, const tensor& T)
97 #   include "createMesh.H"
99     // Read objects in time directory
100     IOobjectList objects(mesh, runTime.timeName());
102     // Read vol fields.
104     PtrList<volScalarField> vsFlds;
105     readAndRotateFields(vsFlds, mesh, T, objects);
107     PtrList<volVectorField> vvFlds;
108     readAndRotateFields(vvFlds, mesh, T, objects);
110     PtrList<volSphericalTensorField> vstFlds;
111     readAndRotateFields(vstFlds, mesh, T, objects);
113     PtrList<volSymmTensorField> vsymtFlds;
114     readAndRotateFields(vsymtFlds, mesh, T, objects);
116     PtrList<volTensorField> vtFlds;
117     readAndRotateFields(vtFlds, mesh, T, objects);
119     // Read surface fields.
121     PtrList<surfaceScalarField> ssFlds;
122     readAndRotateFields(ssFlds, mesh, T, objects);
124     PtrList<surfaceVectorField> svFlds;
125     readAndRotateFields(svFlds, mesh, T, objects);
127     PtrList<surfaceSphericalTensorField> sstFlds;
128     readAndRotateFields(sstFlds, mesh, T, objects);
130     PtrList<surfaceSymmTensorField> ssymtFlds;
131     readAndRotateFields(ssymtFlds, mesh, T, objects);
133     PtrList<surfaceTensorField> stFlds;
134     readAndRotateFields(stFlds, mesh, T, objects);
136     mesh.write();
140 //  Main program:
142 int main(int argc, char *argv[])
144     argList::validOptions.insert("translate", "vector");
145     argList::validOptions.insert("rotate", "(vector vector)");
146     argList::validOptions.insert("rollPitchYaw", "(roll pitch yaw)");
147     argList::validOptions.insert("yawPitchRoll", "(yaw pitch roll)");
148     argList::validOptions.insert("rotateFields", "");
149     argList::validOptions.insert("scale", "vector");
151 #   include "setRootCase.H"
152 #   include "createTime.H"
154     pointIOField points
155     (
156         IOobject
157         (
158             "points",
159             runTime.findInstance(polyMesh::meshSubDir, "points"),
160             polyMesh::meshSubDir,
161             runTime,
162             IOobject::MUST_READ,
163             IOobject::NO_WRITE,
164             false
165         )
166     );
169     if (args.options().empty())
170     {
171         FatalErrorIn(args.executable())
172             << "No options supplied, please use one or more of "
173                "-translate, -rotate or -scale options."
174             << exit(FatalError);
175     }
177     if (args.optionFound("translate"))
178     {
179         vector transVector(args.optionLookup("translate")());
181         Info<< "Translating points by " << transVector << endl;
183         points += transVector;
184     }
186     if (args.optionFound("rotate"))
187     {
188         Pair<vector> n1n2(args.optionLookup("rotate")());
189         n1n2[0] /= mag(n1n2[0]);
190         n1n2[1] /= mag(n1n2[1]);
191         tensor T = rotationTensor(n1n2[0], n1n2[1]);
193         Info<< "Rotating points by " << T << endl;
195         points = transform(T, points);
197         if (args.optionFound("rotateFields"))
198         {
199             rotateFields(runTime, T);
200         }
201     }
202     else if (args.optionFound("rollPitchYaw"))
203     {
204         vector v(args.optionLookup("rollPitchYaw")());
206         Info<< "Rotating points by" << nl
207             << "    roll  " << v.x() << nl
208             << "    pitch " << v.y() << nl
209             << "    yaw   " << v.z() << endl;
212         // Convert to radians
213         v *= pi/180.0;
215         quaternion R(v.x(), v.y(), v.z());
217         Info<< "Rotating points by quaternion " << R << endl;
218         points = transform(R, points);
220         if (args.optionFound("rotateFields"))
221         {
222             rotateFields(runTime, R.R());
223         }
224     }
225     else if (args.optionFound("yawPitchRoll"))
226     {
227         vector v(args.optionLookup("yawPitchRoll")());
229         Info<< "Rotating points by" << nl
230             << "    yaw   " << v.x() << nl
231             << "    pitch " << v.y() << nl
232             << "    roll  " << v.z() << endl;
235         // Convert to radians
236         v *= pi/180.0;
238         scalar yaw = v.x();
239         scalar pitch = v.y();
240         scalar roll = v.z();
242         quaternion R = quaternion(vector(0, 0, 1), yaw);
243         R *= quaternion(vector(0, 1, 0), pitch);
244         R *= quaternion(vector(1, 0, 0), roll);
246         Info<< "Rotating points by quaternion " << R << endl;
247         points = transform(R, points);
249         if (args.optionFound("rotateFields"))
250         {
251             rotateFields(runTime, R.R());
252         }
253     }
255     if (args.optionFound("scale"))
256     {
257         vector scaleVector(args.optionLookup("scale")());
259         Info<< "Scaling points by " << scaleVector << endl;
261         points.replace(vector::X, scaleVector.x()*points.component(vector::X));
262         points.replace(vector::Y, scaleVector.y()*points.component(vector::Y));
263         points.replace(vector::Z, scaleVector.z()*points.component(vector::Z));
264     }
266     // Set the precision of the points data to 10
267     IOstream::defaultPrecision(10);
269     Info << "Writing points into directory " << points.path() << nl << endl;
270     points.write();
272     return 0;
276 // ************************************************************************* //