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 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 "Vector.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class Cmpt>
35:
37{}
38
39
40template<class Cmpt>
41template<class Cmpt2>
43(
44 const VectorSpace<SphericalTensor<Cmpt2>, Cmpt2, 1>& vs
45)
46:
48{}
49
50
51template<class Cmpt>
53{
54 this->v_[II] = stii;
55}
56
57
58template<class Cmpt>
60:
62{}
63
64
65// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66
67template<class Cmpt>
68inline const Cmpt& Foam::SphericalTensor<Cmpt>::ii() const
69{
70 return this->v_[II];
71}
72
73
74template<class Cmpt>
76{
77 return this->v_[II];
78}
79
80
81template<class Cmpt>
84{
85 return *this;
86}
87
88
89// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
90
91namespace Foam
92{
93
94// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
95
96//- Return the trace of a SphericalTensor
97template<class Cmpt>
98inline Cmpt tr(const SphericalTensor<Cmpt>& st)
99{
100 return 3*st.ii();
102
103
104//- Return the spherical part of a SphericalTensor, i.e. itself
105template<class Cmpt>
107{
108 return st;
109}
110
112//- Return the determinant of a SphericalTensor
113template<class Cmpt>
114inline Cmpt det(const SphericalTensor<Cmpt>& st)
115{
116 return st.ii()*st.ii()*st.ii();
117}
118
119
120//- Return the inverse of a SphericalTensor
121template<class Cmpt>
123{
124 #ifdef FULLDEBUG
125 if (mag(st.ii()) < VSMALL)
126 {
128 << "SphericalTensor is not invertible due to the zero determinant:"
129 << "det(SphericalTensor) = " << det(st)
130 << abort(FatalError);
131 }
132 #endif
133
134 return SphericalTensor<Cmpt>(1/st.ii());
135}
136
137
138//- Return the square of Frobenius norm of a SphericalTensor as a Cmpt
139template<class Cmpt>
140inline Cmpt magSqr(const SphericalTensor<Cmpt>& st)
141{
142 return Cmpt(3*mag(st.ii()*st.ii()));
143}
144
145
146//- Return the max component of a SphericalTensor
147template<class Cmpt>
148inline Cmpt cmptMax(const SphericalTensor<Cmpt>& st)
149{
150 return st.ii();
151}
152
153
154//- Return the min component of a SphericalTensor
155template<class Cmpt>
156inline Cmpt cmptMin(const SphericalTensor<Cmpt>& st)
157{
158 return st.ii();
159}
160
161
162//- Return the sum of components of a SphericalTensor
163template<class Cmpt>
164inline Cmpt cmptSum(const SphericalTensor<Cmpt>& st)
165{
166 return 3*st.ii();
167}
168
169
170//- Return the arithmetic average of components of a SphericalTensor
171template<class Cmpt>
172inline Cmpt cmptAv(const SphericalTensor<Cmpt>& st)
173{
174 return st.ii();
175}
176
177
178// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
179
180//- Division of a Cmpt by a SphericalTensor
181template<class Cmpt>
182inline SphericalTensor<Cmpt>
183operator/(const Cmpt s, const SphericalTensor<Cmpt>& st)
184{
185 #ifdef FULLDEBUG
186 if (mag(st.ii()) < VSMALL)
187 {
189 << "Cmpt = " << s
190 << " is not divisible due to a zero element in SphericalTensor:"
191 << "SphericalTensor = " << st
192 << abort(FatalError);
193 }
194 #endif
195
196 return SphericalTensor<Cmpt>(s/st.ii());
197}
198
199
200//- Division of a SphericalTensor by a Cmpt
201template<class Cmpt>
202inline SphericalTensor<Cmpt>
203operator/(const SphericalTensor<Cmpt>& st, const Cmpt s)
204{
205 #ifdef FULLDEBUG
206 if (mag(s) < VSMALL)
207 {
209 << "SphericalTensor = " << st
210 << " is not divisible due to a zero value in Cmpt:"
211 << "Cmpt = " << s
212 << abort(FatalError);
213 }
214 #endif
215
216 return SphericalTensor<Cmpt>(st.ii()/s);
217}
218
219
220//- Inner-product of a SphericalTensor and a SphericalTensor
221template<class Cmpt>
222inline SphericalTensor<Cmpt>
224{
225 return SphericalTensor<Cmpt>(st1.ii()*st2.ii());
226}
227
228
229//- Inner-product of a SphericalTensor and a Vector
230template<class Cmpt>
231inline Vector<Cmpt>
233{
234 return Vector<Cmpt>
235 (
236 st.ii()*v.x(),
237 st.ii()*v.y(),
238 st.ii()*v.z()
239 );
240}
241
242
243//- Inner-product of a Vector and a SphericalTensor
244template<class Cmpt>
245inline Vector<Cmpt>
247{
248 return Vector<Cmpt>
249 (
250 v.x()*st.ii(),
251 v.y()*st.ii(),
252 v.z()*st.ii()
253 );
254}
255
256
257//- Double-inner-product of a SphericalTensor and a SphericalTensor
258template<class Cmpt>
259inline Cmpt
261{
262 return 3*st1.ii()*st2.ii();
263}
264
265
266// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
267
268template<class Cmpt>
270{
271public:
272
274};
275
276template<class Cmpt>
278{
279public:
280
282};
283
284
285template<class Cmpt>
287{
288public:
289
291};
292
293
294template<class Cmpt>
296{
297public:
298
300};
301
302template<class Cmpt>
304{
305public:
306
308};
309
310
311// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
312
313} // End namespace Foam
314
315// ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
A templated (3 x 3) diagonal tensor of objects of <T>, effectively containing 1 element,...
SphericalTensor()=default
Default construct.
const Cmpt & ii() const
const SphericalTensor< Cmpt > & T() const
Return non-Hermitian transpose (no-op)
Templated vector space.
Definition: VectorSpace.H:79
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:65
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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))
Namespace for OpenFOAM.
dimensionedScalar det(const dimensionedSphericalTensor &dt)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Cmpt cmptSum(const SphericalTensor< Cmpt > &st)
Return the sum of components of a SphericalTensor.
tmp< GeometricField< Type, fvPatchField, volMesh > > operator&(const fvMatrix< Type > &, const DimensionedField< Type, volMesh > &)
dimensionedScalar operator/(const scalar s1, const dimensionedScalar &ds2)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
SphericalTensor< Cmpt > sph(const DiagTensor< Cmpt > &dt)
Return the spherical part of a DiagTensor as a SphericalTensor.
Definition: DiagTensorI.H:130
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
error FatalError
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensioned< typename scalarProduct< Type1, Type2 >::type > operator&&(const dimensioned< Type1 > &, const dimensioned< Type2 > &)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh > > cmptAv(const DimensionedField< Type, GeoMesh > &df)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)