SphericalTensorI.H
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) 2020-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "Vector.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class Cmpt>
35 :
37 {}
38 
39 
40 template<class Cmpt>
41 template<class Cmpt2>
43 (
44  const VectorSpace<SphericalTensor<Cmpt2>, Cmpt2, 1>& vs
45 )
46 :
48 {}
49 
50 
51 template<class Cmpt>
53 {
54  this->v_[II] = stii;
55 }
56 
57 
58 template<class Cmpt>
60 :
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
67 template<class Cmpt>
68 inline Foam::scalar Foam::SphericalTensor<Cmpt>::diagSqr() const
69 {
70  return 3*Foam::magSqr(this->ii());
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75 
76 namespace Foam
77 {
78 
79 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
80 
81 //- Return the trace of a SphericalTensor
82 template<class Cmpt>
83 inline Cmpt tr(const SphericalTensor<Cmpt>& st)
84 {
85  return 3*st.ii();
86 }
87 
88 
89 //- Return the spherical part of a SphericalTensor, i.e. itself
90 template<class Cmpt>
92 {
93  return st;
94 }
95 
96 
97 //- Return the determinant of a SphericalTensor
98 template<class Cmpt>
99 inline Cmpt det(const SphericalTensor<Cmpt>& st)
100 {
101  return st.ii()*st.ii()*st.ii();
102 }
103 
104 
105 //- Return the inverse of a SphericalTensor
106 template<class Cmpt>
107 inline SphericalTensor<Cmpt> inv(const SphericalTensor<Cmpt>& st)
108 {
109  #ifdef FULLDEBUG
110  if (mag(st.ii()) < VSMALL)
111  {
113  << "SphericalTensor not invertible, determinant:"
114  << det(st) << " tensor:" << st << nl
115  << abort(FatalError);
116  }
117  #endif
118 
119  return SphericalTensor<Cmpt>(1/st.ii());
120 }
121 
122 
123 //- Return the square of Frobenius norm of a SphericalTensor
124 template<class Cmpt>
125 inline Foam::scalar magSqr(const SphericalTensor<Cmpt>& st)
126 {
127  return st.diagSqr();
128 }
129 
130 
131 //- Return the max component of a SphericalTensor
132 template<class Cmpt>
133 inline Cmpt cmptMax(const SphericalTensor<Cmpt>& st)
134 {
135  return st.ii();
136 }
137 
139 //- Return the min component of a SphericalTensor
140 template<class Cmpt>
141 inline Cmpt cmptMin(const SphericalTensor<Cmpt>& st)
142 {
143  return st.ii();
144 }
145 
146 
147 //- Return the sum of components of a SphericalTensor
148 template<class Cmpt>
149 inline Cmpt cmptSum(const SphericalTensor<Cmpt>& st)
150 {
151  return 3*st.ii();
152 }
153 
154 
155 //- Return the arithmetic average of components of a SphericalTensor
156 template<class Cmpt>
157 inline Cmpt cmptAv(const SphericalTensor<Cmpt>& st)
158 {
159  return st.ii();
160 }
161 
162 
163 //- Linear interpolation of spherical tensors a and b by factor t
164 template<class Cmpt>
165 inline SphericalTensor<Cmpt> lerp
166 (
167  const SphericalTensor<Cmpt>& a,
169  const scalar t
170 )
171 {
172  return SphericalTensor<Cmpt>((1-t)*a.ii() + t*b.ii());
173 }
174 
175 
176 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
177 
178 //- Division of a Cmpt by a SphericalTensor
179 template<class Cmpt>
181 operator/(const Cmpt s, const SphericalTensor<Cmpt>& st)
182 {
183  #ifdef FULLDEBUG
184  if (mag(st.ii()) < VSMALL)
185  {
187  << "Cmpt = " << s
188  << " is not divisible due to a zero element in SphericalTensor:"
189  << "SphericalTensor = " << st
190  << abort(FatalError);
191  }
192  #endif
193 
194  return SphericalTensor<Cmpt>(s/st.ii());
195 }
197 
198 //- Division of a SphericalTensor by a Cmpt
199 template<class Cmpt>
201 operator/(const SphericalTensor<Cmpt>& st, const Cmpt s)
202 {
203  #ifdef FULLDEBUG
204  if (mag(s) < VSMALL)
205  {
207  << "SphericalTensor = " << st
208  << " is not divisible due to a zero value in Cmpt:"
209  << "Cmpt = " << s
210  << abort(FatalError);
211  }
212  #endif
213 
214  return SphericalTensor<Cmpt>(st.ii()/s);
215 }
216 
217 
218 //- Inner-product of a SphericalTensor and a SphericalTensor
219 template<class Cmpt>
222 {
223  return SphericalTensor<Cmpt>(st1.ii()*st2.ii());
224 }
225 
226 
227 //- Inner-product of a SphericalTensor and a Vector
228 template<class Cmpt>
229 inline Vector<Cmpt>
230 operator&(const SphericalTensor<Cmpt>& st, const Vector<Cmpt>& v)
231 {
232  return Vector<Cmpt>
233  (
234  st.ii()*v.x(),
235  st.ii()*v.y(),
236  st.ii()*v.z()
237  );
238 }
239 
241 //- Inner-product of a Vector and a SphericalTensor
242 template<class Cmpt>
243 inline Vector<Cmpt>
244 operator&(const Vector<Cmpt>& v, const SphericalTensor<Cmpt>& st)
245 {
246  return Vector<Cmpt>
247  (
248  v.x()*st.ii(),
249  v.y()*st.ii(),
250  v.z()*st.ii()
251  );
252 }
253 
254 
255 //- Double-inner-product of a SphericalTensor and a SphericalTensor
256 template<class Cmpt>
257 inline Cmpt
258 operator&&(const SphericalTensor<Cmpt>& st1, const SphericalTensor<Cmpt>& st2)
259 {
260  return 3*st1.ii()*st2.ii();
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 
266 template<class Cmpt>
267 class outerProduct<SphericalTensor<Cmpt>, Cmpt>
268 {
269 public:
270 
271  typedef SphericalTensor<Cmpt> type;
272 };
273 
274 template<class Cmpt>
275 class outerProduct<Cmpt, SphericalTensor<Cmpt>>
276 {
277 public:
278 
279  typedef SphericalTensor<Cmpt> type;
280 };
281 
282 
283 template<class Cmpt>
284 class innerProduct<SphericalTensor<Cmpt>, SphericalTensor<Cmpt>>
285 {
286 public:
287 
288  typedef SphericalTensor<Cmpt> type;
289 };
290 
291 
292 template<class Cmpt>
293 class innerProduct<SphericalTensor<Cmpt>, Vector<Cmpt>>
294 {
295 public:
297  typedef Vector<Cmpt> type;
298 };
299 
300 template<class Cmpt>
301 class innerProduct<Vector<Cmpt>, SphericalTensor<Cmpt>>
302 {
303 public:
305  typedef Vector<Cmpt> type;
306 };
307 
308 
309 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310 
311 } // End namespace Foam
312 
313 // ************************************************************************* //
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &f1)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
Cmpt cmptSum(const SphericalTensor< Cmpt > &st)
Return the sum of components of a SphericalTensor.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Templated vector space.
Definition: VectorSpace.H:52
dimensionSet operator &&(const dimensionSet &ds1, const dimensionSet &ds2)
dimensionedScalar operator/(const scalar s1, const dimensionedScalar &ds2)
const Cmpt & y() const noexcept
Access to the vector y component.
Definition: Vector.H:140
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:118
dimensionedScalar det(const dimensionedSphericalTensor &dt)
SphericalTensor< Cmpt > sph(const DiagTensor< Cmpt > &dt)
Return the spherical part of a DiagTensor as a SphericalTensor.
Definition: DiagTensorI.H:110
scalar diagSqr() const
The L2-norm squared of the diagonal.
A templated (3 x 3) diagonal tensor of objects of <T>, effectively containing 1 element, derived from VectorSpace.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
const Cmpt & ii() const noexcept
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross-product operators.
Definition: Vector.H:58
const Cmpt & x() const noexcept
Access to the vector x component.
Definition: Vector.H:135
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
const Cmpt & z() const noexcept
Access to the vector z component.
Definition: Vector.H:145
tmp< GeometricField< Type, faPatchField, areaMesh > > operator &(const faMatrix< Type > &, const DimensionedField< Type, areaMesh > &)
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) - 2 >::type type
Definition: products.H:155
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
SphericalTensor()=default
Default construct.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127