transformField.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) 2018 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
27\*---------------------------------------------------------------------------*/
28
29#include "transformField.H"
30#include "FieldM.H"
31#include "diagTensor.H"
32
33// * * * * * * * * * * * * * * * global functions * * * * * * * * * * * * * //
34
36(
37 vectorField& rtf,
38 const quaternion& q,
39 const vectorField& tf
40)
41{
42 tensor rot = q.R();
44}
45
46
48(
49 const quaternion& q,
50 const vectorField& fld
51)
52{
53 auto tresult = tmp<vectorField>::New(fld.size());
54 transform(tresult.ref(), q, fld);
55 return tresult;
56}
57
58
60(
61 const quaternion& q,
62 const tmp<vectorField>& tfld
63)
64{
65 tmp<vectorField> tresult = New(tfld);
66 transform(tresult.ref(), q, tfld());
67 tfld.clear();
68 return tresult;
69}
70
71
73(
74 vectorField& result,
75 const septernion& tr,
76 const vectorField& fld
77)
78{
79 vector trans = tr.t();
80
81 // Check if any translation
82 if (mag(trans) > VSMALL)
83 {
85 }
86 else
87 {
88 result = fld;
89 }
90
91 // Check if any rotation
92 if (mag(tr.r().R() - I) > SMALL)
93 {
94 transform(result, tr.r(), result);
95 }
96}
97
98
100(
101 const septernion& tr,
102 const vectorField& fld
103)
104{
105 auto tresult = tmp<vectorField>::New(fld.size());
106 transformPoints(tresult.ref(), tr, fld);
107 return tresult;
108}
109
110
112(
113 const septernion& tr,
114 const tmp<vectorField>& tfld
115)
116{
117 tmp<vectorField> tresult = New(tfld);
118 transformPoints(tresult.ref(), tr, tfld());
119 tfld.clear();
120 return tresult;
121}
122
123
124// ************************************************************************* //
High performance macro functions for Field<Type> algebra. These expand using either array element acc...
#define TFOR_ALL_F_OP_F_OP_S(typeF1, f1, OP1, typeF2, f2, OP2, typeS, s)
Definition: FieldM.H:306
#define TFOR_ALL_F_OP_FUNC_S_F(typeF1, f1, OP, FUNC, typeS, s, typeF2, f2)
Definition: FieldM.H:219
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))
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:58
tensor R() const
The rotation tensor corresponding to the quaternion.
Definition: quaternionI.H:358
Septernion class used to perform translations and rotations in 3D space.
Definition: septernion.H:67
A class for managing temporary objects.
Definition: tmp.H:65
void clear() const noexcept
Definition: tmpI.H:287
T & ref() const
Definition: tmpI.H:227
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionSet trans(const dimensionSet &ds)
Check the argument is dimensionless (for transcendental functions)
Definition: dimensionSet.C:486
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
void transformPoints(vectorField &, const septernion &, const vectorField &)
Transform given vectorField of coordinates with the given septernion.
Spatial transformation functions for primitive fields.