48 void Foam::fv::constantbXiIgnition::readCoeffs(
const dictionary&
dict)
61 const word& modelType,
69 start_(
"start",
mesh().time().userUnits(),
dict),
70 duration_(
"duration",
mesh().time().userUnits(),
dict),
84 (curTime > start_ - 0.5*deltaT)
85 && (curTime < start_ +
max(duration_, deltaT))
95 return (curTime > start_ - 0.5*deltaT);
106 if (!igniting())
return;
120 const scalar strength = strength_.value();
121 const scalar duration = duration_.value();
126 const scalar Vc = V[celli];
127 Sp[celli] -= Vc*rhou[celli]*strength/(duration*(
b[celli] + 0.001));
141 XiCorrModel_->XiCorr(Xi,
b, mgb);
151 zone_.topoChange(map);
152 XiCorrModel_->topoChange(map);
159 XiCorrModel_->mapMesh(map);
168 zone_.distribute(map);
169 XiCorrModel_->distribute(map);
176 XiCorrModel_->movePoints();
185 zone_.read(coeffs(
dict));
186 XiCorrModel_->read(coeffs(
dict));
187 readCoeffs(coeffs(
dict));
#define forAll(list, i)
Loop across all elements in list.
Generic GeometricField class.
dimensionedScalar deltaT() const
Return time step.
Base class for ignition kernel flame wrinkling Xi correction.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
void read(const dictionary &, const unitConversion &defaultUnits=NullObjectRef< unitConversion >())
Update the value of dimensioned<Type>
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return the top-level database.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Finite volume model abstract base class.
virtual bool read(const dictionary &dict)
Read source dictionary.
const fvMesh & mesh() const
Return const access to the mesh database.
Abstract base-class for ignition models for the Weller b-Xi combustion models.
Simple constant rate ignition model for the Weller b-Xi combustion models.
virtual bool movePoints()
Update for mesh motion.
dimensionedScalar strength_
Ignition strength.
virtual void XiCorr(volScalarField &Xi, const volScalarField &b, const volScalarField &mgb) const
dimensionedScalar duration_
Ignition duration.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
constantbXiIgnition(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual bool ignited() const
Return true during the combustion duration.
dimensionedScalar start_
Ignition start time.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual void addSup(const volScalarField &rho, const volScalarField &b, fvMatrix< scalar > &eqn) const
Add ignition contribution to b equation.
virtual bool igniting() const
Return true during the ignition duration.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const dimensionSet dimless
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.