43#ifndef thermoSingleLayer_H
44#define thermoSingleLayer_H
55namespace surfaceFilmModels
59class filmViscosityModel;
60class heatTransferModel;
61class phaseChangeModel;
62class filmRadiationModel;
244 const word& modelType,
247 const word& regionType,
315 const scalar massSource,
316 const vector& momentumSource,
317 const scalar pressureSource,
318 const scalar energySource
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Mesh data needed to do the Finite Volume discretisation.
Base class for film radiation models.
Base class for film heat transfer models.
Kinematic form of single-cell layer surface film model.
Base class for surface film phase change models.
const dimensionedVector & g() const
Return the acceleration due to gravity.
Thermodynamic form of single-cell layer surface film model.
volScalarField TPrimary_
Temperature [K].
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
volScalarField kappa_
Thermal conductivity [W/m/K].
autoPtr< heatTransferModel > htcw_
Heat transfer coefficient between wall and film [W/m2/K].
virtual const volScalarField & hs() const
Return the film sensible enthalpy [J/kg].
virtual tmp< fvScalarMatrix > q(volScalarField &hs) const
Return the wall/surface heat transfer term for the enthalpy equation.
virtual void updateSurfaceTemperatures()
Correct the film surface and wall temperatures.
autoPtr< filmRadiationModel > radiation_
Radiation.
virtual void addSources(const label patchi, const label facei, const scalar massSource, const vector &momentumSource, const scalar pressureSource, const scalar energySource)
External hook to add sources to the film.
const SLGThermo & thermo() const
Return const reference to the SLGThermo object.
tmp< scalarField > qconvp(const label patchi) const
Return the convective heat energy from primary region to film.
volScalarField Tw_
Temperature - wall [K].
virtual const volScalarField & kappa() const
Return the film thermal conductivity [W/m/K].
virtual ~thermoSingleLayer()
Destructor.
virtual const volScalarField & T() const
Return the film mean temperature [K].
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
bool hydrophilic_
Activation flag.
volScalarField Cp_
Specific heat capacity [J/kg/K].
virtual void solveEnergy()
Solve energy equation.
volScalarField hsSp_
Energy [J/m2/s].
const PtrList< volScalarField > & YPrimary() const
Specie mass fractions [0-1].
scalar hydrophilicWetScale_
Length scale applied to deltaWet_ to determine when a dry.
virtual const volScalarField & Tw() const
Return the film wall temperature [K].
const heatTransferModel & htcw() const
Return const access to the (wall) heat transfer model.
const volScalarField & hsSp() const
Energy [J/m2/s].
const phaseChangeModel & phaseChange() const
Return const access to the phase change model.
volScalarField hsSpPrimary_
Energy [J/m2/s].
virtual const volScalarField & Cp() const
Return the film specific heat capacity [J/kg/K].
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
scalar deltaWet_
Threshold film thickness beyond which the film is considered 'wet'.
volScalarField T_
Temperature - mean [K].
virtual void preEvolveRegion()
Pre-evolve film hook.
virtual void correctAlpha()
Correct film coverage field.
PtrList< volScalarField > YPrimary_
List of specie mass fractions [0-1].
const heatTransferModel & htcs() const
Return const access to the (surface) heat transfer model.
volScalarField primaryEnergyTrans_
Film energy transfer.
const filmRadiationModel & radiation() const
Return const access to the radiation model.
const volScalarField & TPrimary() const
Temperature [K].
scalar Tmax_
Maximum temperature limit (optional)
scalar Tmin_
Minimum temperature limit (optional)
scalar hydrophilicDryScale_
Length scale applied to deltaWet_ to determine when a wet.
virtual void correctThermoFields()
Correct the thermo fields.
virtual void info()
Provide some feedback.
virtual void updateSubmodels()
Update the film sub-models.
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
TypeName("thermoSingleLayer")
Runtime type information.
const SLGThermo & thermo_
Reference to the SLGThermo.
volScalarField hs_
Sensible enthalpy [J/kg].
const volScalarField & hsSpPrimary() const
Energy [J/m2/s].
autoPtr< phaseChangeModel > phaseChange_
Phase change.
tmp< scalarField > qconvw(const label patchi) const
Return the convective heat energy from film to wall.
virtual tmp< volScalarField::Internal > Srho() const
Return total mass source - Eulerian phase only.
autoPtr< filmViscosityModel > viscosity_
Viscosity model.
autoPtr< heatTransferModel > htcs_
Heat transfer coefficient between film surface and primary.
virtual bool read()
Read control parameters from dictionary.
virtual void correctHsForMappedT()
Correct sensible enthalpy for mapped temperature fields.
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
virtual void evolveRegion()
Evolve the film equations.
volScalarField Ts_
Temperature - surface [K].
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.