DiagTensorI.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 "SphericalTensor.H"
30 #include "SymmTensor.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class Cmpt>
36 :
37  VectorSpace<DiagTensor<Cmpt>, Cmpt, 3>(Zero)
38 {}
39 
40 
41 template<class Cmpt>
42 template<class Cmpt2>
44 (
45  const VectorSpace<DiagTensor<Cmpt2>, Cmpt2, 3>& vs
46 )
47 :
48  VectorSpace<DiagTensor<Cmpt>, Cmpt, 3>(vs)
49 {}
50 
51 
52 template<class Cmpt>
54 {
55  this->v_[XX] = st.ii();
56  this->v_[YY] = st.ii();
57  this->v_[ZZ] = st.ii();
58 }
59 
60 
61 template<class Cmpt>
63 (
64  const Cmpt& vxx,
65  const Cmpt& vyy,
66  const Cmpt& vzz
67 )
68 {
69  this->v_[XX] = vxx;
70  this->v_[YY] = vyy;
71  this->v_[ZZ] = vzz;
72 }
73 
74 
75 template<class Cmpt>
77 :
78  VectorSpace<DiagTensor<Cmpt>, Cmpt, 3>(is)
79 {}
80 
81 
82 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
83 
84 template<class Cmpt>
85 inline Foam::scalar Foam::DiagTensor<Cmpt>::diagSqr() const
86 {
87  return
88  (
89  Foam::magSqr(this->xx())
90  + Foam::magSqr(this->yy())
91  + Foam::magSqr(this->zz())
92  );
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97 
98 namespace Foam
99 {
101 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
102 
103 //- Return the trace of a DiagTensor
104 template<class Cmpt>
105 inline Cmpt tr(const DiagTensor<Cmpt>& dt)
106 {
107  return dt.xx() + dt.yy() + dt.zz();
108 }
109 
111 //- Return the spherical part of a DiagTensor as a SphericalTensor
112 template<class Cmpt>
114 {
115  return SphericalTensor<Cmpt>
116  (
117  (1.0/3.0)*tr(dt)
118  );
119 }
120 
121 
122 //- Return the determinant of a DiagTensor
123 template<class Cmpt>
124 inline Cmpt det(const DiagTensor<Cmpt>& dt)
125 {
126  return dt.xx()*dt.yy()*dt.zz();
127 }
128 
129 
130 //- Return the inverse of a DiagTensor as a DiagTensor
131 template<class Cmpt>
132 inline DiagTensor<Cmpt> inv(const DiagTensor<Cmpt>& dt)
133 {
134  #ifdef FULLDEBUG
135  if (mag(det(dt)) < VSMALL)
136  {
138  << "DiagTensor is not invertible due to the zero determinant:"
139  << "det(DiagTensor) = " << det(dt)
140  << abort(FatalError);
141  }
142  #endif
143 
144  return DiagTensor<Cmpt>(1/dt.xx(), 1/dt.yy(), 1/dt.zz());
145 }
146 
147 
148 //- Return the diagonal of a Tensor as a DiagTensor
149 template<class Cmpt>
150 inline DiagTensor<Cmpt> diag(const Tensor<Cmpt>& t)
151 {
152  return DiagTensor<Cmpt>(t.xx(), t.yy(), t.zz());
153 }
154 
155 
156 //- Return the diagonal of a SymmTensor as a DiagTensor
157 template<class Cmpt>
158 inline DiagTensor<Cmpt> diag(const SymmTensor<Cmpt>& st)
159 {
160  return DiagTensor<Cmpt>(st.xx(), st.yy(), st.zz());
161 }
162 
164 template<class Cmpt>
165 inline Foam::scalar magSqr(const DiagTensor<Cmpt>& t)
166 {
167  return t.diagSqr();
168 }
169 
171 //- Linear interpolation of diagonal tensors a and b by factor t
172 template<class Cmpt>
173 inline DiagTensor<Cmpt> lerp
174 (
175  const DiagTensor<Cmpt>& a,
176  const DiagTensor<Cmpt>& b,
177  const scalar t
178 )
179 {
180  const scalar onet = (1-t);
182  return DiagTensor<Cmpt>
183  (
184  onet*a.xx() + t*b.xx(),
185  onet*a.yy() + t*b.yy(),
186  onet*a.zz() + t*b.zz()
187  );
188 }
189 
190 
191 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
192 
193 //- Sum of a DiagTensor and a Tensor
194 template<class Cmpt>
195 inline Tensor<Cmpt>
196 operator+(const DiagTensor<Cmpt>& dt1, const Tensor<Cmpt>& t2)
197 {
198  return Tensor<Cmpt>
199  (
200  dt1.xx() + t2.xx(), t2.xy(), t2.xz(),
201  t2.yx(), dt1.yy() + t2.yy(), t2.yz(),
202  t2.zx(), t2.zy(), dt1.zz() + t2.zz()
203  );
204 }
206 
207 //- Sum of a Tensor and a DiagTensor
208 template<class Cmpt>
209 inline Tensor<Cmpt>
210 operator+(const Tensor<Cmpt>& t1, const DiagTensor<Cmpt>& dt2)
211 {
212  return Tensor<Cmpt>
213  (
214  t1.xx() + dt2.xx(), t1.xy(), t1.xz(),
215  t1.yx(), t1.yy() + dt2.yy(), t1.yz(),
216  t1.zx(), t1.zy(), t1.zz() + dt2.zz()
217  );
218 }
219 
220 
221 //- Subtract a Tensor from a DiagTensor
222 template<class Cmpt>
223 inline Tensor<Cmpt>
224 operator-(const DiagTensor<Cmpt>& dt1, const Tensor<Cmpt>& t2)
225 {
226  return Tensor<Cmpt>
227  (
228  dt1.xx() - t2.xx(), -t2.xy(), -t2.xz(),
229  -t2.yx(), dt1.yy() - t2.yy(), -t2.yz(),
230  -t2.zx(), -t2.zy(), dt1.zz() - t2.zz()
231  );
232 }
233 
234 
235 //- Subtract a DiagTensor from a Tensor
236 template<class Cmpt>
238 operator-(const Tensor<Cmpt>& t1, const DiagTensor<Cmpt>& dt2)
239 {
240  return Tensor<Cmpt>
241  (
242  t1.xx() - dt2.xx(), t1.xy(), t1.xz(),
243  t1.yx(), t1.yy() - dt2.yy(), t1.yz(),
244  t1.zx(), t1.zy(), t1.zz() - dt2.zz()
245  );
246 }
247 
248 
249 //- Division of a Cmpt by a DiagTensor
250 template<class Cmpt>
251 inline DiagTensor<Cmpt>
252 operator/(const Cmpt s, const DiagTensor<Cmpt>& dt)
253 {
254  #ifdef FULLDEBUG
255  if (mag(det(dt)) < VSMALL)
256  {
258  << "Cmpt = " << s
259  << " is not divisible by the DiagTensor due to a zero element:"
260  << "DiagTensor = " << dt
261  << abort(FatalError);
262  }
263  #endif
264 
265  return DiagTensor<Cmpt>(s/dt.xx(), s/dt.yy(), s/dt.zz());
266 }
267 
268 
269 //- Division of a DiagTensor by a Cmpt
270 template<class Cmpt>
271 inline DiagTensor<Cmpt>
272 operator/(const DiagTensor<Cmpt>& dt, const Cmpt s)
273 {
274  #ifdef FULLDEBUG
275  if (mag(s) < VSMALL)
276  {
278  << "DiagTensor = " << dt
279  << " is not divisible due to a zero value in Cmpt:"
280  << "Cmpt = " << s
281  << abort(FatalError);
282  }
283  #endif
284 
285  return DiagTensor<Cmpt>(dt.xx()/s, dt.yy()/s, dt.zz()/s);
286 }
287 
288 
289 //- Division of a Vector by a DiagTensor
290 template<class Cmpt>
292 operator/(const Vector<Cmpt> v, const DiagTensor<Cmpt>& dt)
293 {
294  #ifdef FULLDEBUG
295  if (mag(det(dt)) < VSMALL)
296  {
298  << "Vector = " << v
299  << " is not divisible by the DiagTensor due to a zero element:"
300  << "DiagTensor = " << dt
301  << abort(FatalError);
302  }
303  #endif
304 
305  return Vector<Cmpt>(v.x()/dt.xx(), v.y()/dt.yy(), v.z()/dt.zz());
306 }
307 
308 
309 //- Inner-product of a DiagTensor and a DiagTensor
310 template<class Cmpt>
311 inline DiagTensor<Cmpt>
312 operator&(const DiagTensor<Cmpt>& dt1, const DiagTensor<Cmpt>& dt2)
313 {
314  return DiagTensor<Cmpt>
315  (
316  dt1.xx()*dt2.xx(),
317  dt1.yy()*dt2.yy(),
318  dt1.zz()*dt2.zz()
319  );
320 }
321 
322 
323 //- Inner-product of a DiagTensor and a Tensor
324 template<class Cmpt>
325 inline Tensor<Cmpt>
326 operator&(const DiagTensor<Cmpt>& dt1, const Tensor<Cmpt>& t2)
327 {
328  return Tensor<Cmpt>
329  (
330  dt1.xx()*t2.xx(),
331  dt1.xx()*t2.xy(),
332  dt1.xx()*t2.xz(),
333 
334  dt1.yy()*t2.yx(),
335  dt1.yy()*t2.yy(),
336  dt1.yy()*t2.yz(),
337 
338  dt1.zz()*t2.zx(),
339  dt1.zz()*t2.zy(),
340  dt1.zz()*t2.zz()
341  );
342 }
343 
344 
345 //- Inner-product of a Tensor and a DiagTensor
346 template<class Cmpt>
347 inline Tensor<Cmpt>
348 operator&(const Tensor<Cmpt>& t1, const DiagTensor<Cmpt>& dt2)
349 {
350  return Tensor<Cmpt>
351  (
352  t1.xx()*dt2.xx(),
353  t1.xy()*dt2.yy(),
354  t1.xz()*dt2.zz(),
355 
356  t1.yx()*dt2.xx(),
357  t1.yy()*dt2.yy(),
358  t1.yz()*dt2.zz(),
359 
360  t1.zx()*dt2.xx(),
361  t1.zy()*dt2.yy(),
362  t1.zz()*dt2.zz()
363  );
364 }
365 
366 
367 //- Inner-product of a DiagTensor and a Vector
368 template<class Cmpt>
369 inline Vector<Cmpt>
370 operator&(const DiagTensor<Cmpt>& dt, const Vector<Cmpt>& v)
371 {
372  return Vector<Cmpt>
373  (
374  dt.xx()*v.x(),
375  dt.yy()*v.y(),
376  dt.zz()*v.z()
377  );
378 }
379 
380 
381 //- Inner-product of a Vector and a DiagTensor
382 template<class Cmpt>
383 inline Vector<Cmpt>
384 operator&(const Vector<Cmpt>& v, const DiagTensor<Cmpt>& dt)
385 {
386  return Vector<Cmpt>
387  (
388  v.x()*dt.xx(),
389  v.y()*dt.yy(),
390  v.z()*dt.zz()
391  );
392 }
393 
394 
395 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
396 
397 } // End namespace Foam
398 
399 // ************************************************************************* //
const Cmpt & zz() const noexcept
Definition: DiagTensor.H:134
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...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const Cmpt & yx() const noexcept
Definition: Tensor.H:196
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
const Cmpt & xy() const noexcept
Definition: Tensor.H:194
Templated vector space.
Definition: VectorSpace.H:52
dimensionedScalar operator/(const scalar s1, const dimensionedScalar &ds2)
const Cmpt & y() const noexcept
Access to the vector y component.
Definition: Vector.H:140
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
tmp< faMatrix< Type > > operator+(const faMatrix< Type > &, const faMatrix< Type > &)
A templated (3 x 3) diagonal tensor of objects of <T>, effectively containing 1 element, derived from VectorSpace.
const Cmpt & yz() const noexcept
Definition: Tensor.H:198
const Cmpt & xx() const noexcept
Definition: Tensor.H:193
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
const Cmpt & yy() const noexcept
Definition: Tensor.H:197
DiagTensor()=default
Default construct.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
const Cmpt & zx() const noexcept
Definition: Tensor.H:199
A templated (3 x 3) diagonal tensor of objects of <T>, effectively containing 3 elements, derived from VectorSpace.
Definition: DiagTensor.H:53
const Cmpt & ii() const noexcept
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:139
tmp< faMatrix< Type > > operator-(const faMatrix< Type > &)
Unary negation.
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 & yy() const noexcept
Definition: DiagTensor.H:133
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 & zy() const noexcept
Definition: Tensor.H:200
const Cmpt & z() const noexcept
Access to the vector z component.
Definition: Vector.H:145
const Cmpt & zz() const noexcept
Definition: Tensor.H:201
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
const Cmpt & xx() const noexcept
Definition: DiagTensor.H:132
scalar diagSqr() const
The L2-norm squared of the diagonal.
Definition: DiagTensorI.H:78
A templated (3 x 3) tensor of objects of <T> derived from MatrixSpace.
Definition: complexI.H:266
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))
const Cmpt & xz() const noexcept
Definition: Tensor.H:195
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127