Go to the documentation of this file.
37 #ifndef surfaceFilmModel_H
38 #define surfaceFilmModel_H
53 namespace regionModels
55 namespace surfaceFilmModels
106 const word& modelType,
109 const word& regionType
111 (modelType,
mesh,
g, regionType)
119 const word& modelType,
122 const word& regionType
133 const word& regionType =
"surfaceFilm"
153 const scalar massSource,
154 const vector& momentumSource,
155 const scalar pressureSource,
156 const scalar energySource
virtual tmp< DimensionedField< scalar, volMesh > > Srho() const
Return total mass source - Eulerian phase only.
virtual const volScalarField & sigma() const =0
Return the film surface tension [N/m].
virtual const volScalarField & Ts() const =0
Return the film surface temperature [K].
virtual scalar CourantNumber() const
Courant number evaluation.
A class for handling words, derived from string.
TypeName("surfaceFilmModel")
Runtime type information.
A class for managing temporary objects.
static autoPtr< surfaceFilmModel > New(const fvMesh &mesh, const dimensionedVector &g, const word ®ionType="surfaceFilm")
Return a reference to the selected surface film model.
virtual const volVectorField & U() const =0
Return the film velocity [m/s].
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
virtual const volScalarField & Tw() const =0
Return the film wall temperature [K].
virtual const volScalarField & delta() const =0
Return the film thickness [m].
virtual ~surfaceFilmModel()
Destructor.
virtual const volVectorField & Us() const =0
Return the film surface velocity [m/s].
const dimensionedVector & g() const
Return the accleration due to gravity.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
virtual const volScalarField & kappa() const =0
Return the film thermal conductivity [W/m/K].
virtual bool read()
Read control parameters from dictionary.
virtual const volScalarField & cloudMassTrans() const =0
Return the film mass available for transfer.
Generic dimensioned Type class.
virtual const volScalarField & cloudDiameterTrans() const =0
Return the parcel diameters originating from film.
Mesh data needed to do the Finite Volume discretisation.
declareRunTimeSelectionTable(autoPtr, surfaceFilmModel, mesh,(const word &modelType, const fvMesh &mesh, const dimensionedVector &g, const word ®ionType),(modelType, mesh, g, regionType))
virtual void addSources(const label patchI, const label faceI, const scalar massSource, const vector &momentumSource, const scalar pressureSource, const scalar energySource)=0
External hook to add sources to the film.
surfaceFilmModel(const surfaceFilmModel &)
Disallow default bitwise copy construct.
const dimensionedVector & g_
Acceleration due to gravity [m/s2].
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
virtual tmp< volScalarField > primaryMassTrans() const =0
Return mass transfer source - Eulerian phase only.
Macros to ease declaration of run-time selection tables.
virtual const volVectorField & Uw() const =0
Return the film wall velocity [m/s].
virtual const volScalarField & T() const =0
Return the film mean temperature [K].
void operator=(const surfaceFilmModel &)
Disallow default bitwise assignment.
virtual tmp< DimensionedField< scalar, volMesh > > Sh() const
Return enthalpy source - Eulerian phase only.
virtual const volScalarField & Cp() const =0
Return the film specific heat capacity [J/kg/K].
Base class for surface film models.
Generic GeometricField class.
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].