31#include "surfaceInterpolate.H"
61template<
class ModelType>
62template<
class GeoField>
63void Foam::BlendedInterfacialModel<ModelType>::correctFixedFluxBCs
68 typename GeoField::Boundary& fieldBf =
field.boundaryFieldRef();
70 forAll(phase1_.phi()().boundaryField(), patchi)
74 isA<fixedValueFvsPatchScalarField>
76 phase1_.phi()().boundaryField()[patchi]
80 fieldBf[patchi] = Zero;
86template<
class ModelType>
90 template<
class>
class PatchField,
95Foam::BlendedInterfacialModel<ModelType>::evaluate
98 (ModelType::*method)(Args ...)
const,
110 if (model_ || model1In2_)
115 blending_.f1(phase1_, phase2_)
119 if (model_ || model2In1_)
124 blending_.f2(phase1_, phase2_)
128 auto x = typeGeoField::New
141 <<
"Cannot treat an interfacial model with no distinction "
142 <<
"between continuous and dispersed phases as signed"
146 x.ref() += (model_().*method)(
args ...)*(scalar(1) - f1() - f2());
151 x.ref() += (model1In2_().*method)(
args ...)*f1;
171 && (model_ || model1In2_ || model2In1_)
174 correctFixedFluxBCs(
x.ref());
183template<
class ModelType>
184Foam::BlendedInterfacialModel<ModelType>::BlendedInterfacialModel
192 const bool correctFixedFluxBCs
208 model1In2_(model1In2),
209 model2In1_(model2In1),
210 correctFixedFluxBCs_(correctFixedFluxBCs)
214template<
class ModelType>
215Foam::BlendedInterfacialModel<ModelType>::BlendedInterfacialModel
222 const bool correctFixedFluxBCs
237 correctFixedFluxBCs_(correctFixedFluxBCs)
251 if (modelTable.
found(pair1In2))
257 modelTable[pair1In2],
263 if (modelTable.
found(pair2In1))
269 modelTable[pair2In1],
279template<
class ModelType>
286template<
class ModelType>
289 const class phaseModel& phase
299template<
class ModelType>
305 return &
phase == &(phase1_) ? model1In2_ : model2In1_;
309template<
class ModelType>
315 return evaluate(
k,
"K", ModelType::dimK,
false);
319template<
class ModelType>
325 return evaluate(
k,
"K", ModelType::dimK,
false, residualAlpha);
329template<
class ModelType>
333 return evaluate(&ModelType::Kf,
"Kf", ModelType::dimK,
false);
337template<
class ModelType>
342 return evaluate(&ModelType::F,
"F", ModelType::dimF,
true);
346template<
class ModelType>
350 return evaluate(&ModelType::Ff,
"Ff", ModelType::dimF*
dimArea,
true);
354template<
class ModelType>
358 return evaluate(&ModelType::D,
"D", ModelType::dimD,
false);
362template<
class ModelType>
366 return evaluate(&ModelType::dmdt,
"dmdt", ModelType::dimDmdt,
false);
370template<
class ModelType>
tmp< surfaceScalarField > Kf() const
Return the face blended force coefficient.
tmp< surfaceScalarField > Ff() const
Return the face blended force.
tmp< volScalarField > D() const
Return the blended diffusivity.
bool hasModel(const phaseModel &phase) const
Return true if a model is specified for the supplied phase.
tmp< GeometricField< Type, fvPatchField, volMesh > > F() const
Return the blended force.
bool writeData(Ostream &os) const
Dummy write for regIOobject.
Foam::tmp< Foam::GeometricField< Type, PatchField, GeoMesh > > evaluate(tmp< GeometricField< Type, PatchField, GeoMesh > >(Foam::dragModel::*method)(Args ...) const, const word &name, const dimensionSet &dims, const bool subtract, Args ... args) const
tmp< volScalarField > K() const
Return the blended force coefficient.
const ModelType & model(const phaseModel &phase) const
Return the model for the supplied phase.
const modelType & phaseModel(const phaseModel &phase) const
Return the model for the supplied phase.
~BlendedInterfacialModel()
Destructor.
tmp< volScalarField > dmdt() const
Return the blended mass transfer rate.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Generic GeometricField class.
bool found(const Key &key) const
Same as contains().
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const Time & time() const noexcept
Return Time associated with the objectRegistry.
const word & name() const noexcept
Return the object name.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
static tmp< GeoField > interpolate(tmp< volScalarField > f)
Convenience function to interpolate blending values. Needs to be.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Generic dimensioned Type class.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Dictionary hash table.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
regIOobject(const IOobject &io, const bool isTimeObject=false)
Construct from IOobject. The optional flag adds special handling if the object is the top-level regIO...
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
void subtract(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
const dimensionSet dimArea(sqr(dimLength))
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.