Go to the documentation of this file.
34 #ifndef twoPhaseSystem_H
35 #define twoPhaseSystem_H
38 #include "phaseModel.H"
39 #include "phasePair.H"
40 #include "orderedPhasePair.H"
43 #include "dragModel.H"
50 class virtualMassModel;
51 class heatTransferModel;
53 class wallLubricationModel;
54 class turbulentDispersionModel;
57 template <
class modelType>
class BlendedInterfacialModel;
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
tmp< surfaceScalarField > & pPrimeByA()
Return non-const access to the dispersion diffusivity.
const fvMesh & mesh() const
Return the mesh.
autoPtr< orderedPhasePair > pair1In2_
Phase pair for phase 1 dispersed in phase 2.
A class for handling words, derived from string.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
A class for managing temporary objects.
tmp< volScalarField > rho() const
Return the mixture density.
tmp< volScalarField > Vm() const
Return the virtual mass coefficient.
autoPtr< BlendedInterfacialModel< virtualMassModel > > virtualMass_
Virtual mass model.
void correct()
Correct two-phase properties other than turbulence.
surfaceScalarField phi_
Total volumetric flux.
autoPtr< BlendedInterfacialModel< heatTransferModel > > heatTransfer_
Heat transfer model.
tmp< volScalarField > Kh() const
Return the heat transfer coefficient.
tmp< volScalarField > sigma() const
Return the surface tension coefficient.
const dimensionedVector & g
autoPtr< BlendedInterfacialModel< dragModel > > drag_
Drag model.
tmp< volScalarField > Kd() const
Return the drag coefficient.
phaseModel phase1_
Phase model 1.
const phaseModel & phase2() const
Constant access phase model 2.
autoPtr< orderedPhasePair > pair2In1_
Phase pair for phase 2 dispersed in phase 1.
autoPtr< phasePair > pair_
Unordered phase pair.
tmp< surfaceScalarField > calcPhi() const
Return the mixture flux.
volScalarField dgdt_
Dilatation term.
const phaseModel & otherPhase(const phaseModel &phase) const
Constant access the phase not given as an argument.
tmp< surfaceScalarField > pPrimeByA_
Optional dispersion diffusivity.
const phaseModel & phase1() const
Constant access phase model 1.
autoPtr< BlendedInterfacialModel< wallLubricationModel > > wallLubrication_
Wall lubrication model.
const surfaceScalarField & phi() const
Return the mixture flux.
Generic dimensioned Type class.
twoPhaseSystem(const fvMesh &)
Construct from fvMesh.
autoPtr< BlendedInterfacialModel< liftModel > > lift_
Lift model.
Mesh data needed to do the Finite Volume discretisation.
phaseModel phase2_
Phase model 2.
phaseModel & phase2_
Phase model 2.
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.
An STL-conforming hash table.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Hashing function class, shared by all the derived classes.
autoPtr< BlendedInterfacialModel< turbulentDispersionModel > > turbulentDispersion_
Wall lubrication model.
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
tmp< volScalarField > D() const
Return the turbulent diffusivity.
bool read()
Read base phaseProperties dictionary.
const volScalarField & dgdt() const
Return the dilatation term.
phaseModel & phase1_
Phase model 1.
void correctTurbulence()
Correct two-phase turbulence.
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.
Generic GeometricField class.
const fvMesh & mesh_
Reference to the mesh.
const dragModel & drag(const phaseModel &phase) const
Return the drag model for the given phase.
virtual ~twoPhaseSystem()
Destructor.
HashTable< autoPtr< blendingMethod >, word, word::hash > blendingMethods_
Blending methods.
tmp< surfaceScalarField > Ff() const
Return the combined face-force (lift + wall-lubrication)