cylindrical.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 "cylindrical.H"
27 #include "axesRotation.H"
28 #include "polyMesh.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
38 }
39 
40 
41 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
42 
44 {
45  vector dir = p - origin_;
46  dir /= mag(dir) + vSmall;
47 
48  const vector axis = axis_/mag(axis_);
49  const vector r = dir - (dir & axis)*axis;
50 
51  if (mag(r) < small)
52  {
53  // If the point is on the axis choose any radial direction
54  return axesRotation(axis, perpendicular(axis)).R();
55  }
56  else
57  {
58  return axesRotation(axis, dir).R();
59  }
60 
61  return tensor(r, axis^r, axis);
62 }
63 
64 
65 void Foam::cylindrical::init(const UList<vector>& points)
66 {
67  Rptr_.reset(new tensorField(points.size()));
68  tensorField& R = Rptr_();
69 
70  forAll(points, i)
71  {
72  R[i] = this->R(points[i]);
73  }
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78 
80 (
81  const vector& axis,
82  const point& origin,
83  const UList<vector>& points
84 )
85 :
86  Rptr_(),
87  origin_(origin),
88  axis_(axis)
89 {
90  init(points);
91 }
92 
93 
95 :
96  Rptr_(),
97  origin_
98  (
99  dict.parent().found("origin")
100  ? dict.parent().lookup("origin")
101  : dict.lookup("origin")
102  ),
103  axis_(dict.lookupBackwardsCompatible<vector>({"axis", "e3"}))
104 {}
105 
106 
108 (
109  const dictionary& dict,
110  const UList<vector>& points
111 )
112 :
114 {
115  init(points);
116 }
117 
118 
119 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
120 
122 {
123  if (Rptr_.valid())
124  {
125  Rptr_().setSize(points.size());
126  }
127  else
128  {
129  Rptr_.reset(new tensorField(points.size()));
130  }
131 
132  tensorField& R = Rptr_();
133 
134  forAll(points, i)
135  {
136  R[i] = this->R(points[i]);
137  }
138 }
139 
140 
142 (
143  const vectorField& vf
144 ) const
145 {
146  if (Rptr_->size() != vf.size())
147  {
149  << "vectorField st has different size to tensorField "
150  << abort(FatalError);
151  }
152 
153  return (Rptr_() & vf);
154 }
155 
156 
158 {
160  return Zero;
161 }
162 
163 
165 (
166  const vector& v,
167  const label cmptI
168 ) const
169 {
170  return (Rptr_()[cmptI] & v);
171 }
172 
173 
175 (
176  const vectorField& vf
177 ) const
178 {
179  if (Rptr_->size() != vf.size())
180  {
182  << "vectorField st has different size to tensorField "
183  << abort(FatalError);
184  }
185 
186  return (Rptr_().T() & vf);
187 }
188 
189 
191 {
193  return Zero;
194 }
195 
196 
198 (
199  const vector& v,
200  const label cmptI
201 ) const
202 {
203  return (Rptr_()[cmptI].T() & v);
204 }
205 
206 
208 (
209  const tensorField& tf
210 ) const
211 {
212  if (Rptr_->size() != tf.size())
213  {
215  << "tensorField st has different size to tensorField Tr"
216  << abort(FatalError);
217  }
218  return (Rptr_() & tf & Rptr_().T());
219 }
220 
221 
223 (
224  const vector& p,
225  const tensor& t
226 ) const
227 {
228  const tensor R(this->R(p));
229  return (R & t & R.T());
230 }
231 
232 
234 (
235  const vectorField& vf
236 ) const
237 {
238  if (Rptr_->size() != vf.size())
239  {
241  << "tensorField vf has different size to tensorField Tr"
242  << abort(FatalError);
243  }
244 
245  tmp<symmTensorField> tfld(new symmTensorField(Rptr_->size()));
246  symmTensorField& fld = tfld.ref();
247 
248  const tensorField& R = Rptr_();
249  forAll(fld, i)
250  {
251  fld[i] = transformVectorDiagTensor(R[i], vf[i]);
252  }
253  return tfld;
254 }
255 
256 
258 (
259  const vector& p,
260  const vector& v
261 ) const
262 {
263  return transformVectorDiagTensor(R(p), v);
264 }
265 
266 
268 {
269  writeEntry(os, "axis", axis());
270 }
271 
272 
273 // ************************************************************************* //
bool found
#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
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
Abstract base class for coordinate rotation.
A local coordinate rotation.
Definition: cylindrical.H:68
virtual void write(Ostream &) const
Write.
Definition: cylindrical.C:267
virtual void updatePoints(const UList< vector > &points)
Update the rotation for a list of points.
Definition: cylindrical.C:121
virtual const tensor & R() const
Return local-to-global transformation tensor.
Definition: cylindrical.H:130
virtual vector transform(const vector &v) const
Transform vector using transformation tensor.
Definition: cylindrical.C:157
virtual const vector axis() const
Return local Cartesian z-axis in global coordinates.
Definition: cylindrical.H:157
virtual symmTensor transformDiagTensor(const vector &p, const vector &v) const
Transform diagTensor masquerading as a vector using transformation.
Definition: cylindrical.C:258
cylindrical(const vector &axis, const point &origin, const UList< vector > &points)
Construct from components for list of points.
Definition: cylindrical.C:80
virtual vector invTransform(const vector &v) const
Inverse transform vector using transformation tensor.
Definition: cylindrical.C:190
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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))
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
errorManip< error > abort(error &err)
Definition: errorManip.H:131
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
Vector< Cmpt > perpendicular(const Vector< Cmpt > &v)
Definition: VectorI.H:166
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.
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
error FatalError
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
dictionary dict
volScalarField & p