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-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 
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),
52  FieldField<Field, Type>(),
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),
73  FieldField<Field, Type>(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 
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 writeOnProc) 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::scopedName(this->name(), "cellValue"),
152  IOobject::NO_REGISTER,
153  mesh_,
155  );
156  auto& cellValue = tcellValue.ref();
157 
159  (
160  IOobject::scopedName(this->name(), "cellGrad"),
161  IOobject::NO_REGISTER,
162  mesh_,
164  );
165  auto& cellGrad = tcellGrad.ref();
166 
168  (
169  IOobject::scopedName(this->name(), "pointValue"),
170  IOobject::NO_REGISTER,
171  pointMesh_,
173  );
174  auto& pointValue = tpointValue.ref();
175 
177  (
178  IOobject::scopedName(this->name(), "pointGrad"),
179  IOobject::NO_REGISTER,
180  pointMesh_,
182  );
183  auto& pointGrad = tpointGrad.ref();
184 
185  // Barycentric coordinates of the tet vertices
187  tetCrds
188  ({
189  barycentric(1, 0, 0, 0),
190  barycentric(0, 1, 0, 0),
191  barycentric(0, 0, 1, 0),
192  barycentric(0, 0, 0, 1)
193  });
194 
195  // tet-volume weighted sums
196  forAll(mesh_.C(), celli)
197  {
198  const List<tetIndices> cellTets =
199  polyMeshTetDecomposition::cellTetIndices(mesh_, celli);
200 
201  forAll(cellTets, tetI)
202  {
203  const tetIndices& tetIs = cellTets[tetI];
204  const triFace triIs = tetIs.faceTriIs(mesh_);
205  const scalar v = tetIs.tet(mesh_).mag();
206 
207  cellValue[celli] += v*interpolate(tetCrds[0], tetIs);
208  cellGrad[celli] += v*interpolateGrad(tetCrds[0], tetIs);
209 
210  forAll(triIs, vertexI)
211  {
212  const label pointi = triIs[vertexI];
213 
214  pointVolume[pointi] += v;
215  pointValue[pointi] += v*interpolate(tetCrds[vertexI], tetIs);
216  pointGrad[pointi] += v*interpolateGrad(tetCrds[vertexI], tetIs);
217  }
218  }
219  }
220 
221  // average
222  cellValue.primitiveFieldRef() /= mesh_.V();
223  cellGrad.primitiveFieldRef() /= mesh_.V();
224  pointValue.primitiveFieldRef() /= pointVolume;
225  pointGrad.primitiveFieldRef() /= pointVolume;
226 
227  // write
228  if (!cellValue.write(writeOnProc)) return false;
229  if (!cellGrad.write(writeOnProc)) return false;
230  if (!pointValue.write(writeOnProc)) return false;
231  if (!pointGrad.write(writeOnProc)) return false;
232 
233  return true;
234 }
235 
236 
237 // ************************************************************************* //
dictionary dict
static autoPtr< AveragingMethod< Type > > New(const IOobject &io, const dictionary &dict, const fvMesh &mesh)
Selector.
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:107
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
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.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Generic dimensioned Type class.
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition: barycentric.H:41
const dimensionSet dimless
Dimensionless.
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:45
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Generic templated field type.
Definition: Field.H:62
A class for handling words, derived from Foam::string.
Definition: word.H:63
face triFace(3)
errorManip< error > abort(error &err)
Definition: errorManip.H:139
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
virtual void updateGrad()
Protected member functions.
virtual bool writeData(Ostream &) const
Dummy write.
rAUs append(new volScalarField(IOobject::groupName("rAU", phase1.name()), 1.0/(U1Eqn.A()+byDt(max(phase1.residualAlpha() - alpha1, scalar(0)) *rho1))))
virtual void average()
Calculate the average.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:68
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Macros to ease declaration of run-time selection tables.
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:645
AveragingMethod(const IOobject &io, const dictionary &dict, const fvMesh &mesh, const labelList &size)
Constructors.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127