Go to the documentation of this file.
41 namespace regionModels
43 namespace pyrolysisModels
115 Qrp =
max(Qrp, scalar(0.0));
123 label localPyrolysisFaceI = 0;
133 const scalar Qr0 = Qrp[faceI];
134 point Cf0 = Cf[faceI];
136 scalar kappaInt = 0.0;
140 const point& Cf1 = cellC[cellI];
141 const scalar
delta =
mag(Cf1 - Cf0);
143 Qr_[cellI] = Qr0*
exp(-kappaInt);
168 label totalFaceId = 0;
179 scalar massInt = 0.0;
183 massInt += RRiGas[cellI]*cellVol[cellI];
184 phiHsGas_[cellI] += massInt*HsiGas[cellI];
187 phiGasp[faceI] += massInt;
191 Info<<
" Gas : " << gasTable[gasI]
192 <<
" on patch : " << patchI
193 <<
" mass produced at face(local) : "
195 <<
" is : " << massInt
196 <<
" [kg/s] " <<
endl;
246 Info<<
"reactingOneDim::solveContinuity()" <<
endl;
278 Info<<
"reactingOneDim::solveSpeciesMass()" <<
endl;
283 for (
label i=0; i<
Ys_.size()-1; i++)
319 Info<<
"reactingOneDim::solveEnergy()" <<
endl;
387 const word& modelType,
389 const word& regionType
394 solidThermo_(solidChemistry_->solidThermo()),
408 Ys_(solidThermo_.composition().Y()),
409 h_(solidThermo_.he()),
473 totalGasMassFlux_(0.0),
477 useChemistrySolvers_(
true)
488 const word& modelType,
491 const word& regionType
496 solidThermo_(solidChemistry_->solidThermo()),
510 Ys_(solidThermo_.composition().Y()),
511 h_(solidThermo_.he()),
573 totalGasMassFlux_(0.0),
577 useChemistrySolvers_(
true)
612 Info<<
"\nPyrolysis region: " <<
type() <<
"added mass : "
613 << massAdded <<
endl;
622 scalar
DiNum = -GREAT;
724 Info<<
"pyrolysis min/max(T) = "
736 Info<<
indent <<
"Total gas mass produced [kg] = "
738 <<
indent <<
"Total solid mass lost [kg] = "
740 <<
indent <<
"Total pyrolysis gases [kg/s] = "
742 <<
indent <<
"Total heat release rate [J/s] = "
Base class for pyrolysis models.
virtual ~reactingOneDim()
Destructor.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
virtual tmp< volScalarField > kappa() const
Return the region thermal conductivity [W/m/k].
Selector class for relaxation factors, solver type and solution.
A class for handling words, derived from string.
virtual void preEvolveRegion()
Pre-evolve region.
virtual const surfaceScalarField & phiGas() const
Return the total gas mass flux to primary region [kg/m2/s].
void solveSpeciesMass()
Solve solid species mass conservation.
scalar maxDiff_
Maximum diffussivity.
dimensionedScalar deltaT() const
Return time step.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
bool useChemistrySolvers_
Use chemistry solvers (ode or sequential)
const dimensionSet dimEnergy
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
bool read(const char *, int32_t &)
void updatePhiGas()
Update enthalpy flux for pyrolysis gases.
virtual const volScalarField & T() const
Temperature [K].
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Calculate the divergence of the given field.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const Time & time() const
Return the reference to the time database.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
void solveContinuity()
Solve continuity equation.
virtual scalar solidRegionDiffNo() const
Mean diffusion number of the solid region.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type & value() const
Return const reference to value.
dimensionedScalar addedGasMass_
Cumulative mass generation of the gas phase [kg].
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar exp(const dimensionedScalar &ds)
volScalarField rho_
Density [kg/m3].
dimensioned< scalar > mag(const dimensioned< Type > &)
bool gasHSource_
Add gas enthalpy source term.
const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
void updateFields()
Update submodels.
bool QrHSource_
Add in depth radiation source term.
labelListList boundaryFaceCells_
Global cell IDs.
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
defineTypeNameAndDebug(noPyrolysis, 0)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
void calculateMassTransfer()
Mass check.
virtual void correct()=0
Update properties.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
tmp< labelField > moveMesh(const scalarList &deltaV, const scalar minDelta=0.0)
Move mesh points according to change in cell volumes.
virtual scalar maxDiff() const
Return max diffusivity allowed in the solid.
A wordList with hashed indices for faster lookup by name.
volScalarField Yt(0.0 *Y[0])
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
label nNonOrthCorr_
Number of non-orthogonal correctors.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
virtual const volScalarField & T() const
Return const temperature [K].
PtrList< volScalarField > & Ys_
List of solid components.
bool read()
Read control parameters from dictionary.
static autoPtr< basicSolidChemistryModel > New(const fvMesh &mesh, const word &phaseName=word::null)
Selector.
const volScalarField & rho() const
Fields.
virtual void preEvolveRegion()
Pre-evolve region.
InternalField & internalField()
Return internal field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read()
Read control parameters.
scalar deltaTValue() const
Return time step value.
scalar minimumDelta_
Minimum delta for combustion.
virtual tmp< volScalarField > kappa() const =0
Thermal diffusivity for temperature of mixture [J/m/s/K].
virtual tmp< volScalarField > kappaRad() const
Return the region absorptivity [1/m].
const dictionary & solution() const
Return the solution dictionary.
dimensionedScalar totalHeatRR_
Total heat release rate [J/s].
const surfaceScalarField & nMagSf() const
Return the face area magnitudes / [m2].
solidReactionThermo & solidThermo_
Reference to solid thermo.
A list of keyword definitions, which are a keyword followed by any number of values (e....
const dictionary & controlDict() const
volScalarField chemistrySh_
Heat release [J/s/m3].
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Switch infoOutput_
Active information output.
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
const double e
Elementary charge.
Reacting, 1-D pyrolysis model.
Ostream & indent(Ostream &os)
Indent stream.
void correctBoundaryConditions()
Correct boundary field.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual void evolveRegion()
Evolve the pyrolysis equations.
void max(const dimensioned< Type > &)
surfaceScalarField phiGas_
Total gas mass flux to the primary region [kg/m2/s].
const vectorField & cellCentres() const
virtual const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
const Time & time_
Reference to the time database.
Calculate the laplacian of the given field.
reactingOneDim(const reactingOneDim &)
Disallow default bitwise copy construct.
volScalarField phiHsGas_
Sensible enthalpy gas flux [J/m2/s].
void updateMesh(const scalarField &mass0)
Update/move mesh based on change in mass.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedScalar lostSolidMass_
Cumulative lost mass of the condensed phase [kg].
label k
Boltzmann constant.
virtual void info()
Provide some feedback.
Volume integrate volField creating a volField.
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
void updateQr()
Update radiative flux in pyrolysis region.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Switch moveMesh_
Flag to allow mesh movement.
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
virtual volScalarField & p()
Pressure [Pa].
virtual scalar addMassSources(const label patchI, const label faceI)
External hook to add mass to the primary region.
const dictionary & coeffs() const
Return the model coefficients dictionary.
void solveEnergy()
Solve energy.
autoPtr< basicSolidChemistryModel > solidChemistry_
Reference to the solid chemistry model.
void readReactingOneDimControls()
Read model controls.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
const dimensionSet dimVolume(pow3(dimLength))
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Generic GeometricField class.
void clear() const
If object pointer points to valid object:
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
addToRunTimeSelectionTable(pyrolysisModel, noPyrolysis, mesh)
scalar totalGasMassFlux_
Total mass gas flux at the pyrolysing walls [kg/s].
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual const tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
const word & name() const
Return reference to name.
labelList primaryPatchIDs_
List of patch IDs on the primary region coupled to this region.
volScalarField Qr_
Coupled region radiative heat flux [W/m2].
const fvMesh & regionMesh() const
Return the region mesh database.