SphericalTensor2DI.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 "Vector2D.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class Cmpt>
35:
37{}
38
39
40template<class Cmpt>
42(
43 const VectorSpace<SphericalTensor2D<Cmpt>, Cmpt, 1>& vs
44)
45:
47{}
48
49
50template<class Cmpt>
52{
53 this->v_[II] = stii;
54}
55
56
57template<class Cmpt>
59:
61{}
62
63
64// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65
66template<class Cmpt>
67inline const Cmpt& Foam::SphericalTensor2D<Cmpt>::ii() const
68{
69 return this->v_[II];
70}
71
72
73template<class Cmpt>
75{
76 return this->v_[II];
77}
78
79
80// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
81
82namespace Foam
83{
84
85// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
86
87//- Return the trace of a SphericalTensor2D
88template<class Cmpt>
89inline Cmpt tr(const SphericalTensor2D<Cmpt>& st)
90{
91 return 2*st.ii();
92}
93
94
95//- Return the spherical part of a SphericalTensor2D, i.e. itself
96template<class Cmpt>
98{
99 return st;
100}
101
102
103//- Return the determinant of a SphericalTensor2D
104template<class Cmpt>
105inline Cmpt det(const SphericalTensor2D<Cmpt>& st)
106{
107 return st.ii()*st.ii();
108}
109
110
111//- Return the inverse of a SphericalTensor2D
112template<class Cmpt>
114{
115 #ifdef FULLDEBUG
116 if (mag(st.ii()) < VSMALL)
117 {
119 << "SphericalTensor2D is not invertible due to zero determinant:"
120 << "det(SphericalTensor2D) = " << det(st)
121 << abort(FatalError);
122 }
123 #endif
124
125 return SphericalTensor2D<Cmpt>(1/st.ii());
126}
127
128
129// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
130
131//- Division of a Cmpt by a SphericalTensor2D
132template<class Cmpt>
133inline SphericalTensor2D<Cmpt>
134operator/(const Cmpt s, const SphericalTensor2D<Cmpt>& st)
135{
136 #ifdef FULLDEBUG
137 if (mag(st.ii()) < VSMALL)
138 {
140 << "Cmpt = " << s
141 << " is not divisible due to a zero element in SphericalTensor2D:"
142 << "SphericalTensor2D = " << st
143 << abort(FatalError);
144 }
145 #endif
146
147 return SphericalTensor2D<Cmpt>(s/st.ii());
148}
149
150
151//- Division of a SphericalTensor2D by a Cmpt
152template<class Cmpt>
153inline SphericalTensor2D<Cmpt>
154operator/(const SphericalTensor2D<Cmpt>& st, const Cmpt s)
155{
156 #ifdef FULLDEBUG
157 if (mag(s) < VSMALL)
158 {
160 << "SphericalTensor2D = " << st
161 << " is not divisible due to a zero value in Cmpt:"
162 << "Cmpt = " << s
163 << abort(FatalError);
164 }
165 #endif
166
167 return SphericalTensor2D<Cmpt>(st.ii()/s);
168}
169
170
171//- Inner-product of a SphericalTensor2D and a SphericalTensor2D
172template<class Cmpt>
173inline SphericalTensor2D<Cmpt>
174operator&
175(
176 const SphericalTensor2D<Cmpt>& st1,
177 const SphericalTensor2D<Cmpt>& st2
178)
179{
180 return SphericalTensor2D<Cmpt>(st1.ii()*st2.ii());
181}
182
183
184//- Inner-product of a SphericalTensor2D and a Vector2D
185template<class Cmpt>
186inline Vector2D<Cmpt>
188{
189 return Vector2D<Cmpt>
190 (
191 st.ii()*v.x(),
192 st.ii()*v.y()
193 );
194}
195
196
197//- Inner-product of a Vector2D and a SphericalTensor2D
198template<class Cmpt>
199inline Vector2D<Cmpt>
201{
202 return Vector2D<Cmpt>
203 (
204 v.x()*st.ii(),
205 v.y()*st.ii()
206 );
207}
208
209
210// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211
212template<class Cmpt>
214{
215public:
216
218};
219
220template<class Cmpt>
222{
223public:
224
226};
227
228
229template<class Cmpt>
231{
232public:
233
235};
236
237
238template<class Cmpt>
240{
241public:
242
244};
245
246template<class Cmpt>
248{
249public:
250
252};
253
254
255// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256
257} // End namespace Foam
258
259// ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
A templated (2 x 2) diagonal tensor of objects of <T>, effectively containing 1 element,...
SphericalTensor2D()=default
Default construct.
const Cmpt & ii() const
Templated 2D Vector derived from VectorSpace adding construction from 2 components,...
Definition: Vector2D.H:58
const Cmpt & y() const
Access to the vector y component.
Definition: Vector2DI.H:73
const Cmpt & x() const
Access to the vector x component.
Definition: Vector2DI.H:66
Templated vector space.
Definition: VectorSpace.H:79
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)
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)