Moment.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) 2013-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "Moment.H"
29
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class Type>
35(
36 const IOobject& io,
37 const dictionary& dict,
38 const fvMesh& mesh
39)
40:
41 AveragingMethod<Type>(io, dict, mesh, labelList(4, mesh.nCells())),
42 data_(FieldField<Field, Type>::operator[](0)),
43 dataX_(FieldField<Field, Type>::operator[](1)),
44 dataY_(FieldField<Field, Type>::operator[](2)),
45 dataZ_(FieldField<Field, Type>::operator[](3)),
46 transform_(mesh.nCells(), Zero),
47 scale_(0.5*cbrt(mesh.V()))
48{
49 scalar a = 1.0/24.0;
50 scalar b = 0.5854101966249685;
51 scalar c = 0.1381966011250105;
52
53 scalarField wQ(4);
54 wQ[0] = a;
55 wQ[1] = a;
56 wQ[2] = a;
57 wQ[3] = a;
58
59 vectorField xQ(4);
60 xQ[0] = vector(b, c, c);
61 xQ[1] = vector(c, b, c);
62 xQ[2] = vector(c, c, b);
63 xQ[3] = vector(c, c, c);
64
65 forAll(mesh.C(), celli)
66 {
67 const List<tetIndices> cellTets =
69
71
72 forAll(cellTets, tetI)
73 {
74 const tetIndices& tetIs = cellTets[tetI];
75 const triFace triIs = tetIs.faceTriIs(mesh);
76
77 const tensor T
78 (
79 tensor
80 (
81 mesh.points()[triIs[0]] - mesh.C()[celli],
82 mesh.points()[triIs[1]] - mesh.C()[celli],
83 mesh.points()[triIs[2]] - mesh.C()[celli]
84 ).T()
85 );
86
87 const vectorField d((T & xQ)/scale_[celli]);
88
89 const scalar v(6.0*tetIs.tet(mesh).mag()/mesh.V()[celli]);
90
91 A += v*sum(wQ*sqr(d));
92 }
93
94 transform_[celli] = inv(A);
95 }
96}
97
98
99template<class Type>
101(
102 const Moment<Type>& am
103)
104:
105 AveragingMethod<Type>(am),
106 data_(FieldField<Field, Type>::operator[](0)),
107 dataX_(FieldField<Field, Type>::operator[](1)),
108 dataY_(FieldField<Field, Type>::operator[](2)),
109 dataZ_(FieldField<Field, Type>::operator[](3)),
110 transform_(am.transform_)
111{}
112
113
114// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
115
116template<class Type>
118{}
119
120
121// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
122
123template<class Type>
125{}
126
127
128// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129
130template<class Type>
132(
134 const tetIndices& tetIs,
135 const Type& value
136)
137{
138 const label celli = tetIs.cell();
139 const triFace triIs = tetIs.faceTriIs(this->mesh_);
140
141 const point delta =
142 (coordinates[0] - 1)*this->mesh_.C()[celli]
143 + coordinates[1]*this->mesh_.points()[triIs[0]]
144 + coordinates[2]*this->mesh_.points()[triIs[1]]
145 + coordinates[3]*this->mesh_.points()[triIs[2]];
146
147 const Type v = value/this->mesh_.V()[celli];
148 const TypeGrad dv = transform_[celli] & (v*delta/scale_[celli]);
149
150 data_[celli] += v;
151 dataX_[celli] += v + dv.x();
152 dataY_[celli] += v + dv.y();
153 dataZ_[celli] += v + dv.z();
154}
155
156
157template<class Type>
159(
161 const tetIndices& tetIs
162) const
163{
164 const label celli = tetIs.cell();
165 const triFace triIs = tetIs.faceTriIs(this->mesh_);
166
167 const point delta =
168 (coordinates[0] - 1)*this->mesh_.C()[celli]
169 + coordinates[1]*this->mesh_.points()[triIs[0]]
170 + coordinates[2]*this->mesh_.points()[triIs[1]]
171 + coordinates[3]*this->mesh_.points()[triIs[2]];
172
173 return
174 data_[celli]
175 + (
177 (
178 dataX_[celli] - data_[celli],
179 dataY_[celli] - data_[celli],
180 dataZ_[celli] - data_[celli]
181 )
182 & delta/scale_[celli]
183 );
184}
185
186
187template<class Type>
190(
192 const tetIndices& tetIs
193) const
194{
195 const label celli(tetIs.cell());
196
197 return
199 (
200 dataX_[celli] - data_[celli],
201 dataY_[celli] - data_[celli],
202 dataZ_[celli] - data_[celli]
203 )/scale_[celli];
204}
205
206
207template<class Type>
210{
211 return tmp<Field<Type>>(data_);
212}
213
214
215// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
scalar delta
Base class for lagrangian averaging methods.
Moment lagrangian averaging procedure.
Definition: Moment.H:69
tmp< Field< Type > > primitiveField() const
Return an internal field of the average.
Definition: Moment.C:209
virtual ~Moment()
Destructor.
Definition: Moment.C:117
TypeGrad interpolateGrad(const barycentric &coordinates, const tetIndices &tetIs) const
Interpolate gradient.
Definition: Moment.C:190
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
tmp< FieldField< Field, Type > > T() const
Return the field transpose (only defined for second rank tensors)
Definition: FieldField.C:286
Generic templated field type.
Definition: Field.H:82
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Tensor< Cmpt > T() const
Return non-Hermitian transpose.
Definition: TensorI.H:526
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Sums a given list of (at least two or more) fields and outputs the result into a new field,...
Definition: add.H:161
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const volVectorField & C() const
Return cell centres as volVectorField.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
static List< tetIndices > cellTetIndices(const polyMesh &mesh, label cI)
Return the tet decomposition of the given cell, see.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
bool interpolate() const noexcept
Same as isPointData()
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:84
triFace faceTriIs(const polyMesh &mesh, const bool warn=true) const
Return the indices corresponding to the tri on the face for.
Definition: tetIndicesI.H:55
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet.
Definition: tetIndicesI.H:134
label cell() const noexcept
Return the cell index.
Definition: tetIndices.H:133
scalar mag() const
Return volume.
Definition: tetrahedronI.H:172
A class for managing temporary objects.
Definition: tmp.H:65
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:72
type
Volume classification types.
Definition: volumeType.H:66
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensionedScalar cbrt(const dimensionedScalar &ds)
Vector< scalar > vector
Definition: vector.H:61
dictionary dict
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333