Go to the documentation of this file.
37 #ifndef pyrolysisModel_H
38 #define pyrolysisModel_H
53 namespace regionModels
55 namespace pyrolysisModels
71 void constructMeshObjects();
74 void readPyrolysisControls();
108 const word& modelType,
110 const word& regionType
112 (modelType,
mesh, regionType)
121 const word& modelType,
124 const word& regionType
136 const word& regionType
142 const word& modelType,
144 const word& regionType
150 const word& modelType,
153 const word& regionType
170 const word& regionType =
"pyrolysis"
178 const word& regionType =
"pyrolysis"
227 virtual scalar
maxDiff()
const;
Base class for pyrolysis models.
A class for handling words, derived from Foam::string.
virtual const volScalarField & T() const =0
Return const temperature [K].
A class for managing temporary objects.
virtual const surfaceScalarField & phiGas() const =0
Return the total gas mass flux to primary region [kg/m2/s].
virtual tmp< volScalarField > kappaRad() const =0
Return the region absorptivity [1/m].
virtual ~pyrolysisModel()
Destructor.
virtual const volScalarField & rho() const =0
Return density [kg/m3].
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
autoPtr< pyrolysisModel > clone() const
Return clone.
static autoPtr< pyrolysisModel > New(const fvMesh &mesh, const word ®ionType="pyrolysis")
Return a reference to the selected pyrolysis model.
virtual bool read()
Read control parameters.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
TypeName("pyrolysisModel")
Runtime type information.
Mesh data needed to do the Finite Volume discretisation.
virtual const tmp< volScalarField > Cp() const =0
Return specific heat capacity [J/kg/K].
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
declareRunTimeSelectionTable(autoPtr, pyrolysisModel, mesh,(const word &modelType, const fvMesh &mesh, const word ®ionType),(modelType, mesh, regionType))
Macros to ease declaration of run-time selection tables.
virtual scalar addMassSources(const label patchi, const label facei)
External hook to add mass to the primary region.
virtual tmp< volScalarField > kappa() const =0
Return the region thermal conductivity [W/m/k].
Base class for 1-D region models.
virtual scalar solidRegionDiffNo() const
Mean diffusion number of the solid region.
virtual scalar maxDiff() const
Return max diffusivity allowed in the solid.