CloudAverageField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2025 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 "CloudAverageField.H"
27 #include "LagrangianFields.H"
28 #include "LagrangianSubFields.H"
29 #include "volMesh.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class Type>
34 template<class MethodWeightSum, class MethodNoWeightSum, class ... Args>
36 (
37  MethodWeightSum mws,
38  MethodNoWeightSum mnws,
39  const LagrangianSubMesh& subMesh,
40  const Args& ... args
41 ) const
42 {
43  tcellWeightSum_.valid()
44  ? (
45  notNull(weightPsiOrPsiState_)
46  ? (psiAveragePtr_().*mws)
47  (
48  subMesh.sub(weightPsiOrPsiState_),
49  args ...
50  )
51  // isNull(weightPsiOrPsiState_)
52  : (psiAveragePtr_().*mws)
53  (
54  toSubField(weightPsiOrPsiDerived_(subMesh)),
55  args ...
56  )
57  )
58  : (
59  notNull(weightState_) && notNull(weightPsiOrPsiState_)
60  ? (psiAveragePtr_().*mnws)
61  (
62  subMesh.sub(weightState_),
63  subMesh.sub(weightPsiOrPsiState_),
64  args ...
65  )
66  : notNull(weightState_) && isNull(weightPsiOrPsiState_)
67  ? (psiAveragePtr_().*mnws)
68  (
69  subMesh.sub(weightState_),
70  toSubField(weightPsiOrPsiDerived_(subMesh)),
71  args ...
72  )
73  : isNull(weightState_) && notNull(weightPsiOrPsiState_)
74  ? (psiAveragePtr_().*mnws)
75  (
76  toSubField(weightDerived_(subMesh)),
77  subMesh.sub(weightPsiOrPsiState_),
78  args ...
79  )
80  // isNull(weightState_) && isNull(weightPsiOrPsiState_)
81  : (psiAveragePtr_().*mnws)
82  (
83  toSubField(weightDerived_(subMesh)),
84  toSubField(weightPsiOrPsiDerived_(subMesh)),
85  args ...
86  )
87  );
88 }
89 
90 
91 template<class Type>
94 (
95  const LagrangianModelRef&,
96  const LagrangianSubMesh& subMesh
97 ) const
98 {
99  const LagrangianMesh& mesh = subMesh.mesh();
100 
101  if (psiAverageState_ == psiAverageState::removed)
102  {
104  << "Cannot interpolate an average that has removed elements"
105  << exit(FatalError);
106  }
107 
108  if (!psiAveragePtr_.valid())
109  {
110  psiAveragePtr_ =
111  tcellWeightSum_.valid()
112  ? (
113  notNull(weightPsiOrPsiState_)
115  (
116  word(mesh.schemes().averaging(name_)),
117  "average(" + name_ + ')',
118  tcellWeightSum_(),
119  weightPsiOrPsiState_(mesh)
120  )
121  // isNull(weightPsiOrPsiState_)
122  : LagrangianAverage<Type>::New
123  (
124  word(mesh.schemes().averaging(name_)),
125  "average(" + name_ + ')',
126  tcellWeightSum_(),
127  weightPsiOrPsiDerived_(mesh)()
128  )
129  )
130  : (
131  notNull(weightState_) && notNull(weightPsiOrPsiState_)
133  (
134  word(mesh.schemes().averaging(name_)),
135  "average(" + name_ + ')',
136  weightState_(mesh),
137  weightPsiOrPsiState_(mesh)
138  )
139  : notNull(weightState_) && isNull(weightPsiOrPsiState_)
140  ? LagrangianAverage<Type>::New
141  (
142  word(mesh.schemes().averaging(name_)),
143  "average(" + name_ + ')',
144  weightState_(mesh),
145  weightPsiOrPsiDerived_(mesh)()
146  )
147  : isNull(weightState_) && notNull(weightPsiOrPsiState_)
148  ? LagrangianAverage<Type>::New
149  (
150  word(mesh.schemes().averaging(name_)),
151  "average(" + name_ + ')',
152  weightDerived_(mesh)(),
153  weightPsiOrPsiState_(mesh)
154  )
155  // isNull(weightState_) && isNull(weightPsiOrPsiState_)
156  : LagrangianAverage<Type>::New
157  (
158  word(mesh.schemes().averaging(name_)),
159  "average(" + name_ + ')',
160  weightDerived_(mesh)(),
161  weightPsiOrPsiDerived_(mesh)()
162  )
163  );
164 
165  // Pointers to average manipulation methods
166  void (LagrangianAverage<Type>::*removeWeightSum)
167  (
170  void (LagrangianAverage<Type>::*removeNoWeightSum)
171  (
175  void (LagrangianAverage<Type>::*addWeightSum)
176  (
178  const bool
180  void (LagrangianAverage<Type>::*addNoWeightSum)
181  (
184  const bool
186 
187  switch (psiAverageState_)
188  {
189  case psiAverageState::removed:
190  op(removeWeightSum, removeNoWeightSum, subMesh);
191  break;
192  case psiAverageState::complete:
193  break;
194  case psiAverageState::cached:
195  op(removeWeightSum, removeNoWeightSum, subMesh);
196  op(addWeightSum, addNoWeightSum, subMesh, true);
197  break;
198  }
199  }
200 
201  return psiAveragePtr_->interpolate(subMesh);
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
206 
207 template<class Type>
209 (
210  const word& name,
211  const DimensionedField<scalar, volMesh>& cellWeightSum,
212  const CloudDerivedField<Type>& weightPsi
213 )
214 :
216  name_(name),
217  tcellWeightSum_(cellWeightSum),
218  weightState_(NullObjectRef<CloudStateField<scalar>>()),
219  weightDerived_(NullObjectRef<CloudDerivedField<scalar>>()),
220  weightPsiOrPsiState_(NullObjectRef<CloudStateField<Type>>()),
221  weightPsiOrPsiDerived_(weightPsi),
222  psiAveragePtr_(nullptr),
223  psiAverageState_(psiAverageState::complete)
224 {}
225 
226 
227 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
228 
229 template<class Type>
231 {
232  psiAverageState_ = psiAverageState::removed;
233 
234  if (!psiAveragePtr_.valid()) return;
235 
237 
238  // Pointers to average manipulation methods
239  void (LagrangianAverage<Type>::*removeWeightSum)
240  (
243  void (LagrangianAverage<Type>::*removeNoWeightSum)
244  (
248 
249  op(removeWeightSum, removeNoWeightSum, subMesh);
250 }
251 
252 
253 template<class Type>
255 (
256  const LagrangianSubMesh& subMesh,
257  const bool cache
258 )
259 {
260  psiAverageState_ =
261  cache ? psiAverageState::cached : psiAverageState::complete;
262 
263  if (!psiAveragePtr_.valid()) return;
264 
266 
267  // Pointers to average manipulation methods
268  void (LagrangianAverage<Type>::*addWeightSum)
269  (
271  const bool
273  void (LagrangianAverage<Type>::*addNoWeightSum)
274  (
277  const bool
279 
280  op(addWeightSum, addNoWeightSum, subMesh, cache);
281 }
282 
283 
284 template<class Type>
286 (
287  const LagrangianSubMesh& subMesh,
288  const bool cache
289 )
290 {
291  psiAverageState_ =
292  cache ? psiAverageState::cached : psiAverageState::complete;
293 
294  if (!psiAveragePtr_.valid()) return;
295 
297 
298  // Pointers to average manipulation methods
299  void (LagrangianAverage<Type>::*correctWeightSum)
300  (
302  const bool
304  void (LagrangianAverage<Type>::*correctNoWeightSum)
305  (
308  const bool
310 
311  op(correctWeightSum, correctNoWeightSum, subMesh, cache);
312 }
313 
314 
315 template<class Type>
317 {
319 
320  psiAveragePtr_.clear();
321  psiAverageState_ = psiAverageState::complete;
322 }
323 
324 
325 // ************************************************************************* //
A cloud field that is averaged and then interpolated back to the cloud. Uses CloudDerivedField to pro...
void remove(const LagrangianSubMesh &subMesh)
Remove this sub-mesh from the average.
CloudAverageField(const word &name, const DimensionedField< scalar, volMesh > &cellWeightSum, const CloudDerivedField< Type > &weightPsi)
Construct from a name, a cell weight sum and a derived field.
void add(const LagrangianSubMesh &subMesh, const bool cache)
Add this sub-mesh to the average.
void correct(const LagrangianSubMesh &subMesh, const bool cache)
Correct the average with values from the sub-mesh.
A field derived from other state fields of the cloud. Stores and virtualises a function or a method w...
void clear(const bool final)
Clear.
A field which is stored as part of the state of the cloud. This is a light wrapper around a dynamic L...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Basic Lagrangian averaging process in which values are averaged in the cells and then interpolated as...
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1799
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:443
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:64
const T & NullObjectRef()
Return const reference to the nullObject of type T.
Definition: nullObjectI.H:27
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:58
error FatalError
tmp< DimensionedField< Type, GeoMesh, SubField > > toSubField(const DimensionedField< Type, GeoMesh, Field > &)
Return a temporary sub-field from a reference to a field.
Foam::argList args(argc, argv)