axesRotation.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 \*---------------------------------------------------------------------------*/
25 
26 #include "axesRotation.H"
27 #include "dictionary.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
37 }
38 
39 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 
41 void Foam::axesRotation::calcTransform
42 (
43  const vector& axis1,
44  const vector& axis2,
45  const axisOrder& order
46 )
47 {
48  vector a = axis1/mag(axis1);
49  vector b = axis2;
50 
51  b = b - (b & a)*a;
52 
53  if (mag(b) < small)
54  {
56  << "axis1, axis2 appear co-linear: "
57  << axis1 << ", " << axis2 << endl
58  << abort(FatalError);
59  }
60 
61  b = b/mag(b);
62  vector c = a^b;
63 
64  tensor Rtr;
65  switch (order)
66  {
67  case e1e2:
68  {
69  Rtr = tensor(a, b, c);
70  break;
71  }
72  case e2e3:
73  {
74  Rtr = tensor(c, a, b);
75  break;
76  }
77  case e3e1:
78  {
79  Rtr = tensor(b, c, a);
80  break;
81  }
82  default:
83  {
85  << "Unhandled axes specification" << endl
86  << abort(FatalError);
87 
88  Rtr = Zero;
89  break;
90  }
91  }
92 
93  // Global->local transformation
94  Rtr_ = Rtr;
95 
96  // Local->global transformation
97  R_ = Rtr.T();
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
102 
104 (
105  const vector& axis,
106  const vector& dir
107 )
108 :
109  R_(sphericalTensor::I),
110  Rtr_(R_)
111 {
112  calcTransform(axis, dir, e3e1);
113 }
114 
115 
117 :
118  R_(R),
119  Rtr_(R_.T())
120 {}
121 
122 
124 (
125  const dictionary& dict
126 )
127 :
128  R_(sphericalTensor::I),
129  Rtr_(R_)
130 {
131  operator=(dict);
132 }
133 
134 
136 (
137  const dictionary& dict,
138  const UList<vector>& points
139 )
140 :
142 {}
143 
144 
145 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
146 
148 (
149  const vectorField& vf
150 ) const
151 {
152  return (R_ & vf);
153 }
154 
155 
157 {
158  return (R_ & v);
159 }
160 
161 
163 (
164  const vectorField& vf
165 ) const
166 {
167  return (Rtr_ & vf);
168 }
169 
170 
172 {
173  return (Rtr_ & v);
174 }
175 
176 
178 (
179  const tensorField& tf
180 ) const
181 {
183  return tmp<tensorField>(nullptr);
184 }
185 
186 
188 (
189  const vector& p,
190  const tensor& t
191 ) const
192 {
193  return (R_ & t & Rtr_);
194 }
195 
196 
198 (
199  const vectorField& vf
200 ) const
201 {
202  tmp<symmTensorField> tfld(new symmTensorField(vf.size()));
203  symmTensorField& fld = tfld.ref();
204 
205  forAll(fld, i)
206  {
207  fld[i] = transformVectorDiagTensor(R_, vf[i]);
208  }
209  return tfld;
210 }
211 
212 
214 (
215  const vector& p,
216  const vector& v
217 ) const
218 {
219  return transformVectorDiagTensor(R_, v);
220 }
221 
222 
224 {
225  writeEntry(os, "e1", e1());
226  writeEntry(os, "e2", e2());
227  writeEntry(os, "e3", e3());
228 }
229 
230 
231 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
232 
234 {
235  if (debug)
236  {
237  Pout<< "axesRotation::operator=(const dictionary&) : "
238  << "assign from " << dict << endl;
239  }
240 
241  vector axis1, axis2;
242  axisOrder order(e3e1);
243 
244  if (dict.readIfPresent("e1", axis1) && dict.readIfPresent("e2", axis2))
245  {
246  order = e1e2;
247  }
248  else if (dict.readIfPresent("e2", axis1)&& dict.readIfPresent("e3", axis2))
249  {
250  order = e2e3;
251  }
252  else if (dict.readIfPresent("e3", axis1)&& dict.readIfPresent("e1", axis2))
253  {
254  order = e3e1;
255  }
256  else if (dict.found("axis") || dict.found("direction"))
257  {
258  // Both "axis" and "direction" are required
259  // If one is missing the appropriate error message will be generated
260  order = e3e1;
261  dict.lookup("axis") >> axis1;
262  dict.lookup("direction") >> axis2;
263  }
264  else
265  {
267  << "not entry of the type (e1, e2) or (e2, e3) or (e3, e1) "
268  << "found "
269  << exit(FatalError);
270  }
271 
272  calcTransform(axis1, axis2, order);
273 }
274 
275 
276 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Templated 3D SphericalTensor derived from VectorSpace adding construction from 1 component,...
Tensor< Cmpt > T() const
Return transpose.
Definition: TensorI.H:331
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:74
A coordinate rotation specified using global axis.
Definition: axesRotation.H:67
virtual void write(Ostream &) const
Write.
Definition: axesRotation.C:223
axesRotation(const vector &axis, const vector &dir)
Construct from 2 axes.
Definition: axesRotation.C:104
virtual vector transform(const vector &v) const
Transform vector using transformation tensor.
Definition: axesRotation.C:156
virtual symmTensor transformDiagTensor(const vector &p, const vector &v) const
Transform diagTensor masquerading as a vector using transformation.
Definition: axesRotation.C:214
void operator=(const dictionary &)
Assign from dictionary.
Definition: axesRotation.C:233
virtual vector invTransform(const vector &v) const
Inverse transform vector using transformation tensor.
Definition: axesRotation.C:171
Abstract base class for coordinate rotation.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const tensorField & tf
const pointField & points
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:27
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
errorManip< error > abort(error &err)
Definition: errorManip.H:131
static const Identity< scalar > I
Definition: Identity.H:93
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
dictionary dict
volScalarField & p