36#ifndef twoPhaseSystem_H
37#define twoPhaseSystem_H
40#include "phaseModel.H"
42#include "orderedPhasePair.H"
52class virtualMassModel;
53class heatTransferModel;
55class wallLubricationModel;
56class turbulentDispersionModel;
59template<
class modelType>
class BlendedInterfacialModel;
87 tmp<surfaceScalarField> pPrimeByA_;
90 autoPtr<phasePair> pair_;
93 autoPtr<orderedPhasePair> pair1In2_;
96 autoPtr<orderedPhasePair> pair2In1_;
99 HashTable<autoPtr<blendingMethod>> blendingMethods_;
102 autoPtr<BlendedInterfacialModel<dragModel>> drag_;
105 autoPtr<BlendedInterfacialModel<virtualMassModel>> virtualMass_;
108 autoPtr<BlendedInterfacialModel<heatTransferModel>> heatTransfer_;
111 autoPtr<BlendedInterfacialModel<liftModel>> lift_;
114 autoPtr<BlendedInterfacialModel<wallLubricationModel>>
118 autoPtr<BlendedInterfacialModel<turbulentDispersionModel>>
119 turbulentDispersion_;
125 tmp<surfaceScalarField> calcPhi()
const;
188 template<
class modelType>
192 template<
class modelType>
const uniformDimensionedVectorField & g
Mesh data needed to do the Finite Volume discretisation.
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...
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
A class for managing temporary objects.
tmp< surfaceScalarField > Ff() const
Return the combined face-force (lift + wall-lubrication)
tmp< volVectorField > F() const
Return the combined force (lift + wall-lubrication)
const phaseModel & phase2() const
Return phase model 2.
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument.
tmp< volScalarField > D() const
Return the turbulent diffusivity.
virtual ~twoPhaseSystem()
Destructor.
const surfaceScalarField & phi() const
Return the mixture flux.
tmp< volScalarField > Kd() const
Return the drag coefficient.
void correct()
Correct two-phase properties other than turbulence.
virtual ~twoPhaseSystem()
Destructor.
phaseModel & phase1_
Phase model 1.
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
tmp< surfaceScalarField > Vmf() const
Return the face virtual mass coefficient.
tmp< volScalarField > Vm() const
Return the virtual mass coefficient.
tmp< volVectorField > U() const
Return the mixture velocity.
const volScalarField & dgdt() const
Return the dilatation term.
const modelType & lookupSubModel(const phasePair &key) const
Access a sub model between a phase pair.
phaseModel & phase2_
Phase model 2.
tmp< surfaceScalarField > & pPrimeByA()
Return non-const access to the dispersion diffusivity.
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
tmp< volScalarField > Kd() const
Return the drag coefficient.
void correctTurbulence()
Correct two-phase turbulence.
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument.
phaseModel & phase1()
Return non-const access to phase model 1.
phaseModel & phase2()
Return non-const access to phase model 2.
const dimensionedScalar & sigma() const
Return the surface tension coefficient.
const phaseModel & phase1() const
Return phase model 1.
const fvMesh & mesh() const
Return the mesh.
tmp< volScalarField > sigma() const
Return the surface tension coefficient.
tmp< volScalarField > rho() const
Return the mixture density.
const modelType & lookupSubModel(const phaseModel &dispersed, const phaseModel &continuous) const
Access a sub model between two phases.
twoPhaseSystem(const fvMesh &)
Construct from fvMesh.
const phaseModel & phase1() const
Return phase model 1.
const phaseModel & phase2() const
Return phase model 2.
tmp< volScalarField > Kh() const
Return the heat transfer coefficient.
void solve()
Solve for the two-phase-fractions.
bool read()
Read base phaseProperties dictionary.
tmp< volScalarField > Vm() const
Return the virtual mass coefficient.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.