Go to the documentation of this file.
37 #ifndef surfaceFilmRegionModel_H
38 #define surfaceFilmRegionModel_H
47 namespace regionModels
49 namespace surfaceFilmModels
95 const word& modelType,
98 const word& regionType
118 const scalar massSource,
119 const vector& momentumSource,
120 const scalar pressureSource,
121 const scalar energySource
virtual const volScalarField & Ts() const =0
Return the film surface temperature [K].
Base class for single layer region models.
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
virtual const volScalarField & Cp() const =0
Return the film specific heat capacity [J/kg/K].
virtual const volScalarField & delta() const =0
Return the film thickness [m].
virtual const volVectorField & U() const =0
Return the film velocity [m/s].
virtual const volScalarField & kappa() const =0
Return the film thermal conductivity [W/m/K].
virtual ~surfaceFilmRegionModel()
Destructor.
const dimensionedVector & g() const
Return the acceleration due to gravity.
Base class for surface film models.
virtual const volScalarField & cloudDiameterTrans() const =0
Return the parcel diameters originating from film.
virtual void evolve()
Main driver routing to evolve the region - calls other evolves.
virtual const volScalarField & sigma() const =0
Return the film surface tension [N/m].
virtual const volVectorField & Us() const =0
Return the film surface velocity [m/s].
virtual bool read()
Read control parameters from dictionary.
virtual const volScalarField & cloudMassTrans() const =0
Return the film mass available for transfer.
Base class for surface film models.
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].
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.
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
virtual const volScalarField & hs() const =0
Return the film surface temperature [J/kg].
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
Mesh data needed to do the Finite Volume discretisation.
TypeName("surfaceFilmRegionModel")
Runtime type information.
virtual const volScalarField & Tw() const =0
Return the film wall temperature [K].
virtual tmp< volScalarField > primaryMassTrans() const =0
Return mass transfer source - Eulerian phase only.
const dimensionedVector & g_
Acceleration due to gravity [m/s2].