Go to the documentation of this file.
45 #ifndef multiphaseSystem_H
46 #define multiphaseSystem_H
50 #include "phaseModel.H"
54 #include "dragModel.H"
82 public Hash<interfacePair>
97 public Hash<interfacePair>
129 friend bool operator==
137 ((a.first() ==
b.first()) && (a.
second() ==
b.second()))
138 || ((a.first() ==
b.second()) && (a.
second() ==
b.first()))
142 friend bool operator!=
212 void correctContactAngle
216 surfaceVectorField::Boundary& nHatb
const T & second() const noexcept
Return second element, which is also the last element.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
virtual ~multiphaseSystem()=default
Destructor.
A class for handling words, derived from Foam::string.
void solve()
Solve for the mixture phase-fractions.
tmp< volScalarField > rho() const
Return the mixture density.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
A class for managing temporary objects.
const volScalarField & alpha2
Template dictionary class which manages the storage associated with it.
Incompressible multi-phase mixture with built in solution for the phase fractions with interface comp...
interfacePair(const word &alpha1Name, const word &alpha2Name)
HashPtrTable< dragModel, interfacePair, interfacePair::symmHash > dragModelTable
const volScalarField & alpha1
Dimension set for the base types.
void correct()
Dummy correct.
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
const word & name() const
interfacePair(const phaseModel &alpha1, const phaseModel &alpha2)
HashPtrTable< volScalarField, interfacePair, interfacePair::symmHash > dragCoeffFields
CGAL::Exact_predicates_exact_constructions_kernel K
Hash function class. The default definition is for primitives, non-primitives used to hash entries on...
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
label operator()(const interfacePair &key) const
Hashing function for string and derived string classes.
autoPtr< dragCoeffFields > dragCoeffs() const
Return the drag coefficients for all of the interfaces.
bool read()
Read base transportProperties dictionary.
Mesh data needed to do the Finite Volume discretisation.
const PtrDictionary< phaseModel > & phases() const
Return the phases.
Base-class for all transport models used by the incompressible turbulence models.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
An ordered pair of two objects of type <T> with first() and second() elements.
PtrDictionary< phaseModel > & phases()
Return the phases.
label operator()(const interfacePair &key) const
tmp< volScalarField > dragCoeff(const phaseModel &phase, const dragCoeffFields &dragCoeffs) const
Return the sum of the drag coefficients for the given phase.
tmp< volVectorField > Svm(const phaseModel &phase) const
Return the virtual-mass source for the given phase.
const dragModelTable & dragModels() const
Return the table of drag models.
tmp< volScalarField > nu() const
Return the mixture laminar viscosity.
tmp< volScalarField > Cvm(const phaseModel &phase) const
Return the virtual-mass coefficient for the given phase.
multiphaseSystem(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.