38 #ifndef phaseInterface_H
39 #define phaseInterface_H
51 class phaseInterfaceKey;
99 const word& oldSeparator,
278 template<
class ModelType,
class Derived>
281 if (!isA<Derived>(*
this))
284 <<
"Constructing " << ModelType::typeName
285 <<
" for interface " <<
name()
286 <<
" which is not of the required type "
290 return refCast<const Derived>(*
this);
Generic GeometricField class.
An STL-conforming hash table.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
A 2-tuple for storing two objects of different types.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Mesh data needed to do the Finite Volume discretisation.
Word-pair based class used for keying interface models in hash tables.
const phaseModel & operator*() const
bool operator!=(const const_iterator &) const
const_iterator & operator++()
const phaseModel & otherPhase() const
const phaseModel & operator()() const
bool operator==(const const_iterator &) const
label index() const
Return the current index.
Class used for construction of PtrLists of phaseInterfaces.
iNew(const phaseSystem &fluid)
autoPtr< phaseInterface > operator()(Istream &is) const
Class to represent an interface between phases. Derivations can further specify the configuration of ...
static word oldNamePartsToName(const phaseSystem &fluid, const wordList &oldNameParts)
Convert old-format interface name parts to an interface name. Used.
static const phaseModel & getPhase1(const phaseModel &phase1, const phaseModel &phase2)
Get a reference to phase1 after sorting the phases by index.
static Tuple2< const phaseModel &, const phaseModel & > identifyPhases(const phaseSystem &fluid, const word &name, const wordList &separators)
Return references to the phases associated with a given name, and a.
const uniformDimensionedVectorField & g() const
Return gravitational acceleration.
const_iterator end() const
const_iterator set to beyond the end of the pair
const phaseSystem & fluid() const
Return the phase system.
TypeName("phaseInterface")
Runtime type information.
const phaseModel & otherPhase(const phaseModel &phase) const
Return the other phase relative to the given phase.
static word separator()
Return the separator that delimits this interface's name.
virtual autoPtr< phaseInterface > clone() const
Clone function.
virtual word name() const
Name.
static wordList nameToNameParts(const phaseSystem &fluid, const word &name)
Split an interface name and return all its parts.
static word nameToTypeName(const phaseSystem &fluid, const word &name)
Convert an interface name into a type name. Essentially just.
virtual const rhoFluidThermo & thermo2() const
Return the thermo for phase 2.
const_iterator begin() const
const_iterator set to the beginning of the pair
virtual const rhoFluidThermo & thermo1() const
Return the thermo for phase 1.
const_iterator cend() const
const_iterator set to beyond the end of the pair
static bool addHeadSeparator(const word &separator)
Add a head separator to the list.
virtual ~phaseInterface()
Destructor.
tmp< volScalarField > magUr() const
Relative velocity magnitude.
static word separatorsToTypeName(const wordList &separators)
Convert a list of separators into a type name.
static const phaseModel & getPhase2(const phaseModel &phase1, const phaseModel &phase2)
Get a reference to phase2 after sorting the phases by index.
virtual const volScalarField & rho2() const
Return the density of phase 2.
declareRunTimeSelectionTable(autoPtr, phaseInterface, word,(const phaseSystem &fluid, const word &name),(fluid, name))
tmp< volScalarField > sigma() const
Surface tension coefficient.
static wordList nameToSeparators(const phaseSystem &fluid, const word &name)
Split an interface name and return its separators.
phaseInterface(const phaseModel &phase1, const phaseModel &phase2)
Construct from phases.
label index(const phaseModel &phase) const
Return the index of the given phase. Generates a FatalError if.
static bool addOldSeparatorToSeparator(const word &oldSeparator, const word &separator)
Add a old separator to separator to the table.
const_iterator cbegin() const
const_iterator set to the beginning of the pair
const fvMesh & mesh() const
Return the mesh.
bool contains(const phaseModel &phase) const
Return true if this phaseInterface contains the given phase.
tmp< volScalarField > rho() const
Average density.
virtual const volScalarField & rho1() const
Return the density of phase 1.
const Derived & modelCast() const
Cast to derived type for use in a model.
static word namePartsToName(const phaseSystem &fluid, const wordList &nameParts)
Convert interface name parts to an interface name.
virtual const volScalarField & alpha2() const
Return the volume fraction of phase 2.
static autoPtr< phaseInterface > New(const phaseSystem &fluid, const word &name)
Select given fluid and name.
const phaseModel & phase1() const
Return phase 1.
const phaseModel & phase2() const
Return phase 2.
virtual const volScalarField & alpha1() const
Return the volume fraction of phase 1.
Class to represent a system of phases and model interfacial transfers between them.
Base-class for fluid thermodynamic properties based on density.
A class for managing temporary objects.
A class for handling words, derived from string.
static const word null
An empty word.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
errorManipArg< error, int > exit(error &err, const int errNo=1)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Macros to ease declaration of run-time selection tables.