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 | Copyright (C) 2013-2016 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 "Moment.H"
27 
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const IOobject& io,
35  const dictionary& dict,
36  const fvMesh& mesh
37 )
38 :
39  AveragingMethod<Type>(io, dict, mesh, labelList(4, mesh.nCells())),
44  transform_(mesh.nCells(), Zero),
45  scale_(0.5*pow(mesh.V(), 1.0/3.0))
46 {
47  scalar a = 1.0/24.0;
48  scalar b = 0.5854101966249685;
49  scalar c = 0.1381966011250105;
50 
51  scalarField wQ(4);
52  wQ[0] = a;
53  wQ[1] = a;
54  wQ[2] = a;
55  wQ[3] = a;
56 
57  vectorField xQ(4);
58  xQ[0] = vector(b, c, c);
59  xQ[1] = vector(c, b, c);
60  xQ[2] = vector(c, c, b);
61  xQ[3] = vector(c, c, c);
62 
63  forAll(mesh.C(), celli)
64  {
65  const List<tetIndices> cellTets =
66  polyMeshTetDecomposition::cellTetIndices(mesh, celli);
67 
68  symmTensor A(Zero);
69 
70  forAll(cellTets, tetI)
71  {
72  const tetIndices& tetIs = cellTets[tetI];
73  const label facei = tetIs.face();
74  const face& f = mesh.faces()[facei];
75 
76  const tensor T
77  (
78  tensor
79  (
80  mesh.points()[f[tetIs.faceBasePt()]] - mesh.C()[celli],
81  mesh.points()[f[tetIs.facePtA()]] - mesh.C()[celli],
82  mesh.points()[f[tetIs.facePtB()]] - mesh.C()[celli]
83  ).T()
84  );
85 
86  const vectorField d((T & xQ)/scale_[celli]);
87 
88  const scalar v(6.0*tetIs.tet(mesh).mag()/mesh.V()[celli]);
89 
90  A += v*sum(wQ*sqr(d));
91  }
92 
93  transform_[celli] = inv(A);
94  }
95 }
96 
97 
98 template<class Type>
100 (
101  const Moment<Type>& am
102 )
103 :
109  transform_(am.transform_)
110 {}
111 
112 
113 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
114 
115 template<class Type>
117 {}
118 
119 
120 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
121 
122 template<class Type>
124 {}
125 
126 
127 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128 
129 template<class Type>
131 (
132  const point position,
133  const tetIndices& tetIs,
134  const Type& value
135 )
136 {
137  const label celli = tetIs.cell();
138 
139  const Type v = value/this->mesh_.V()[celli];
140  const TypeGrad dv =
141  transform_[celli]
142  & (
143  v
144  * (position - this->mesh_.C()[celli])
145  / scale_[celli]
146  );
147 
148  data_[celli] += v;
149  dataX_[celli] += v + dv.x();
150  dataY_[celli] += v + dv.y();
151  dataZ_[celli] += v + dv.z();
152 }
153 
154 
155 template<class Type>
157 (
158  const point position,
159  const tetIndices& tetIs
160 ) const
161 {
162  const label celli = tetIs.cell();
163 
164  return
165  data_[celli]
166  + (
167  TypeGrad
168  (
169  dataX_[celli] - data_[celli],
170  dataY_[celli] - data_[celli],
171  dataZ_[celli] - data_[celli]
172  )
173  & (position - this->mesh_.C()[celli])
174  / scale_[celli]
175  );
176 }
177 
178 
179 template<class Type>
182 (
183  const point position,
184  const tetIndices& tetIs
185 ) const
186 {
187  const label celli(tetIs.cell());
188 
189  return
190  TypeGrad
191  (
192  dataX_[celli] - data_[celli],
193  dataY_[celli] - data_[celli],
194  dataZ_[celli] - data_[celli]
195  )/scale_[celli];
196 }
197 
198 
199 template<class Type>
202 {
203  return tmp<Field<Type>>(data_);
204 }
205 
206 
207 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
label cell() const
Return the cell.
Definition: tetIndicesI.H:30
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
label facePtA() const
Return face point A.
Definition: tetIndicesI.H:48
TypeGrad interpolateGrad(const point position, const tetIndices &tetIs) const
Interpolate gradient.
Definition: Moment.C:182
tmp< Field< Type > > primitiveField() const
Return an internal field of the average.
Definition: Moment.C:201
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
void add(const point position, const tetIndices &tetIs, const Type &value)
Member Functions.
Definition: Moment.C:131
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
AveragingMethod< Type >::TypeGrad TypeGrad
Public typedefs.
Definition: Moment.H:70
outerProduct< vector, Type >::type TypeGrad
Protected typedefs.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Generic field type.
Definition: FieldField.H:51
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
label faceBasePt() const
Return the face base point.
Definition: tetIndicesI.H:42
Moment lagrangian averaging procedure.
Definition: Moment.H:61
dynamicFvMesh & mesh
label facePtB() const
Return face point B.
Definition: tetIndicesI.H:54
Type interpolate(const point position, const tetIndices &tetIs) const
Interpolate.
Definition: Moment.C:157
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Base class for lagrangian averaging methods.
Definition: MPPICParcel.H:55
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:91
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
label face() const
Return the face.
Definition: tetIndicesI.H:36
scalar mag() const
Return volume.
Definition: tetrahedronI.H:171
labelList f(nPoints)
const volScalarField & T
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionedScalar c
Speed of light in a vacuum.
Moment(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Constructors.
Definition: Moment.C:33
virtual ~Moment()
Destructor.
Definition: Moment.C:116
A class for managing temporary objects.
Definition: PtrList.H:54
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet from the.
Definition: tetIndicesI.H:66