36#ifndef multiphaseInterSystem_H
37#define multiphaseInterSystem_H
41#include "phaseModel.H"
43#include "orderedPhasePair.H"
59namespace multiphaseInter
61class surfaceTensionModel;
188 template<
class modelType>
201 template<
class modelType>
215 template<
class modelType>
218 const word& modelName,
229 template<
class modelType>
232 const word& modelName,
244 template<
class modelType>
247 const word& modelName,
601 const word speciesName
611 virtual void solve() = 0;
649 template <
class modelType>
653 template <
class modelType>
CGAL::Exact_predicates_exact_constructions_kernel K
compressible::turbulenceModel & turb
const volScalarField & alpha1
const volScalarField & alpha2
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
A HashTable similar to std::unordered_map.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Abstract base-class for fluid and solid thermodynamic properties.
virtual const volScalarField & T() const
Temperature [K].
virtual volScalarField & p()
Pressure [Pa].
Base-class for all transport models used by the compressible turbulence models.
Abstract base class for turbulence models (RAS, LES and laminar).
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Mesh data needed to do the Finite Volume discretisation.
virtual tmp< fvScalarMatrix > volTransfer(const volScalarField &p)=0
Return the volumetric rate transfer matrix.
void createSubModels(const dictTable &modelDicts, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables using pair keys.
void generatePairsTable()
Generate pair table.
tmp< surfaceScalarField > surfaceTensionForce() const
Calculate surface tension of the mixture.
tmp< volScalarField > nearInterface() const
Near Interface of alpha'n.
virtual ~multiphaseInterSystem()
Destructor.
compressibleTurbulenceModel * turb_
Turbulence model.
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol] of the mixture.
dimensionedScalar Prt_
Turbulent Prandt number.
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const =0
Return interfacial source mass rate per phase pair.
static const word phasePropertiesName
Default name of the phase properties dictionary.
const dimensionedScalar & Prt() const
Return Prandt number.
virtual void alphaTransfer(SuSpTable &Su, SuSpTable &Sp)=0
Calculate mass transfer for alpha's.
virtual tmp< volScalarField > Cp() const
Return Cp of the mixture.
void generatePairs(const dictTable &modelDicts)
Generate pairs.
tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
virtual word thermoName() const
Return the name of the thermo physics.
const fvMesh & mesh_
Reference to the mesh.
virtual tmp< volScalarField > surfaceTensionCoeff(const phasePairKey &key) const
Return the surface tension coefficient.
HashTable< autoPtr< multiphaseInter::surfaceTensionModel >, phasePairKey, phasePairKey::hash > surfaceTensionModelTable
const surfaceScalarField & phi() const
Constant access to the total flux.
virtual void correct()
Correct the mixture thermos.
wordList phaseNames_
Phase names.
tmp< volVectorField > U() const
Mixture U.
const phaseModelTable & phases() const
Constant access the phases.
TypeName("multiphaseInterSystem")
Runtime type information.
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
virtual tmp< volScalarField > mu() const
Dynamic viscosity of mixture [kg/m/s].
compressibleTurbulenceModel * turbulence() const
Return pointer to turbulence model.
const modelType & lookupSubModel(const phasePair &key) const
Access a sub model between a phase pair.
tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
interfacePorousModelTable interfacePorousModelTable_
Interface porous models.
virtual tmp< fvScalarMatrix > heatTransfer(const volScalarField &T)=0
Return the heat transfer matrices.
virtual tmp< scalarField > Cp(const scalarField &p, const scalarField &T, const labelList &cells) const
Heat capacity using pressure and temperature.
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
virtual void correctMassSources(const volScalarField &T)=0
Correct mass sources.
const surfaceScalarField & rhoPhi() const
Constant access to the mixture mass flux.
surfaceTensionModelTable surfaceTensionModels_
Surface tension models.
const phasePairTable & totalPhasePairs() const
Constant access the total phase pairs.
virtual bool includeVolChange()=0
Add volume change in pEq.
virtual void massSpeciesTransfer(const multiphaseInter::phaseModel &phase, volScalarField::Internal &Su, volScalarField::Internal &Sp, const word speciesName)=0
Calculate mass transfer for species.
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
virtual void correctTurbulence()
Correct the turbulence.
virtual tmp< volScalarField > coeffs(const word &key) const
Return coefficients (1/rho)
tmp< volScalarField > alphaEff() const
Effective thermal turbulent diffusivity of mixture [kg/m/s].
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
HashTable< autoPtr< multiphaseInter::phaseModel > > generatePhaseModels(const wordList &names) const
Generate the phases.
phasePairTable phasePairs_
Phase pairs.
surfaceScalarField phi_
Mixture total volumetric flux.
HashTable< autoPtr< multiphaseInter::phaseModel > > phaseModelTable
HashTable< autoPtr< porousModel >, phasePairKey, phasePairKey::hash > interfacePorousModelTable
HashTable< volScalarField::Internal > SuSpTable
virtual tmp< volScalarField > hc() const
Chemical enthalpy of the mixture [J/kg].
tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
virtual void solve()=0
Solve for the phase transport equations.
void calcMu()
Calculate and return the laminar viscosity.
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
virtual volScalarField & he()
Return access to the internal energy field [J/Kg].
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [J/m/s/K].
virtual tmp< scalarField > THE(const scalarField &h, const scalarField &p, const scalarField &T0, const labelList &cells) const
Temperature from enthalpy/internal energy for cell-set.
volScalarField mu_
Dynamic viscocity.
void addInterfacePorosity(fvVectorMatrix &UEqn)
Add interface porosity on phasePair.
const fvMesh & mesh() const
Return mesh.
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
Interface normal surface vector.
virtual tmp< volScalarField > rho() const
Return the mixture density.
tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
surfaceScalarField rhoPhi_
Mixture total mass flux.
tmp< surfaceVectorField > nHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Interface normal surface vector.
virtual tmp< scalarField > rhoEoS(const scalarField &p, const scalarField &T, const labelList &cells) const
Density from pressure and temperature.
tmp< volScalarField > kappaEff() const
Effective thermal turbulent diffusivity for temperature.
virtual tmp< volScalarField > Cv() const
Return Cv of the mixture.
tmp< volVectorField > nVolHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Interface normal volField vector.
void setTurbulence(compressibleTurbulenceModel &turb)
Set turbulence model.
phaseModelTable phaseModels_
Phase models.
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
virtual bool read()
Read base phaseProperties dictionary.
virtual bool isochoric() const
Return true if the equation of state is isochoric for all phasses.
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
virtual bool incompressible() const
Return true if the equation of state is incompressible for all.
tmp< surfaceScalarField > generatePhi(const HashTable< autoPtr< multiphaseInter::phaseModel > > &phaseModels) const
Generate the mixture flux.
phasePairTable totalPhasePairs_
Total ordered phase pairs in the system.
An ordered or unorder pair of phase names. Typically specified as follows.
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.
A class for handling words, derived from Foam::string.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Forward declarations of fvMatrix specializations.
List< label > labelList
A List of labels.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Hashing functor for phasePairKey.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.