Loading...
Searching...
No Matches
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\*---------------------------------------------------------------------------*/
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)
42 const tensor rot = q.R();
43 // std::transform
45}
46
47
49(
50 const quaternion& q,
51 const vectorField& fld
52)
54 auto tresult = tmp<vectorField>::New(fld.size());
55 transform(tresult.ref(), q, fld);
56 return tresult;
57}
58
59
61(
62 const quaternion& q,
63 const tmp<vectorField>& tfld
64)
65{
66 tmp<vectorField> tresult = New(tfld);
67 transform(tresult.ref(), q, tfld());
68 tfld.clear();
69 return tresult;
70}
71
72
74(
75 vectorField& result,
76 const septernion& tr,
77 const vectorField& fld
78)
79{
80 const vector trans = tr.t();
81
82 // Check if any translation
83 if (mag(trans) > VSMALL)
84 {
85 // std::transform
86 TFOR_ALL_F_OP_F_OP_S(vector, result, =, vector, fld, -, vector, trans);
87 }
88 else
89 {
90 result = fld;
91 }
92
93 // Check if any rotation
94 if (mag(tr.r().R() - I) > SMALL)
95 {
96 transform(result, tr.r(), result);
97 }
98}
99
100
102(
103 const septernion& tr,
104 const vectorField& fld
105)
107 auto tresult = tmp<vectorField>::New(fld.size());
108 transformPoints(tresult.ref(), tr, fld);
109 return tresult;
110}
111
112
114(
115 const septernion& tr,
116 const tmp<vectorField>& tfld
117)
118{
119 tmp<vectorField> tresult = New(tfld);
120 transformPoints(tresult.ref(), tr, tfld());
121 tfld.clear();
122 return tresult;
123}
124
125
126// ************************************************************************* //
Declaration macros for Field<Type> algebra.
#define TFOR_ALL_F_OP_F_OP_S(typeF1, f1, OP1, typeF2, f2, OP2, typeS, s)
Definition FieldM.H:436
#define TFOR_ALL_F_OP_FUNC_S_F(typeF1, f1, OP, FUNC, typeS, s, typeF2, f2)
Definition FieldM.H:277
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-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:54
tensor R() const
The rotation tensor corresponding to the quaternion.
Septernion class used to perform translations and rotations in 3D space.
Definition septernion.H:63
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
static const Identity< scalar > I
Definition Identity.H:100
Tensor< scalar > tensor
Definition symmTensor.H:57
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionSet trans(const dimensionSet &ds)
Check the argument is dimensionless (for transcendental functions).
Vector< scalar > vector
Definition vector.H:57
void transformPoints(vectorField &, const septernion &, const vectorField &)
Transform given vectorField of coordinates with the given septernion.
Spatial transformation functions for primitive fields.