Go to the documentation of this file.
26 #include "twoPhaseSystem.H"
28 #include "BlendedInterfacialModel.H"
29 #include "virtualMassModel.H"
30 #include "heatTransferModel.H"
31 #include "liftModel.H"
32 #include "wallLubricationModel.H"
33 #include "turbulentDispersionModel.H"
47 #include "blendingMethod.H"
66 IOobject::MUST_READ_IF_MODIFIED,
107 IOobject::READ_IF_PRESENT,
114 phase2_.volScalarField::operator=(scalar(1) - phase1_);
120 blendingMethods_.insert
177 new BlendedInterfacialModel<dragModel>
181 blendingMethods_.found(
"drag")
182 ? blendingMethods_[
"drag"]
183 : blendingMethods_[
"default"]
194 new BlendedInterfacialModel<virtualMassModel>
198 blendingMethods_.found(
"virtualMass")
199 ? blendingMethods_[
"virtualMass"]
200 : blendingMethods_[
"default"]
210 new BlendedInterfacialModel<heatTransferModel>
214 blendingMethods_.found(
"heatTransfer")
215 ? blendingMethods_[
"heatTransfer"]
216 : blendingMethods_[
"default"]
226 new BlendedInterfacialModel<liftModel>
230 blendingMethods_.found(
"lift")
231 ? blendingMethods_[
"lift"]
232 : blendingMethods_[
"default"]
242 new BlendedInterfacialModel<wallLubricationModel>
244 lookup(
"wallLubrication"),
246 blendingMethods_.found(
"wallLubrication")
247 ? blendingMethods_[
"wallLubrication"]
248 : blendingMethods_[
"default"]
256 turbulentDispersion_.set
258 new BlendedInterfacialModel<turbulentDispersionModel>
260 lookup(
"turbulentDispersion"),
262 blendingMethods_.found(
"turbulentDispersion")
263 ? blendingMethods_[
"turbulentDispersion"]
264 : blendingMethods_[
"default"]
284 return phase1_*phase1_.thermo().rho() + phase2_*phase2_.thermo().rho();
290 return phase1_*phase1_.U() + phase2_*phase2_.U();
316 return virtualMass_->K();
322 return virtualMass_->Kf();
328 return heatTransfer_->K();
334 return lift_->F<
vector>() + wallLubrication_->F<
vector>();
340 return lift_->Ff() + wallLubrication_->Ff();
346 return turbulentDispersion_->D();
352 const Time& runTime = mesh_.time();
368 word alphaScheme(
"div(phi," +
alpha1.name() +
')');
371 alpha1.correctBoundaryConditions();
377 tmp<surfaceScalarField> alpha1alpha2f;
379 if (pPrimeByA_.valid())
422 if (dgdt_[celli] > 0.0)
424 Sp[celli] -= dgdt_[celli]/
max(1.0 -
alpha1[celli], 1
e-4);
425 Su[celli] += dgdt_[celli]/
max(1.0 -
alpha1[celli], 1
e-4);
427 else if (dgdt_[celli] < 0.0)
453 alphaPhic1.boundaryField()[
patchi];
455 if (!alphaPhic1p.coupled())
460 forAll(alphaPhic1p, facei)
462 if (phi1p[facei] < 0)
464 alphaPhic1p[facei] = alpha1p[facei]*phi1p[facei];
475 !(++alphaSubCycle).end();
486 (alphaSubCycle.index()*
Sp)(),
487 (
Su - (alphaSubCycle.index() - 1)*
Sp*
alpha1)(),
492 if (alphaSubCycle.index() == 1)
494 phase1_.alphaPhi() = alphaPhic10;
498 phase1_.alphaPhi() += alphaPhic10;
518 phase1_.alphaPhi() = alphaPhic1;
521 if (pPrimeByA_.valid())
532 phase1_.alphaPhi() += alpha1Eqn.flux();
535 phase1_.alphaRhoPhi() =
538 phase2_.alphaPhi() = phi_ - phase1_.alphaPhi();
540 phase2_.alphaRhoPhi() =
544 <<
alpha1.weightedAverage(mesh_.V()).value()
561 phase1_.turbulence().correct();
562 phase2_.turbulence().correct();
572 readOK &= phase1_.read(*
this);
573 readOK &= phase2_.read(*
this);
588 return drag_->phaseModel(phase);
595 return virtualMass_->phaseModel(phase);
601 return pair_->sigma();
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
surfaceScalarField & phi1
tmp< volScalarField > rho() const
Return the mixture density.
static autoPtr< blendingMethod > New(const dictionary &dict, const wordList &phaseNames)
tmp< volScalarField > Vm() const
Return the virtual mass coefficient.
fvsPatchField< scalar > fvsPatchScalarField
Calculate the snGrad of the given volField.
void correct()
Correct two-phase properties other than turbulence.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
virtual bool read()
Read object.
Calculate the divergence of the given field.
tmp< volScalarField > Kh() const
Return the heat transfer coefficient.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< volScalarField > sigma() const
Return the surface tension coefficient.
const dimensionedVector & g
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type & value() const
Return const reference to value.
surfaceScalarField phir(phic *interface.nHatf())
surfaceScalarField & phi2
tmp< volScalarField > Kd() const
Return the drag coefficient.
Constant dispersed-phase particle diameter model.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Calculate the curl of the given volField by constructing the Hodge-dual of the symmetric part of the ...
tmp< surfaceScalarField > calcPhi() const
Return the mixture flux.
List< word > wordList
A List of words.
const word dictName("particleTrackDict")
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
fvMatrix< scalar > fvScalarMatrix
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Calculate the matrix for the laplacian of the field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Generic dimensioned Type class.
twoPhaseSystem(const fvMesh &)
Construct from fvMesh.
HashTable< scalar, phasePairKey, phasePairKey::hash > scalarTable
Scalar hash table.
const double e
Elementary charge.
Vector< scalar > vector
A scalar version of the templated Vector.
tmp< surfaceScalarField > Vmf() const
Return the face virtual mass coefficient.
const virtualMassModel & virtualMass(const phaseModel &phase) const
Return the virtual mass model for the given phase.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
word alpharScheme("div(phirb,alpha)")
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
tmp< volScalarField > D() const
Return the turbulent diffusivity.
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Dictionary hash table.
bool read()
Read base phaseProperties dictionary.
Calculate the face-flux of the given field.
label readLabel(Istream &is)
void correctTurbulence()
Correct two-phase turbulence.
Calculate the first temporal derivative.
DimensionedField< Type, GeoMesh > DimensionedInternalField
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
tmp< volVectorField > F() const
Return the combined force (lift + wall-lubrication)
virtual void solve()
Solve for the phase fractions.
tmp< volVectorField > U() const
Return the mixture velocity.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
MULES: Multidimensional universal limiter for explicit solution.
const dragModel & drag(const phaseModel &phase) const
Return the drag model for the given phase.
virtual ~twoPhaseSystem()
Destructor.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Calulate the matrix for the first temporal derivative.
stressControl lookup("compactNormalStress") >> compactNormalStress
tmp< surfaceScalarField > Ff() const
Return the combined face-force (lift + wall-lubrication)