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 -------------------------------------------------------------------------------
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  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 
66 template<class Cmpt>
68 :
69  VectorSpace<DiagTensor<Cmpt>, Cmpt, 3>(is)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74 
75 namespace Foam
76 {
77 
78 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
79 
80 //- Return the trace of a DiagTensor
81 template<class Cmpt>
82 inline Cmpt tr(const DiagTensor<Cmpt>& dt)
83 {
84  return dt.xx() + dt.yy() + dt.zz();
85 }
86 
87 
88 //- Return the spherical part of a DiagTensor as a SphericalTensor
89 template<class Cmpt>
91 {
93  (
94  (1.0/3.0)*tr(dt)
95  );
96 }
97 
98 
99 //- Return the determinant of a DiagTensor
100 template<class Cmpt>
101 inline Cmpt det(const DiagTensor<Cmpt>& dt)
102 {
103  return dt.xx()*dt.yy()*dt.zz();
104 }
105 
106 
107 //- Return the inverse of a DiagTensor as a DiagTensor
108 template<class Cmpt>
109 inline DiagTensor<Cmpt> inv(const DiagTensor<Cmpt>& dt)
110 {
111  #ifdef FULLDEBUG
112  if (mag(det(dt)) < VSMALL)
113  {
115  << "DiagTensor is not invertible due to the zero determinant:"
116  << "det(DiagTensor) = " << det(dt)
117  << abort(FatalError);
118  }
119  #endif
120 
121  return DiagTensor<Cmpt>(1/dt.xx(), 1/dt.yy(), 1/dt.zz());
122 }
123 
124 
125 //- Return the diagonal of a Tensor as a DiagTensor
126 template<class Cmpt>
127 inline DiagTensor<Cmpt> diag(const Tensor<Cmpt>& t)
128 {
129  return DiagTensor<Cmpt>(t.xx(), t.yy(), t.zz());
130 }
131 
132 
133 //- Return the diagonal of a SymmTensor as a DiagTensor
134 template<class Cmpt>
135 inline DiagTensor<Cmpt> diag(const SymmTensor<Cmpt>& st)
136 {
137  return DiagTensor<Cmpt>(st.xx(), st.yy(), st.zz());
138 }
139 
141 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
142 
143 //- Sum of a DiagTensor and a Tensor
144 template<class Cmpt>
145 inline Tensor<Cmpt>
146 operator+(const DiagTensor<Cmpt>& dt1, const Tensor<Cmpt>& t2)
147 {
148  return Tensor<Cmpt>
149  (
150  dt1.xx() + t2.xx(), t2.xy(), t2.xz(),
151  t2.yx(), dt1.yy() + t2.yy(), t2.yz(),
152  t2.zx(), t2.zy(), dt1.zz() + t2.zz()
153  );
154 }
155 
156 
157 //- Sum of a Tensor and a DiagTensor
158 template<class Cmpt>
159 inline Tensor<Cmpt>
160 operator+(const Tensor<Cmpt>& t1, const DiagTensor<Cmpt>& dt2)
161 {
162  return Tensor<Cmpt>
163  (
164  t1.xx() + dt2.xx(), t1.xy(), t1.xz(),
165  t1.yx(), t1.yy() + dt2.yy(), t1.yz(),
166  t1.zx(), t1.zy(), t1.zz() + dt2.zz()
167  );
168 }
170 
171 //- Subtract a Tensor from a DiagTensor
172 template<class Cmpt>
173 inline Tensor<Cmpt>
174 operator-(const DiagTensor<Cmpt>& dt1, const Tensor<Cmpt>& t2)
175 {
176  return Tensor<Cmpt>
177  (
178  dt1.xx() - t2.xx(), -t2.xy(), -t2.xz(),
179  -t2.yx(), dt1.yy() - t2.yy(), -t2.yz(),
180  -t2.zx(), -t2.zy(), dt1.zz() - t2.zz()
181  );
182 }
183 
184 
185 //- Subtract a DiagTensor from a Tensor
186 template<class Cmpt>
187 inline Tensor<Cmpt>
188 operator-(const Tensor<Cmpt>& t1, const DiagTensor<Cmpt>& dt2)
189 {
190  return Tensor<Cmpt>
191  (
192  t1.xx() - dt2.xx(), t1.xy(), t1.xz(),
193  t1.yx(), t1.yy() - dt2.yy(), t1.yz(),
194  t1.zx(), t1.zy(), t1.zz() - dt2.zz()
195  );
196 }
197 
198 
199 //- Division of a Cmpt by a DiagTensor
200 template<class Cmpt>
202 operator/(const Cmpt s, const DiagTensor<Cmpt>& dt)
203 {
204  #ifdef FULLDEBUG
205  if (mag(det(dt)) < VSMALL)
206  {
208  << "Cmpt = " << s
209  << " is not divisible by the DiagTensor due to a zero element:"
210  << "DiagTensor = " << dt
211  << abort(FatalError);
212  }
213  #endif
214 
215  return DiagTensor<Cmpt>(s/dt.xx(), s/dt.yy(), s/dt.zz());
216 }
218 
219 //- Division of a DiagTensor by a Cmpt
220 template<class Cmpt>
221 inline DiagTensor<Cmpt>
222 operator/(const DiagTensor<Cmpt>& dt, const Cmpt s)
223 {
224  #ifdef FULLDEBUG
225  if (mag(s) < VSMALL)
226  {
228  << "DiagTensor = " << dt
229  << " is not divisible due to a zero value in Cmpt:"
230  << "Cmpt = " << s
231  << abort(FatalError);
232  }
233  #endif
234 
235  return DiagTensor<Cmpt>(dt.xx()/s, dt.yy()/s, dt.zz()/s);
236 }
237 
238 
239 //- Division of a Vector by a DiagTensor
240 template<class Cmpt>
241 inline Vector<Cmpt>
242 operator/(const Vector<Cmpt> v, const DiagTensor<Cmpt>& dt)
243 {
244  #ifdef FULLDEBUG
245  if (mag(det(dt)) < VSMALL)
246  {
248  << "Vector = " << v
249  << " is not divisible by the DiagTensor due to a zero element:"
250  << "DiagTensor = " << dt
251  << abort(FatalError);
252  }
253  #endif
254 
255  return Vector<Cmpt>(v.x()/dt.xx(), v.y()/dt.yy(), v.z()/dt.zz());
256 }
257 
258 
259 //- Inner-product of a DiagTensor and a DiagTensor
260 template<class Cmpt>
262 operator&(const DiagTensor<Cmpt>& dt1, const DiagTensor<Cmpt>& dt2)
263 {
264  return DiagTensor<Cmpt>
265  (
266  dt1.xx()*dt2.xx(),
267  dt1.yy()*dt2.yy(),
268  dt1.zz()*dt2.zz()
269  );
270 }
271 
272 
273 //- Inner-product of a DiagTensor and a Tensor
274 template<class Cmpt>
275 inline Tensor<Cmpt>
276 operator&(const DiagTensor<Cmpt>& dt1, const Tensor<Cmpt>& t2)
277 {
278  return Tensor<Cmpt>
279  (
280  dt1.xx()*t2.xx(),
281  dt1.xx()*t2.xy(),
282  dt1.xx()*t2.xz(),
284  dt1.yy()*t2.yx(),
285  dt1.yy()*t2.yy(),
286  dt1.yy()*t2.yz(),
287 
288  dt1.zz()*t2.zx(),
289  dt1.zz()*t2.zy(),
290  dt1.zz()*t2.zz()
291  );
292 }
293 
294 
295 //- Inner-product of a Tensor and a DiagTensor
296 template<class Cmpt>
297 inline Tensor<Cmpt>
298 operator&(const Tensor<Cmpt>& t1, const DiagTensor<Cmpt>& dt2)
299 {
300  return Tensor<Cmpt>
301  (
302  t1.xx()*dt2.xx(),
303  t1.xy()*dt2.yy(),
304  t1.xz()*dt2.zz(),
305 
306  t1.yx()*dt2.xx(),
307  t1.yy()*dt2.yy(),
308  t1.yz()*dt2.zz(),
309 
310  t1.zx()*dt2.xx(),
311  t1.zy()*dt2.yy(),
312  t1.zz()*dt2.zz()
313  );
314 }
315 
316 
317 //- Inner-product of a DiagTensor and a Vector
318 template<class Cmpt>
319 inline Vector<Cmpt>
320 operator&(const DiagTensor<Cmpt>& dt, const Vector<Cmpt>& v)
321 {
322  return Vector<Cmpt>
323  (
324  dt.xx()*v.x(),
325  dt.yy()*v.y(),
326  dt.zz()*v.z()
327  );
328 }
329 
330 
331 //- Inner-product of a Vector and a DiagTensor
332 template<class Cmpt>
333 inline Vector<Cmpt>
334 operator&(const Vector<Cmpt>& v, const DiagTensor<Cmpt>& dt)
335 {
336  return Vector<Cmpt>
337  (
338  v.x()*dt.xx(),
339  v.y()*dt.yy(),
340  v.z()*dt.zz()
341  );
342 }
343 
344 
345 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
346 
347 } // End namespace Foam
348 
349 // ************************************************************************* //
const Cmpt & zz() const noexcept
Definition: DiagTensor.H:128
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:578
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:87
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 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:52
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:127
const Cmpt & x() const noexcept
Access to the vector x component.
Definition: Vector.H:135
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:58
const Cmpt & xx() const noexcept
Definition: DiagTensor.H:126
A templated (3 x 3) tensor of objects of <T> derived from MatrixSpace.
Definition: complexI.H:268
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
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157