Loading...
Searching...
No Matches
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-2025 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\*---------------------------------------------------------------------------*/
29#include "SphericalTensor.H"
30#include "SymmTensor.H"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
34template<class Cmpt>
37 VectorSpace<DiagTensor<Cmpt>, Cmpt, 3>(Foam::zero{})
38{}
39
40
41template<class Cmpt>
42template<class Cmpt2>
44(
45 const VectorSpace<DiagTensor<Cmpt2>, Cmpt2, 3>& vs
47:
48 VectorSpace<DiagTensor<Cmpt>, Cmpt, 3>(vs)
49{}
50
51
52template<class Cmpt>
54{
55 this->v_[XX] = st.ii();
56 this->v_[YY] = st.ii();
57 this->v_[ZZ] = st.ii();
58}
59
60
61template<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
75template<class Cmpt>
77:
78 VectorSpace<DiagTensor<Cmpt>, Cmpt, 3>(is)
79{}
80
81
82// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
83
84template<class Cmpt>
86{
87 return Vector<Cmpt>(this->xx(), this->yy(), this->zz());
88}
89
90
91template<class Cmpt>
92inline Foam::scalar Foam::DiagTensor<Cmpt>::diagSqr() const
93{
94 return
95 (
96 Foam::magSqr(this->xx())
97 + Foam::magSqr(this->yy())
98 + Foam::magSqr(this->zz())
99 );
100}
101
102
103// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104
105namespace Foam
106{
108// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
109
110//- Return the trace of a DiagTensor
111template<class Cmpt>
112inline Cmpt tr(const DiagTensor<Cmpt>& dt)
113{
114 return dt.xx() + dt.yy() + dt.zz();
115}
116
118//- Return the spherical part of a DiagTensor as a SphericalTensor
119template<class Cmpt>
121{
123 (
124 (1.0/3.0)*tr(dt)
125 );
126}
127
128
129//- Return the determinant of a DiagTensor
130template<class Cmpt>
131inline Cmpt det(const DiagTensor<Cmpt>& dt)
132{
133 return dt.xx()*dt.yy()*dt.zz();
134}
135
136
137//- Return the inverse of a DiagTensor as a DiagTensor
138template<class Cmpt>
139inline DiagTensor<Cmpt> inv(const DiagTensor<Cmpt>& dt)
141 #ifdef FULLDEBUG
142 if (mag(det(dt)) < VSMALL)
143 {
145 << "DiagTensor is not invertible due to the zero determinant:"
146 << "det(DiagTensor) = " << det(dt)
147 << abort(FatalError);
148 }
149 #endif
150
151 return DiagTensor<Cmpt>(1/dt.xx(), 1/dt.yy(), 1/dt.zz());
152}
153
154
155//- Return the diagonal of a Tensor as a DiagTensor
156template<class Cmpt>
157inline DiagTensor<Cmpt> diag(const Tensor<Cmpt>& t)
158{
159 return DiagTensor<Cmpt>(t.xx(), t.yy(), t.zz());
161
162
163//- Return the diagonal of a SymmTensor as a DiagTensor
164template<class Cmpt>
165inline DiagTensor<Cmpt> diag(const SymmTensor<Cmpt>& st)
166{
167 return DiagTensor<Cmpt>(st.xx(), st.yy(), st.zz());
168}
169
171template<class Cmpt>
172inline Foam::scalar magSqr(const DiagTensor<Cmpt>& t)
173{
174 return t.diagSqr();
175}
176
178//- Linear interpolation of diagonal tensors a and b by factor t
179template<class Cmpt>
181(
182 const DiagTensor<Cmpt>& a,
183 const DiagTensor<Cmpt>& b,
184 const scalar t
185)
186{
187 const scalar onet = (1-t);
188
189 return DiagTensor<Cmpt>
190 (
191 onet*a.xx() + t*b.xx(),
192 onet*a.yy() + t*b.yy(),
193 onet*a.zz() + t*b.zz()
194 );
195}
196
197
198// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
199
200//- Sum of a DiagTensor and a Tensor
201template<class Cmpt>
202inline Tensor<Cmpt>
203operator+(const DiagTensor<Cmpt>& dt1, const Tensor<Cmpt>& t2)
204{
205 return Tensor<Cmpt>
206 (
207 dt1.xx() + t2.xx(), t2.xy(), t2.xz(),
208 t2.yx(), dt1.yy() + t2.yy(), t2.yz(),
209 t2.zx(), t2.zy(), dt1.zz() + t2.zz()
210 );
211}
213
214//- Sum of a Tensor and a DiagTensor
215template<class Cmpt>
216inline Tensor<Cmpt>
217operator+(const Tensor<Cmpt>& t1, const DiagTensor<Cmpt>& dt2)
218{
219 return Tensor<Cmpt>
220 (
221 t1.xx() + dt2.xx(), t1.xy(), t1.xz(),
222 t1.yx(), t1.yy() + dt2.yy(), t1.yz(),
223 t1.zx(), t1.zy(), t1.zz() + dt2.zz()
224 );
225}
226
227
228//- Subtract a Tensor from a DiagTensor
229template<class Cmpt>
230inline Tensor<Cmpt>
231operator-(const DiagTensor<Cmpt>& dt1, const Tensor<Cmpt>& t2)
232{
233 return Tensor<Cmpt>
234 (
235 dt1.xx() - t2.xx(), -t2.xy(), -t2.xz(),
236 -t2.yx(), dt1.yy() - t2.yy(), -t2.yz(),
237 -t2.zx(), -t2.zy(), dt1.zz() - t2.zz()
238 );
239}
240
241
242//- Subtract a DiagTensor from a Tensor
243template<class Cmpt>
245operator-(const Tensor<Cmpt>& t1, const DiagTensor<Cmpt>& dt2)
246{
247 return Tensor<Cmpt>
248 (
249 t1.xx() - dt2.xx(), t1.xy(), t1.xz(),
250 t1.yx(), t1.yy() - dt2.yy(), t1.yz(),
251 t1.zx(), t1.zy(), t1.zz() - dt2.zz()
252 );
253}
254
255
256//- Division of a Cmpt by a DiagTensor
257template<class Cmpt>
258inline DiagTensor<Cmpt>
259operator/(const Cmpt s, const DiagTensor<Cmpt>& dt)
261 #ifdef FULLDEBUG
262 if (mag(det(dt)) < VSMALL)
263 {
265 << "Cmpt = " << s
266 << " is not divisible by the DiagTensor due to a zero element:"
267 << "DiagTensor = " << dt
268 << abort(FatalError);
269 }
270 #endif
271
272 return DiagTensor<Cmpt>(s/dt.xx(), s/dt.yy(), s/dt.zz());
273}
274
275
276//- Division of a DiagTensor by a Cmpt
277template<class Cmpt>
278inline DiagTensor<Cmpt>
279operator/(const DiagTensor<Cmpt>& dt, const Cmpt s)
280{
281 #ifdef FULLDEBUG
282 if (mag(s) < VSMALL)
283 {
285 << "DiagTensor = " << dt
286 << " is not divisible due to a zero value in Cmpt:"
287 << "Cmpt = " << s
288 << abort(FatalError);
289 }
290 #endif
291
292 return DiagTensor<Cmpt>(dt.xx()/s, dt.yy()/s, dt.zz()/s);
293}
294
295
296//- Division of a Vector by a DiagTensor
297template<class Cmpt>
299operator/(const Vector<Cmpt> v, const DiagTensor<Cmpt>& dt)
300{
301 #ifdef FULLDEBUG
302 if (mag(det(dt)) < VSMALL)
303 {
305 << "Vector = " << v
306 << " is not divisible by the DiagTensor due to a zero element:"
307 << "DiagTensor = " << dt
308 << abort(FatalError);
309 }
310 #endif
311
312 return Vector<Cmpt>(v.x()/dt.xx(), v.y()/dt.yy(), v.z()/dt.zz());
313}
314
315
316//- Inner-product of a DiagTensor and a DiagTensor
317template<class Cmpt>
318inline DiagTensor<Cmpt>
319operator&(const DiagTensor<Cmpt>& dt1, const DiagTensor<Cmpt>& dt2)
321 return DiagTensor<Cmpt>
322 (
323 dt1.xx()*dt2.xx(),
324 dt1.yy()*dt2.yy(),
325 dt1.zz()*dt2.zz()
326 );
327}
328
329
330//- Inner-product of a DiagTensor and a Tensor
331template<class Cmpt>
332inline Tensor<Cmpt>
333operator&(const DiagTensor<Cmpt>& dt1, const Tensor<Cmpt>& t2)
334{
335 return Tensor<Cmpt>
336 (
337 dt1.xx()*t2.xx(),
338 dt1.xx()*t2.xy(),
339 dt1.xx()*t2.xz(),
340
341 dt1.yy()*t2.yx(),
342 dt1.yy()*t2.yy(),
343 dt1.yy()*t2.yz(),
344
345 dt1.zz()*t2.zx(),
346 dt1.zz()*t2.zy(),
347 dt1.zz()*t2.zz()
348 );
349}
350
351
352//- Inner-product of a Tensor and a DiagTensor
353template<class Cmpt>
354inline Tensor<Cmpt>
355operator&(const Tensor<Cmpt>& t1, const DiagTensor<Cmpt>& dt2)
356{
357 return Tensor<Cmpt>
359 t1.xx()*dt2.xx(),
360 t1.xy()*dt2.yy(),
361 t1.xz()*dt2.zz(),
362
363 t1.yx()*dt2.xx(),
364 t1.yy()*dt2.yy(),
365 t1.yz()*dt2.zz(),
366
367 t1.zx()*dt2.xx(),
368 t1.zy()*dt2.yy(),
369 t1.zz()*dt2.zz()
370 );
371}
372
373
374//- Inner-product of a DiagTensor and a Vector
375template<class Cmpt>
376inline Vector<Cmpt>
377operator&(const DiagTensor<Cmpt>& dt, const Vector<Cmpt>& v)
378{
379 return Vector<Cmpt>
380 (
381 dt.xx()*v.x(),
382 dt.yy()*v.y(),
383 dt.zz()*v.z()
384 );
385}
386
387
388//- Inner-product of a Vector and a DiagTensor
389template<class Cmpt>
390inline Vector<Cmpt>
391operator&(const Vector<Cmpt>& v, const DiagTensor<Cmpt>& dt)
392{
393 return Vector<Cmpt>
394 (
395 v.x()*dt.xx(),
396 v.y()*dt.yy(),
397 v.z()*dt.zz()
398 );
399}
400
401
402// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
403
404} // End namespace Foam
405
406// ************************************************************************* //
A templated (3 x 3) diagonal tensor of objects of <T>, effectively containing 3 elements,...
Definition DiagTensor.H:56
scalar diagSqr() const
The L2-norm squared of the diagonal.
Definition DiagTensorI.H:85
const Cmpt & yy() const noexcept
Definition DiagTensor.H:133
DiagTensor()=default
Default construct.
const Cmpt & xx() const noexcept
Definition DiagTensor.H:132
Vector< Cmpt > diag() const
Extract the diagonal as a vector.
Definition DiagTensorI.H:78
const Cmpt & zz() const noexcept
Definition DiagTensor.H:134
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
A templated (3 x 3) diagonal tensor of objects of <T>, effectively containing 1 element,...
const Cmpt & ii() const noexcept
A templated (3 x 3) symmetric tensor of objects of <T>, effectively containing 6 elements,...
Definition SymmTensor.H:53
const Cmpt & yy() const noexcept
Definition SymmTensor.H:154
const Cmpt & xx() const noexcept
Definition SymmTensor.H:150
const Cmpt & zz() const noexcept
Definition SymmTensor.H:158
A templated (3 x 3) tensor of objects of <T> derived from MatrixSpace.
Definition Tensor.H:60
const Cmpt & yy() const noexcept
Definition Tensor.H:197
const Cmpt & zy() const noexcept
Definition Tensor.H:200
const Cmpt & yx() const noexcept
Definition Tensor.H:196
const Cmpt & xx() const noexcept
Definition Tensor.H:193
const Cmpt & yz() const noexcept
Definition Tensor.H:198
const Cmpt & zz() const noexcept
Definition Tensor.H:201
const Cmpt & xy() const noexcept
Definition Tensor.H:194
const Cmpt & xz() const noexcept
Definition Tensor.H:195
const Cmpt & zx() const noexcept
Definition Tensor.H:199
Templated vector space.
Definition VectorSpace.H:75
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
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.
tmp< faMatrix< Type > > operator-(const faMatrix< Type > &)
Unary negation.
dimensionedScalar det(const dimensionedSphericalTensor &dt)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
tmp< faMatrix< Type > > operator+(const faMatrix< Type > &, const faMatrix< Type > &)
dimensionedScalar operator/(const scalar s1, const dimensionedScalar &ds2)
tmp< GeometricField< Type, faPatchField, areaMesh > > operator&(const faMatrix< Type > &, const DimensionedField< Type, areaMesh > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
SphericalTensor< Cmpt > sph(const DiagTensor< Cmpt > &dt)
Return the spherical part of a DiagTensor as a SphericalTensor.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
volScalarField & b