GeometricTensorField.C
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-2013 OpenFOAM Foundation
9  Copyright (C) 2019-2023 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 "GeometricTensorField.H"
30 #include "tensorFieldField.H"
31 
32 #define TEMPLATE template<template<class> class PatchField, class GeoMesh>
34 
35 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
36 
37 template<class Cmpt, template<class> class PatchField, class GeoMesh>
38 void Foam::zip
39 (
40  GeometricField<Tensor<Cmpt>, PatchField, GeoMesh>& result,
50 )
51 {
52  Foam::zip
53  (
54  result.primitiveFieldRef(),
58  );
59 
60  Foam::zip
61  (
62  result.boundaryFieldRef(),
63  xx.boundaryField(), xy.boundaryField(), xz.boundaryField(),
64  yx.boundaryField(), yy.boundaryField(), yz.boundaryField(),
66  );
67 }
68 
69 
70 template<class Cmpt, template<class> class PatchField, class GeoMesh>
71 void Foam::unzip
72 (
73  const GeometricField<Tensor<Cmpt>, PatchField, GeoMesh>& input,
74  GeometricField<Cmpt, PatchField, GeoMesh>& xx,
75  GeometricField<Cmpt, PatchField, GeoMesh>& xy,
76  GeometricField<Cmpt, PatchField, GeoMesh>& xz,
77  GeometricField<Cmpt, PatchField, GeoMesh>& yx,
78  GeometricField<Cmpt, PatchField, GeoMesh>& yy,
79  GeometricField<Cmpt, PatchField, GeoMesh>& yz,
80  GeometricField<Cmpt, PatchField, GeoMesh>& zx,
81  GeometricField<Cmpt, PatchField, GeoMesh>& zy,
82  GeometricField<Cmpt, PatchField, GeoMesh>& zz
83 )
84 {
86  (
87  input.primitiveField(),
88  xx.primitiveFieldRef(), xy.primitiveFieldRef(), xz.primitiveFieldRef(),
89  yx.primitiveFieldRef(), yy.primitiveFieldRef(), yz.primitiveFieldRef(),
90  zx.primitiveFieldRef(), zy.primitiveFieldRef(), zz.primitiveFieldRef()
91  );
92 
94  (
95  input.boundaryField(),
96  xx.boundaryFieldRef(), xy.boundaryFieldRef(), xz.boundaryFieldRef(),
97  yx.boundaryFieldRef(), yy.boundaryFieldRef(), yz.boundaryFieldRef(),
98  zx.boundaryFieldRef(), zy.boundaryFieldRef(), zz.boundaryFieldRef()
99  );
100 }
101 
102 
103 template<class Cmpt, template<class> class PatchField, class GeoMesh>
104 void Foam::zipRows
105 (
106  GeometricField<Tensor<Cmpt>, PatchField, GeoMesh>& result,
107  const GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& x,
108  const GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& y,
109  const GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& z
110 )
111 {
113  (
114  result.primitiveFieldRef(),
115  x.primitiveField(),
116  y.primitiveField(),
117  z.primitiveField()
118  );
119 
121  (
122  result.boundaryFieldRef(),
123  x.boundaryField(),
124  y.boundaryField(),
125  z.boundaryField()
126  );
127 }
128 
129 
130 template<class Cmpt, template<class> class PatchField, class GeoMesh>
131 void Foam::zipCols
132 (
133  GeometricField<Tensor<Cmpt>, PatchField, GeoMesh>& result,
134  const GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& x,
135  const GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& y,
136  const GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& z
137 )
138 {
140  (
141  result.primitiveFieldRef(),
142  x.primitiveField(),
143  y.primitiveField(),
144  z.primitiveField()
145  );
146 
148  (
149  result.boundaryFieldRef(),
150  x.boundaryField(),
151  y.boundaryField(),
152  z.boundaryField()
153  );
154 }
155 
156 
157 template<class Cmpt, template<class> class PatchField, class GeoMesh>
158 void Foam::unzipRows
159 (
160  const GeometricField<Tensor<Cmpt>, PatchField, GeoMesh>& input,
161  GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& x,
162  GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& y,
163  GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& z
164 )
165 {
167  (
168  input.primitiveField(),
169  x.primitiveFieldRef(),
170  y.primitiveFieldRef(),
171  z.primitiveFieldRef()
172  );
173 
175  (
176  input.boundaryField(),
177  x.boundaryFieldRef(),
178  y.boundaryFieldRef(),
179  z.boundaryFieldRef()
180  );
181 }
182 
183 
184 template<class Cmpt, template<class> class PatchField, class GeoMesh>
185 void Foam::unzipCols
186 (
187  const GeometricField<Tensor<Cmpt>, PatchField, GeoMesh>& input,
188  GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& x,
189  GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& y,
190  GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& z
191 )
192 {
194  (
195  input.primitiveField(),
196  x.primitiveFieldRef(),
197  y.primitiveFieldRef(),
198  z.primitiveFieldRef()
199  );
200 
202  (
203  input.boundaryField(),
204  x.boundaryFieldRef(),
205  y.boundaryFieldRef(),
206  z.boundaryFieldRef()
207  );
208 }
209 
210 
211 template<class Cmpt, template<class> class PatchField, class GeoMesh>
212 void Foam::unzipRow
213 (
214  const GeometricField<Tensor<Cmpt>, PatchField, GeoMesh>& input,
215  const direction idx,
216  GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& result
217 )
218 {
219  Foam::unzipRow(input.primitiveField(), idx, result.primitiveFieldRef());
221  Foam::unzipRow(input.boundaryField(), idx, result.boundaryFieldRef());
222 }
223 
224 
225 template<class Cmpt, template<class> class PatchField, class GeoMesh>
226 void Foam::unzipCol
227 (
228  const GeometricField<Tensor<Cmpt>, PatchField, GeoMesh>& input,
229  const direction idx,
230  GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& result
231 )
232 {
233  Foam::unzipCol(input.primitiveField(), idx, result.primitiveFieldRef());
235  Foam::unzipCol(input.boundaryField(), idx, result.boundaryFieldRef());
236 }
237 
238 
239 template<class Cmpt, template<class> class PatchField, class GeoMesh>
240 void Foam::unzipDiag
241 (
242  const GeometricField<Tensor<Cmpt>, PatchField, GeoMesh>& input,
243  GeometricField<Vector<Cmpt>, PatchField, GeoMesh>& result
244 )
245 {
246  Foam::unzipDiag(input.primitiveField(), result.primitiveFieldRef());
247 
248  Foam::unzipDiag(input.boundaryField(), result.boundaryFieldRef());
249 }
250 
251 
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 
254 namespace Foam
255 {
256 
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258 
269 UNARY_FUNCTION(scalar, tensor, det, pow3)
272 
275 
276 
277 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
278 
281 
284 
285 
286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 
288 } // End namespace Foam
289 
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 
292 #include "undefFieldFunctionsM.H"
293 
294 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
void divide(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
uint8_t direction
Definition: direction.H:46
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
Specialisation of FieldField<T> for tensor.
void unzipCols(const FieldField< Field, SymmTensor< Cmpt >> &input, FieldField< Field, Vector< Cmpt >> &x, FieldField< Field, Vector< Cmpt >> &y, FieldField< Field, Vector< Cmpt >> &z)
Extract symmTensor field field columns.
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
dimensionedTensor skew(const dimensionedTensor &dt)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
void zip(FieldField< Field, SphericalTensor< Cmpt >> &result, const FieldField< Field, Cmpt > &ii)
Zip together sphericalTensor field field from components.
void zipCols(FieldField< Field, SymmTensor< Cmpt >> &result, const FieldField< Field, Vector< Cmpt >> &x, const FieldField< Field, Vector< Cmpt >> &y, const FieldField< Field, Vector< Cmpt >> &z)
Zip together symmTensor field from column components.
Tensor< scalar > tensor
Definition: symmTensor.H:57
void unzipRow(const FieldField< Field, SymmTensor< Cmpt >> &input, const direction idx, FieldField< Field, Vector< Cmpt >> &result)
Extract a symmTensor field field row (x,y,z) == (0,1,2)
SymmTensor< Cmpt > devSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:481
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
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:110
void unzipRows(const FieldField< Field, SymmTensor< Cmpt >> &input, FieldField< Field, Vector< Cmpt >> &x, FieldField< Field, Vector< Cmpt >> &y, FieldField< Field, Vector< Cmpt >> &z)
Extract symmTensor field field rows.
void zipRows(FieldField< Field, SymmTensor< Cmpt >> &result, const FieldField< Field, Vector< Cmpt >> &x, const FieldField< Field, Vector< Cmpt >> &y, const FieldField< Field, Vector< Cmpt >> &z)
Zip together symmTensor field field from row components.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
dimensionedTensor eigenVectors(const dimensionedSymmTensor &dt)
void unzip(const FieldField< Field, SphericalTensor< Cmpt >> &input, FieldField< Field, Cmpt > &ii)
Unzip sphericalTensor field field into components.
void unzipDiag(const FieldField< Field, SymmTensor< Cmpt >> &input, FieldField< Field, Vector< Cmpt >> &result)
Extract a symmTensor field field diagonal.
scalar y
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:55
static Istream & input(Istream &is, IntRange< T > &range)
Definition: IntRanges.C:33
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Vector< scalar > vector
Definition: vector.H:57
dimensionSet pow2(const dimensionSet &ds)
Definition: dimensionSet.C:352
void hdual(pointPatchField< vector > &, const pointPatchField< tensor > &)
Tensor specific part of the implementation of GeometricField.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void unzipCol(const FieldField< Field, SymmTensor< Cmpt >> &input, const direction idx, FieldField< Field, Vector< Cmpt >> &result)
Extract a symmTensor field field column (x,y,z) == (0,1,2)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedVector eigenValues(const dimensionedSymmTensor &dt)
#define UNARY_OPERATOR(ReturnType, Type1, Op, OpFunc, Dfunc)
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:42
A templated (3 x 3) tensor of objects of <T> derived from MatrixSpace.
Definition: complexI.H:266
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
dimensionedSymmTensor cof(const dimensionedSymmTensor &dt)
Namespace for OpenFOAM.
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491