42#ifndef PhaseChangeModel_H
43#define PhaseChangeModel_H
59template<
class CloudType>
213#define makePhaseChangeModel(CloudType) \
215 typedef Foam::CloudType::reactingCloudType reactingCloudType; \
216 defineNamedTemplateTypeNameAndDebug \
218 Foam::PhaseChangeModel<reactingCloudType>, \
223 defineTemplateRunTimeSelectionTable \
225 PhaseChangeModel<reactingCloudType>, \
231#define makePhaseChangeModelType(SS, CloudType) \
233 typedef Foam::CloudType::reactingCloudType reactingCloudType; \
234 defineNamedTemplateTypeNameAndDebug(Foam::SS<reactingCloudType>, 0); \
236 Foam::PhaseChangeModel<reactingCloudType>:: \
237 adddictionaryConstructorToTable<Foam::SS<reactingCloudType>> \
238 add##SS##CloudType##reactingCloudType##ConstructorToTable_;
Base class for cloud sub-models.
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Templated phase change model class.
virtual void calculate(const scalar dt, const label celli, const scalar Re, const scalar Pr, const scalar d, const scalar nu, const scalar rho, const scalar T, const scalar Ts, const scalar pc, const scalar Tc, const scalarField &X, const scalarField &solMass, const scalarField &liqMass, scalarField &dMassPC) const =0
Update model.
scalar dMass_
Mass of lagrangian phase converted.
static autoPtr< PhaseChangeModel< CloudType > > New(const dictionary &dict, CloudType &owner)
Selector.
virtual autoPtr< PhaseChangeModel< CloudType > > clone() const =0
Construct and return a clone.
TypeName("phaseChangeModel")
Runtime type information.
virtual scalar dh(const label idc, const label idl, const scalar p, const scalar T) const
Return the enthalpy per unit mass.
const enthalpyTransferType & enthalpyTransfer() const
Return the enthalpy transfer type enumeration.
virtual void info(Ostream &os)
Write injection info to stream.
enthalpyTransferType enthalpyTransfer_
Enthalpy transfer type enumeration.
void addToPhaseChangeMass(const scalar dMass)
Add to phase change mass.
enthalpyTransferType
Enthalpy transfer type.
virtual scalar Tvap(const scalarField &X) const
Return vapourisation temperature.
scalar Sh() const
Sherwood number.
enthalpyTransferType wordToEnthalpyTransfer(const word &etName) const
Convert word to enthalpy transfer type.
virtual ~PhaseChangeModel()=default
Destructor.
virtual scalar TMax(const scalar p, const scalarField &X) const
Return maximum/limiting temperature.
static const wordList enthalpyTransferTypeNames
Name representations of enthalpy transfer types.
declareRunTimeSelectionTable(autoPtr, PhaseChangeModel, dictionary,(const dictionary &dict, CloudType &owner),(dict, owner))
Declare runtime constructor selection table.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const dictionary & dict() const
Return const access to the cloud dictionary.
A class for handling words, derived from Foam::string.
OBJstream os(runTime.globalPath()/outputName)
scalarField Re(const UList< complex > &cf)
Extract real component.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
dimensionedScalar Pr("Pr", dimless, laminarTransport)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.