AveragingMethod.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  Copyright (C) 2019-2021 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 "AveragingMethod.H"
30 #include "runTimeSelectionTables.H"
31 #include "pointMesh.H"
32 
33 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
34 
35 template<class Type>
37 {}
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class Type>
44 (
45  const IOobject& io,
46  const dictionary& dict,
47  const fvMesh& mesh,
48  const labelList& size
49 )
50 :
51  regIOobject(io),
53  dict_(dict),
54  mesh_(mesh)
55 {
56  forAll(size, i)
57  {
59  (
60  new Field<Type>(size[i], Zero)
61  );
62  }
63 }
64 
65 
66 template<class Type>
68 (
69  const AveragingMethod<Type>& am
70 )
71 :
72  regIOobject(am),
74  dict_(am.dict_),
75  mesh_(am.mesh_)
76 {}
77 
78 
79 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
80 
81 template<class Type>
84 (
85  const IOobject& io,
86  const dictionary& dict,
87  const fvMesh& mesh
88 )
89 {
90  const word modelType
91  (
92  dict.template getOrDefault<word>(typeName, "basic")
93  );
94 
95  auto* ctorPtr = dictionaryConstructorTable(modelType);
96 
97  if (!ctorPtr)
98  {
100  (
101  dict,
102  "averaging limiter",
103  modelType,
104  *dictionaryConstructorTablePtr_
105  ) << abort(FatalIOError);
106  }
107 
108  return autoPtr<AveragingMethod<Type>>(ctorPtr(io, dict, mesh));
109 }
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
114 template<class Type>
116 {
117  updateGrad();
118 }
119 
120 
121 template<class Type>
123 (
124  const AveragingMethod<scalar>& weight
125 )
126 {
127  updateGrad();
128 
129  *this /= max(weight, SMALL);
130 }
131 
132 
133 template<class Type>
135 {
136  return os.good();
137 }
138 
139 
140 template<class Type>
141 bool Foam::AveragingMethod<Type>::write(const bool valid) const
142 {
143  const pointMesh pointMesh_(mesh_);
144 
145  // point volumes
146  Field<scalar> pointVolume(mesh_.nPoints(), Zero);
147 
148  // output fields
150  (
151  IOobject
152  (
153  this->name() + ":cellValue",
154  this->time().timeName(),
155  mesh_
156  ),
157  mesh_,
159  );
161  (
162  IOobject
163  (
164  this->name() + ":cellGrad",
165  this->time().timeName(),
166  mesh_
167  ),
168  mesh_,
170  );
172  (
173  IOobject
174  (
175  this->name() + ":pointValue",
176  this->time().timeName(),
177  mesh_
178  ),
179  pointMesh_,
181  );
183  (
184  IOobject
185  (
186  this->name() + ":pointGrad",
187  this->time().timeName(),
188  mesh_
189  ),
190  pointMesh_,
192  );
193 
194  // Barycentric coordinates of the tet vertices
196  tetCrds
197  ({
198  barycentric(1, 0, 0, 0),
199  barycentric(0, 1, 0, 0),
200  barycentric(0, 0, 1, 0),
201  barycentric(0, 0, 0, 1)
202  });
203 
204  // tet-volume weighted sums
205  forAll(mesh_.C(), celli)
206  {
207  const List<tetIndices> cellTets =
208  polyMeshTetDecomposition::cellTetIndices(mesh_, celli);
209 
210  forAll(cellTets, tetI)
211  {
212  const tetIndices& tetIs = cellTets[tetI];
213  const triFace triIs = tetIs.faceTriIs(mesh_);
214  const scalar v = tetIs.tet(mesh_).mag();
215 
216  cellValue[celli] += v*interpolate(tetCrds[0], tetIs);
217  cellGrad[celli] += v*interpolateGrad(tetCrds[0], tetIs);
218 
219  forAll(triIs, vertexI)
220  {
221  const label pointi = triIs[vertexI];
222 
223  pointVolume[pointi] += v;
224  pointValue[pointi] += v*interpolate(tetCrds[vertexI], tetIs);
225  pointGrad[pointi] += v*interpolateGrad(tetCrds[vertexI], tetIs);
226  }
227  }
228  }
229 
230  // average
231  cellValue.primitiveFieldRef() /= mesh_.V();
232  cellGrad.primitiveFieldRef() /= mesh_.V();
233  pointValue.primitiveFieldRef() /= pointVolume;
234  pointGrad.primitiveFieldRef() /= pointVolume;
235 
236  // write
237  if (!cellValue.write(valid)) return false;
238  if (!cellGrad.write(valid)) return false;
239  if (!pointValue.write(valid)) return false;
240  if (!pointGrad.write(valid)) return false;
241 
242  return true;
243 }
244 
245 
246 // ************************************************************************* //
Foam::AveragingMethod
Base class for lagrangian averaging methods.
Definition: KinematicParcel.H:69
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::tetIndices::tet
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet.
Definition: tetIndicesI.H:155
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::AveragingMethod::write
virtual bool write(const bool valid=true) const
Write using setting from DB.
Definition: AveragingMethod.C:141
Foam::FatalIOError
IOerror FatalIOError
AveragingMethod.H
Foam::AveragingMethod::New
static autoPtr< AveragingMethod< Type > > New(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Selector.
Definition: AveragingMethod.C:84
append
rAUs append(new volScalarField(IOobject::groupName("rAU", phase1.name()), 1.0/(U1Eqn.A()+byDt(max(phase1.residualAlpha() - alpha1, scalar(0)) *rho1))))
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::AveragingMethod::dict_
const dictionary & dict_
Protected data.
Definition: AveragingMethod.H:74
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::AveragingMethod::updateGrad
virtual void updateGrad()
Protected member functions.
Definition: AveragingMethod.C:36
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam::tetIndices::faceTriIs
triFace faceTriIs(const polyMesh &mesh, const bool warn=true) const
Return the indices corresponding to the tri on the face for.
Definition: tetIndicesI.H:68
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:51
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::AveragingMethod::AveragingMethod
AveragingMethod(const IOobject &io, const dictionary &dict, const fvMesh &mesh, const labelList &size)
Constructors.
Definition: AveragingMethod.C:44
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:83
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::regIOobject
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:73
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:69
runTimeSelectionTables.H
Macros to ease declaration of run-time selection tables.
Foam::List< label >
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::AveragingMethod::mesh_
const fvMesh & mesh_
The mesh on which the averaging is to be done.
Definition: AveragingMethod.H:77
Foam::AveragingMethod::writeData
virtual bool writeData(Ostream &) const
Dummy write.
Definition: AveragingMethod.C:134
Foam::barycentric
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition: barycentric.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::tetrahedron::mag
scalar mag() const
Return volume.
Definition: tetrahedronI.H:172
Foam::AveragingMethod::average
virtual void average()
Calculate the average.
Definition: AveragingMethod.C:115
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
pointMesh.H
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189