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-2022 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
30#include "tensorFieldField.H"
31
32#define TEMPLATE template<template<class> class PatchField, class GeoMesh>
34
35// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
36
37template<class Cmpt, template<class> class PatchField, class GeoMesh>
39(
40 GeometricField<Tensor<Cmpt>, PatchField, GeoMesh>& result,
50)
51{
53 (
54 result.primitiveFieldRef(),
58 );
59
61 (
62 result.boundaryFieldRef(),
66 );
67}
68
69
70template<class Cmpt, template<class> class PatchField, class GeoMesh>
72(
73 const GeometricField<Tensor<Cmpt>, PatchField, GeoMesh>& input,
83)
84{
86 (
87 input.primitiveField(),
91 );
92
94 (
95 input.boundaryField(),
99 );
100}
101
102
103template<class Cmpt, template<class> class PatchField, class GeoMesh>
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
130template<class Cmpt, template<class> class PatchField, class GeoMesh>
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
157template<class Cmpt, template<class> class PatchField, class GeoMesh>
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
184template<class Cmpt, template<class> class PatchField, class GeoMesh>
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
211template<class Cmpt, template<class> class PatchField, class GeoMesh>
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());
220
221 Foam::unzipRow(input.boundaryField(), idx, result.boundaryFieldRef());
222}
223
224
225template<class Cmpt, template<class> class PatchField, class GeoMesh>
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());
234
235 Foam::unzipCol(input.boundaryField(), idx, result.boundaryFieldRef());
236}
237
238
239template<class Cmpt, template<class> class PatchField, class GeoMesh>
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
254namespace Foam
255{
256
257// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258
267UNARY_FUNCTION(scalar, tensor, det, pow3)
270
273
274
275// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
276
279
282
283
284// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
285
286} // End namespace Foam
287
288// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
289
290#include "undefFieldFunctionsM.H"
291
292// ************************************************************************* //
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define UNARY_OPERATOR(ReturnType, Type1, Op, OpFunc, Dfunc)
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
Tensor specific part of the implementation of GeometricField.
scalar y
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:49
Generic GeometricField class.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A templated (3 x 3) tensor of objects of <T> derived from MatrixSpace.
Definition: Tensor.H:64
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:65
Tensor of scalars, i.e. Tensor<scalar>.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const volScalarField & T
Namespace for OpenFOAM.
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar det(const dimensionedSphericalTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
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.
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
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.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
void unzipDiag(const FieldField< Field, SymmTensor< Cmpt > > &input, FieldField< Field, Vector< Cmpt > > &result)
Extract a symmTensor field field diagonal.
dimensionedTensor eigenVectors(const dimensionedSymmTensor &dt)
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedVector eigenValues(const dimensionedSymmTensor &dt)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:59
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)
dimensionedSymmTensor cof(const dimensionedSymmTensor &dt)
void hdual(pointPatchField< vector > &, const pointPatchField< tensor > &)
SphericalTensor< Cmpt > sph(const DiagTensor< Cmpt > &dt)
Return the spherical part of a DiagTensor as a SphericalTensor.
Definition: DiagTensorI.H:130
uint8_t direction
Definition: direction.H:56
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 unzip(const FieldField< Field, SphericalTensor< Cmpt > > &input, FieldField< Field, Cmpt > &ii)
Unzip sphericalTensor field field into components.
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
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)
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.
dimensionSet pow2(const dimensionSet &ds)
Definition: dimensionSet.C:369
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
dimensionedTensor skew(const dimensionedTensor &dt)
void zip(FieldField< Field, SphericalTensor< Cmpt > > &result, const FieldField< Field, Cmpt > &ii)
Zip together sphericalTensor field field from components.
Specialisation of FieldField<T> for tensor.