34 template<
class MethodWeightSum,
class MethodNoWeightSum,
class ... Args>
38 MethodNoWeightSum mnws,
39 const LagrangianSubMesh& subMesh,
43 tcellWeightSum_.valid()
46 ? (psiAveragePtr_().*mws)
48 subMesh.sub(weightPsiOrPsiState_),
52 : (psiAveragePtr_().*mws)
60 ? (psiAveragePtr_().*mnws)
62 subMesh.sub(weightState_),
63 subMesh.sub(weightPsiOrPsiState_),
67 ? (psiAveragePtr_().*mnws)
69 subMesh.sub(weightState_),
74 ? (psiAveragePtr_().*mnws)
77 subMesh.sub(weightPsiOrPsiState_),
81 : (psiAveragePtr_().*mnws)
95 const LagrangianModelRef&,
96 const LagrangianSubMesh& subMesh
99 const LagrangianMesh&
mesh = subMesh.
mesh();
101 if (psiAverageState_ == psiAverageState::removed)
104 <<
"Cannot interpolate an average that has removed elements"
108 if (!psiAveragePtr_.valid())
111 tcellWeightSum_.valid()
117 "average(" + name_ +
')',
119 weightPsiOrPsiState_(
mesh)
124 word(
mesh.schemes().averaging(name_)),
125 "average(" + name_ +
')',
127 weightPsiOrPsiDerived_(
mesh)()
134 word(
mesh.schemes().averaging(name_)),
135 "average(" + name_ +
')',
137 weightPsiOrPsiState_(
mesh)
142 word(
mesh.schemes().averaging(name_)),
143 "average(" + name_ +
')',
145 weightPsiOrPsiDerived_(
mesh)()
150 word(
mesh.schemes().averaging(name_)),
151 "average(" + name_ +
')',
152 weightDerived_(
mesh)(),
153 weightPsiOrPsiState_(
mesh)
158 word(
mesh.schemes().averaging(name_)),
159 "average(" + name_ +
')',
160 weightDerived_(
mesh)(),
161 weightPsiOrPsiDerived_(
mesh)()
187 switch (psiAverageState_)
189 case psiAverageState::removed:
190 op(removeWeightSum, removeNoWeightSum, subMesh);
192 case psiAverageState::complete:
194 case psiAverageState::cached:
195 op(removeWeightSum, removeNoWeightSum, subMesh);
196 op(addWeightSum, addNoWeightSum, subMesh,
true);
201 return psiAveragePtr_->interpolate(subMesh);
217 tcellWeightSum_(cellWeightSum),
221 weightPsiOrPsiDerived_(weightPsi),
222 psiAveragePtr_(nullptr),
223 psiAverageState_(psiAverageState::
complete)
232 psiAverageState_ = psiAverageState::removed;
234 if (!psiAveragePtr_.valid())
return;
249 op(removeWeightSum, removeNoWeightSum, subMesh);
261 cache ? psiAverageState::cached : psiAverageState::complete;
263 if (!psiAveragePtr_.valid())
return;
280 op(addWeightSum, addNoWeightSum, subMesh, cache);
292 cache ? psiAverageState::cached : psiAverageState::complete;
294 if (!psiAveragePtr_.valid())
return;
311 op(correctWeightSum, correctNoWeightSum, subMesh, cache);
320 psiAveragePtr_.clear();
321 psiAverageState_ = psiAverageState::complete;
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.
const polyMesh & mesh() const
Return reference to polyMesh.
A class for managing temporary objects.
A class for handling words, derived from string.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
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)
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
const T & NullObjectRef()
Return const reference to the nullObject of type T.
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.
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)