SphericalTensorI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "Vector.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class Cmpt>
32 {}
33 
34 
35 template<class Cmpt>
37 :
39 {}
40 
41 
42 template<class Cmpt>
43 template<class Cmpt2>
45 (
46  const VectorSpace<SphericalTensor<Cmpt2>, Cmpt2, 1>& vs
47 )
48 :
50 {}
51 
52 
53 template<class Cmpt>
55 {
56  this->v_[II] = stii;
57 }
58 
59 
60 template<class Cmpt>
62 :
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
69 template<class Cmpt>
70 inline const Cmpt& Foam::SphericalTensor<Cmpt>::ii() const
71 {
72  return this->v_[II];
73 }
74 
75 
76 template<class Cmpt>
78 {
79  return this->v_[II];
80 }
81 
82 
83 template<class Cmpt>
84 inline const Foam::SphericalTensor<Cmpt>&
86 {
87  return *this;
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
92 
93 namespace Foam
94 {
95 
96 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
97 
98 //- Inner-product between two spherical tensors
99 template<class Cmpt>
102 {
103  return SphericalTensor<Cmpt>(st1.ii()*st2.ii());
104 }
105 
106 
107 //- Inner-product between a spherical tensor and a vector
108 template<class Cmpt>
109 inline Vector<Cmpt>
111 {
112  return Vector<Cmpt>
113  (
114  st.ii()*v.x(),
115  st.ii()*v.y(),
116  st.ii()*v.z()
117  );
118 }
119 
120 
121 //- Inner-product between a vector and a spherical tensor
122 template<class Cmpt>
123 inline Vector<Cmpt>
125 {
126  return Vector<Cmpt>
127  (
128  v.x()*st.ii(),
129  v.y()*st.ii(),
130  v.z()*st.ii()
131  );
132 }
133 
134 
135 //- Double-dot-product between a spherical tensor and a spherical tensor
136 template<class Cmpt>
137 inline Cmpt
139 {
140  return 3*st1.ii()*st2.ii();
141 }
142 
143 
144 //- Division of a scalar by a sphericalTensor
145 template<class Cmpt>
146 inline SphericalTensor<Cmpt>
147 operator/(const scalar s, const SphericalTensor<Cmpt>& st)
148 {
149  return SphericalTensor<Cmpt>(s/st.ii());
150 }
151 
152 
153 template<class Cmpt>
154 inline Cmpt magSqr(const SphericalTensor<Cmpt>& st)
155 {
156  return 3*magSqr(st.ii());
157 }
158 
159 
160 //- Return the trace of a spherical tensor
161 template<class Cmpt>
162 inline Cmpt tr(const SphericalTensor<Cmpt>& st)
163 {
164  return 3*st.ii();
165 }
166 
167 
168 //- Return the spherical part of a spherical tensor, i.e. itself
169 template<class Cmpt>
171 {
172  return st;
173 }
174 
175 
176 //- Return the determinant of a spherical tensor
177 template<class Cmpt>
178 inline Cmpt det(const SphericalTensor<Cmpt>& st)
179 {
180  return st.ii()*st.ii()*st.ii();
181 }
182 
183 
184 //- Return the inverse of a spherical tensor
185 template<class Cmpt>
187 {
188  return SphericalTensor<Cmpt>(1.0/st.ii());
189 }
190 
191 
192 template<class Cmpt>
193 class outerProduct<SphericalTensor<Cmpt>, Cmpt>
194 {
195 public:
196 
198 };
199 
200 template<class Cmpt>
201 class outerProduct<Cmpt, SphericalTensor<Cmpt>>
202 {
203 public:
204 
206 };
207 
208 
209 template<class Cmpt>
211 {
212 public:
213 
215 };
216 
217 
218 template<class Cmpt>
220 {
221 public:
222 
224 };
225 
226 template<class Cmpt>
228 {
229 public:
230 
232 };
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 } // End namespace Foam
238 
239 // ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
Templated 3D SphericalTensor derived from VectorSpace adding construction from 1 component,...
const Cmpt & ii() const
const SphericalTensor< Cmpt > & T() const
Transpose.
SphericalTensor()
Construct null.
Templated vector space.
Definition: VectorSpace.H:85
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:60
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:50
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
dimensionedScalar det(const dimensionedSphericalTensor &dt)
static const zero Zero
Definition: zero.H:97
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
SphericalTensor< Cmpt > sph(const DiagTensor< Cmpt > &dt)
Return the spherical part of a diagonal tensor.
Definition: DiagTensorI.H:353
tmp< VolField< Type > > operator&(const fvMatrix< Type > &, const DimensionedField< Type, volMesh > &)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
tmp< fvMatrix< Type > > operator/(const fvMatrix< Type > &, const volScalarField::Internal &)
dimensioned< typename scalarProduct< Type1, Type2 >::type > operator&&(const dimensioned< Type1 > &, const dimensioned< Type2 > &)