Go to the documentation of this file.
37 #ifndef interfaceCompositionModel_H
38 #define interfaceCompositionModel_H
73 TypeName(
"interfaceCompositionModel");
132 const word& speciesName,
139 const word& speciesName,
146 const word& speciesName,
153 const word& speciesName
159 const word& speciesName,
A class for handling words, derived from string.
A class for managing temporary objects.
virtual tmp< volScalarField > D(const word &speciesName) const =0
Mass diffusivity.
declareRunTimeSelectionTable(autoPtr, interfaceCompositionModel, dictionary,(const dictionary &dict, const phasePair &pair),(dict, pair))
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const =0
Mass fraction difference between the interface and the field.
Generic base class for interface composition models. These models describe the composition in phase 1...
CGAL::Exact_predicates_exact_constructions_kernel K
static autoPtr< interfaceCompositionModel > New(const dictionary &dict, const phasePair &pair)
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const =0
Latent heat.
TypeName("interfaceCompositionModel")
Runtime type information.
A wordList with hashed indices for faster lookup by name.
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const =0
The interface mass fraction derivative w.r.t. temperature.
virtual void update(const volScalarField &Tf)=0
Update the composition.
const phasePair & pair_
Phase pair.
A list of keyword definitions, which are a keyword followed by any number of values (e....
const hashedWordList & species() const
Return the transferring species names.
bool transports(word &speciesName) const
Returns whether the species is transported by the model and.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
const hashedWordList speciesNames_
Names of the transferring species.
interfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.
Macros to ease declaration of run-time selection tables.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const =0
Interface mass fraction.
virtual void addMDotL(const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const =0
Add latent heat flow rate to total.
Generic GeometricField class.
virtual ~interfaceCompositionModel()
Destructor.