31#include "surfaceInterpolate.H"
35template<
class modelType>
36template<
class GeometricField>
37void Foam::BlendedInterfacialModel<modelType>::correctFixedFluxBCs
42 auto& fieldBf =
field.boundaryFieldRef();
44 forAll(pair_.phase1().phi().boundaryField(), patchi)
48 isA<fixedValueFvsPatchScalarField>
50 pair_.phase1().phi().boundaryField()[patchi]
54 fieldBf[patchi] = Zero;
62template<
class modelType>
63Foam::BlendedInterfacialModel<modelType>::BlendedInterfacialModel
65 const phasePair::dictTable& modelTable,
66 const blendingMethod& blending,
67 const phasePair& pair,
68 const orderedPhasePair& pair1In2,
69 const orderedPhasePair& pair2In1,
70 const bool correctFixedFluxBCs
77 correctFixedFluxBCs_(correctFixedFluxBCs)
79 if (modelTable.found(pair_))
91 if (modelTable.found(pair1In2_))
97 modelTable[pair1In2_],
103 if (modelTable.found(pair2In1_))
109 modelTable[pair2In1_],
119template<
class modelType>
126template<
class modelType>
127Foam::tmp<Foam::volScalarField>
132 if (model_ || model1In2_)
134 f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
137 if (model_ || model2In1_)
139 f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
146 pair_.phase1().mesh(),
152 x.ref() += model_->K()*(f1() - f2());
157 x.ref() += model1In2_->K()*(1 - f1);
162 x.ref() += model2In1_->K()*f2;
168 && (model_ || model1In2_ || model2In1_)
171 correctFixedFluxBCs(
x.ref());
178template<
class modelType>
179Foam::tmp<Foam::surfaceScalarField>
184 if (model_ || model1In2_)
188 blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed())
192 if (model_ || model2In1_)
196 blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed())
203 pair_.phase1().mesh(),
209 x.ref() += model_->Kf()*(f1() - f2());
214 x.ref() += model1In2_->Kf()*(1 - f1);
219 x.ref() += model2In1_->Kf()*f2;
225 && (model_ || model1In2_ || model2In1_)
228 correctFixedFluxBCs(
x.ref());
235template<
class modelType>
237Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh>>
242 if (model_ || model1In2_)
244 f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
247 if (model_ || model2In1_)
249 f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
256 pair_.phase1().mesh(),
262 x.ref() += model_->F()*(f1() - f2());
267 x.ref() += model1In2_->F()*(1 - f1);
272 x.ref() -= model2In1_->F()*f2;
278 && (model_ || model1In2_ || model2In1_)
281 correctFixedFluxBCs(
x.ref());
288template<
class modelType>
289Foam::tmp<Foam::surfaceScalarField>
294 if (model_ || model1In2_)
298 blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed())
302 if (model_ || model2In1_)
306 blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed())
314 pair_.phase1().mesh(),
318 x.ref().setOriented();
322 x.ref() += model_->Ff()*(f1() - f2());
327 x.ref() += model1In2_->Ff()*(1 - f1);
332 x.ref() -= model2In1_->Ff()*f2;
338 && (model_ || model1In2_ || model2In1_)
341 correctFixedFluxBCs(
x.ref());
348template<
class modelType>
349Foam::tmp<Foam::volScalarField>
354 if (model_ || model1In2_)
356 f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
359 if (model_ || model2In1_)
361 f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
368 pair_.phase1().mesh(),
374 x.ref() += model_->D()*(f1() - f2());
379 x.ref() += model1In2_->D()*(1 - f1);
384 x.ref() += model2In1_->D()*f2;
390 && (model_ || model1In2_ || model2In1_)
393 correctFixedFluxBCs(
x.ref());
400template<
class modelType>
408 &
phase == &(pair_.phase1())
415template<
class modelType>
421 return &
phase == &(pair_.phase1()) ? *model1In2_ : *model2In1_;
if(maxValue - minValue< SMALL)
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.
tmp< volScalarField > K() const
Return the blended force coefficient.
const modelType & phaseModel(const phaseModel &phase) const
Return the model for the supplied phase.
~BlendedInterfacialModel()
Destructor.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Generic dimensioned Type class.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
A class for managing temporary objects.
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.
const dimensionSet dimArea(sqr(dimLength))
static constexpr const zero Zero
Global zero (0).
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
#define forAll(list, i)
Loop across all elements in list.