26 #include "BlendedInterfacialModel.H" 28 #include "surfaceInterpolate.H" 38 inline tmp<Foam::volScalarField>
46 inline tmp<Foam::surfaceScalarField>
58 template<
class ModelType>
59 template<
class GeoField>
65 typename GeoField::Boundary& fieldBf = field.boundaryFieldRef();
71 isA<fixedValueFvsPatchScalarField>
73 phase1_.phi()().boundaryField()[
patchi]
83 template<
class ModelType>
87 template<
class>
class PatchField,
94 tmp<GeometricField<Type, PatchField, GeoMesh>>
95 (ModelType::*method)(Args ...)
const,
97 const dimensionSet& dims,
102 typedef GeometricField<scalar, PatchField, GeoMesh> scalarGeoField;
103 typedef GeometricField<Type, PatchField, GeoMesh> typeGeoField;
105 tmp<scalarGeoField>
f1, f2;
107 if (model_.valid() || model1In2_.valid())
110 blendedInterfacialModel::interpolate<scalarGeoField>
112 blending_.f1(phase1_, phase2_)
116 if (model_.valid() || model2In1_.valid())
119 blendedInterfacialModel::interpolate<scalarGeoField>
121 blending_.f2(phase1_, phase2_)
131 ModelType::typeName +
":" + name,
132 phase1_.mesh().time().timeName(),
139 dimensioned<Type>(
"zero", dims,
Zero)
148 <<
"Cannot treat an interfacial model with no distinction " 149 <<
"between continuous and dispersed phases as signed" 153 x.ref() += (model_().*method)(args ...)*(scalar(1) -
f1() - f2());
156 if (model1In2_.valid())
158 x.ref() += (model1In2_().*method)(args ...)*
f1;
161 if (model2In1_.valid())
163 tmp<typeGeoField> dx = (model2In1_().*method)(args ...)*f2;
178 && (model_.valid() || model1In2_.valid() || model2In1_.valid())
181 correctFixedFluxBCs(
x.ref());
190 template<
class ModelType>
195 const blendingMethod& blending,
196 autoPtr<ModelType> model,
197 autoPtr<ModelType> model1In2,
198 autoPtr<ModelType> model2In1,
199 const bool correctFixedFluxBCs
206 IOobject::groupName(typeName, phasePair(phase1, phase2).name()),
215 model1In2_(model1In2),
216 model2In1_(model2In1),
217 correctFixedFluxBCs_(correctFixedFluxBCs)
221 template<
class ModelType>
225 const blendingMethod& blending,
226 const phasePair& pair,
227 const orderedPhasePair& pair1In2,
228 const orderedPhasePair& pair2In1,
229 const bool correctFixedFluxBCs
236 IOobject::groupName(typeName, pair.name()),
244 correctFixedFluxBCs_(correctFixedFluxBCs)
246 if (modelTable.found(pair))
258 if (modelTable.found(pair1In2))
264 modelTable[pair1In2],
270 if (modelTable.found(pair2In1))
276 modelTable[pair2In1],
286 template<
class ModelType>
293 template<
class ModelType>
296 const class phaseModel& phase
302 : model2In1_.valid();
306 template<
class ModelType>
309 const class phaseModel& phase
312 return &phase == &(phase1_) ? model1In2_ : model2In1_;
316 template<
class ModelType>
322 return evaluate(
k,
"K", ModelType::dimK,
false);
326 template<
class ModelType>
330 tmp<volScalarField> (ModelType::*
k)(
const scalar)
const = &
ModelType::K;
332 return evaluate(
k,
"K", ModelType::dimK,
false, residualAlpha);
336 template<
class ModelType>
340 return evaluate(&ModelType::Kf,
"Kf", ModelType::dimK,
false);
344 template<
class ModelType>
349 return evaluate(&
ModelType::F,
"F", ModelType::dimF,
true);
353 template<
class ModelType>
361 template<
class ModelType>
365 return evaluate(&ModelType::D,
"D", ModelType::dimD,
false);
369 template<
class ModelType>
373 return evaluate(&ModelType::dmdt,
"dmdt", ModelType::dimDmdt,
false);
377 template<
class ModelType>
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Dictionary hash table.
#define forAll(list, i)
Loop across all elements in list.
bool hasModel(const phaseModel &phase) const
Return true if a model is specified for the supplied phase.
errorManipArg< error, int > exit(error &err, const int errNo=1)
bool writeData(Ostream &os) const
Dummy write for regIOobject.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
~BlendedInterfacialModel()
Destructor.
tmp< volScalarField > K() const
Return the blended force coefficient.
tmp< surfaceScalarField > Ff() const
Return the face blended force.
tmp< surfaceScalarField > Kf() const
Return the face blended force coefficient.
tmp< volScalarField > dmdt() const
Return the blended mass transfer rate.
label k
Boltzmann constant.
tmp< GeometricField< Type, fvPatchField, volMesh > > F() const
Return the blended force.
CGAL::Exact_predicates_exact_constructions_kernel K
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
const ModelType & model(const phaseModel &phase) const
Return the model for the supplied phase.
static tmp< GeoField > interpolate(tmp< volScalarField > f)
Convenience function to interpolate blending values. Needs to be.
volVectorField F(fluid.F())
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.
tmp< volScalarField > D() const
Return the blended diffusivity.
A class for managing temporary objects.
surfaceScalarField Ff(fluid.Ff())
const dimensionSet dimArea(sqr(dimLength))