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