00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089 #include <OpenFOAM/argList.H>
00090 #include <OpenFOAM/Time.H>
00091 #include <finiteVolume/fvMesh.H>
00092 #include <finiteVolume/volFields.H>
00093 #include <finiteVolume/surfaceFields.H>
00094 #include <OpenFOAM/ReadFields.H>
00095 #include <OpenFOAM/pointFields.H>
00096 #include <OpenFOAM/transformField.H>
00097 #include <OpenFOAM/transformGeometricField.H>
00098 #include <OpenFOAM/IStringStream.H>
00099
00100 using namespace Foam;
00101 using namespace Foam::mathematicalConstant;
00102
00103
00104
00105 template<class GeoField>
00106 void readAndRotateFields
00107 (
00108 PtrList<GeoField>& flds,
00109 const fvMesh& mesh,
00110 const tensor& T,
00111 const IOobjectList& objects
00112 )
00113 {
00114 ReadFields(mesh, objects, flds);
00115 forAll(flds, i)
00116 {
00117 Info<< "Transforming " << flds[i].name() << endl;
00118 dimensionedTensor dimT("t", flds[i].dimensions(), T);
00119 transform(flds[i], dimT, flds[i]);
00120 }
00121 }
00122
00123
00124 void rotateFields(const argList& args, const Time& runTime, const tensor& T)
00125 {
00126 # include <OpenFOAM/createNamedMesh.H>
00127
00128
00129 IOobjectList objects(mesh, runTime.timeName());
00130
00131
00132
00133 PtrList<volScalarField> vsFlds;
00134 readAndRotateFields(vsFlds, mesh, T, objects);
00135
00136 PtrList<volVectorField> vvFlds;
00137 readAndRotateFields(vvFlds, mesh, T, objects);
00138
00139 PtrList<volSphericalTensorField> vstFlds;
00140 readAndRotateFields(vstFlds, mesh, T, objects);
00141
00142 PtrList<volSymmTensorField> vsymtFlds;
00143 readAndRotateFields(vsymtFlds, mesh, T, objects);
00144
00145 PtrList<volTensorField> vtFlds;
00146 readAndRotateFields(vtFlds, mesh, T, objects);
00147
00148
00149
00150 PtrList<surfaceScalarField> ssFlds;
00151 readAndRotateFields(ssFlds, mesh, T, objects);
00152
00153 PtrList<surfaceVectorField> svFlds;
00154 readAndRotateFields(svFlds, mesh, T, objects);
00155
00156 PtrList<surfaceSphericalTensorField> sstFlds;
00157 readAndRotateFields(sstFlds, mesh, T, objects);
00158
00159 PtrList<surfaceSymmTensorField> ssymtFlds;
00160 readAndRotateFields(ssymtFlds, mesh, T, objects);
00161
00162 PtrList<surfaceTensorField> stFlds;
00163 readAndRotateFields(stFlds, mesh, T, objects);
00164
00165 mesh.write();
00166 }
00167
00168
00169
00170
00171 int main(int argc, char *argv[])
00172 {
00173 # include <OpenFOAM/addRegionOption.H>
00174 argList::validOptions.insert("translate", "vector");
00175 argList::validOptions.insert("rotate", "(vector vector)");
00176 argList::validOptions.insert("rollPitchYaw", "(roll pitch yaw)");
00177 argList::validOptions.insert("yawPitchRoll", "(yaw pitch roll)");
00178 argList::validOptions.insert("rotateFields", "");
00179 argList::validOptions.insert("scale", "vector");
00180
00181 # include <OpenFOAM/setRootCase.H>
00182 # include <OpenFOAM/createTime.H>
00183
00184 word regionName = polyMesh::defaultRegion;
00185 fileName meshDir;
00186
00187 if (args.optionReadIfPresent("region", regionName))
00188 {
00189 meshDir = regionName/polyMesh::meshSubDir;
00190 }
00191 else
00192 {
00193 meshDir = polyMesh::meshSubDir;
00194 }
00195
00196 pointIOField points
00197 (
00198 IOobject
00199 (
00200 "points",
00201 runTime.findInstance(meshDir, "points"),
00202 meshDir,
00203 runTime,
00204 IOobject::MUST_READ,
00205 IOobject::NO_WRITE,
00206 false
00207 )
00208 );
00209
00210
00211 if (args.options().empty())
00212 {
00213 FatalErrorIn(args.executable())
00214 << "No options supplied, please use one or more of "
00215 "-translate, -rotate or -scale options."
00216 << exit(FatalError);
00217 }
00218
00219 if (args.optionFound("translate"))
00220 {
00221 vector transVector(args.optionLookup("translate")());
00222
00223 Info<< "Translating points by " << transVector << endl;
00224
00225 points += transVector;
00226 }
00227
00228 if (args.optionFound("rotate"))
00229 {
00230 Pair<vector> n1n2(args.optionLookup("rotate")());
00231 n1n2[0] /= mag(n1n2[0]);
00232 n1n2[1] /= mag(n1n2[1]);
00233 tensor T = rotationTensor(n1n2[0], n1n2[1]);
00234
00235 Info<< "Rotating points by " << T << endl;
00236
00237 points = transform(T, points);
00238
00239 if (args.optionFound("rotateFields"))
00240 {
00241 rotateFields(args, runTime, T);
00242 }
00243 }
00244 else if (args.optionFound("rollPitchYaw"))
00245 {
00246 vector v(args.optionLookup("rollPitchYaw")());
00247
00248 Info<< "Rotating points by" << nl
00249 << " roll " << v.x() << nl
00250 << " pitch " << v.y() << nl
00251 << " yaw " << v.z() << endl;
00252
00253
00254
00255 v *= pi/180.0;
00256
00257 quaternion R(v.x(), v.y(), v.z());
00258
00259 Info<< "Rotating points by quaternion " << R << endl;
00260 points = transform(R, points);
00261
00262 if (args.optionFound("rotateFields"))
00263 {
00264 rotateFields(args, runTime, R.R());
00265 }
00266 }
00267 else if (args.optionFound("yawPitchRoll"))
00268 {
00269 vector v(args.optionLookup("yawPitchRoll")());
00270
00271 Info<< "Rotating points by" << nl
00272 << " yaw " << v.x() << nl
00273 << " pitch " << v.y() << nl
00274 << " roll " << v.z() << endl;
00275
00276
00277
00278 v *= pi/180.0;
00279
00280 scalar yaw = v.x();
00281 scalar pitch = v.y();
00282 scalar roll = v.z();
00283
00284 quaternion R = quaternion(vector(0, 0, 1), yaw);
00285 R *= quaternion(vector(0, 1, 0), pitch);
00286 R *= quaternion(vector(1, 0, 0), roll);
00287
00288 Info<< "Rotating points by quaternion " << R << endl;
00289 points = transform(R, points);
00290
00291 if (args.optionFound("rotateFields"))
00292 {
00293 rotateFields(args, runTime, R.R());
00294 }
00295 }
00296
00297 if (args.optionFound("scale"))
00298 {
00299 vector scaleVector(args.optionLookup("scale")());
00300
00301 Info<< "Scaling points by " << scaleVector << endl;
00302
00303 points.replace(vector::X, scaleVector.x()*points.component(vector::X));
00304 points.replace(vector::Y, scaleVector.y()*points.component(vector::Y));
00305 points.replace(vector::Z, scaleVector.z()*points.component(vector::Z));
00306 }
00307
00308
00309 IOstream::defaultPrecision(10);
00310
00311 Info << "Writing points into directory " << points.path() << nl << endl;
00312 points.write();
00313
00314 return 0;
00315 }
00316
00317
00318