42template<
class ReactionThermo,
class ThermoType>
43FSD<ReactionThermo, ThermoType>::FSD
45 const word& modelType,
48 const word& combustionProperties
51 singleStepCombustion<ReactionThermo, ThermoType>
58 reactionRateFlameArea_
71 this->
thermo().phasePropertyName(
"ft"),
94template<
class ReactionThermo,
class ThermoType>
101template<
class ReactionThermo,
class ThermoType>
102void FSD<ReactionThermo, ThermoType>::calculateSourceNorm()
104 this->singleMixturePtr_->fresCorrect();
106 const label fuelI = this->singleMixturePtr_->fuelIndex();
115 (
s*YFuel - (YO2 - YO2OxiStream_))/(
s*YFuelFuelStream_ + YO2OxiStream_);
127 (ft_*cAux*mgft)().weightedAverage(this->
mesh().V())
128 /((ft_*cAux)().weightedAverage(this->
mesh().V()) + SMALL)
142 reactionRateFlameArea_->correct(
sigma);
144 const volScalarField& omegaFuel = reactionRateFlameArea_->omega();
147 const scalar ftStoich =
148 YO2OxiStream_.value()
150 s.value()*YFuelFuelStream_.value() + YO2OxiStream_.value()
155 this->
thermo().phasePropertyName(
"Pc"),
160 auto& pc = tPc.ref();
164 this->
thermo().phasePropertyName(
"omegaFuelBar"),
169 auto& omegaFuelBar = tomegaFuel.ref();
192 scalar deltaFt = 1.0/ftDim_;
196 if (ft_[celli] > ftMin_ && ft_[celli] < ftMax_)
198 scalar ftCell = ft_[celli];
200 if (ftVar[celli] > ftVarMin_)
202 scalar ftVarc = ftVar[celli];
204 max(ftCell*(ftCell*(1.0 - ftCell)/ftVarc - 1.0), 0.0);
205 scalar
b =
max(a/ftCell - a, 0.0);
207 for (
int i=1; i<ftDim_; i++)
209 scalar ft = i*deltaFt;
210 pc[celli] +=
pow(ft, a-1.0)*
pow(1.0 - ft,
b - 1.0)*deltaFt;
213 for (
int i=1; i<ftDim_; i++)
215 scalar ft = i*deltaFt;
216 omegaFuelBar[celli] +=
217 omegaFuel[celli]/omegaF[celli]
221 /(2.0*
sqr(0.01*omegaF[celli]))
224 *
pow(1.0 - ft,
b - 1.0)
227 omegaFuelBar[celli] /=
max(pc[celli], 1
e-4);
231 omegaFuelBar[celli] =
232 omegaFuel[celli]/omegaF[celli]
233 *
exp(-
sqr(ftCell - ftStoich)/(2.0*
sqr(0.01*omegaF[celli])));
238 omegaFuelBar[celli] = 0.0;
245 List<label> productsIndex(2, label(-1));
248 forAll(this->singleMixturePtr_->specieProd(), specieI)
250 if (this->singleMixturePtr_->specieProd()[specieI] < 0)
252 productsIndex[i] = specieI;
260 scalar YprodTotal = 0;
263 YprodTotal += this->singleMixturePtr_->Yprod0()[productsIndex[j]];
268 if (ft_[celli] < ftStoich)
270 pc[celli] = ft_[celli]*(YprodTotal/ftStoich);
274 pc[celli] = (1.0 - ft_[celli])*(YprodTotal/(1.0 - ftStoich));
280 this->
thermo().phasePropertyName(
"products"),
285 auto& products = tproducts.ref();
289 label specieI = productsIndex[j];
296 max(scalar(1) - products/
max(pc, scalar(1
e-5)), scalar(0))
299 pc =
min(C_*c, scalar(1));
303 this->
wFuel_ == mgft*pc*omegaFuelBar;
307template<
class ReactionThermo,
class ThermoType>
310 this->wFuel_ ==
Zero;
314 calculateSourceNorm();
319template<
class ReactionThermo,
class ThermoType>
324 this->coeffs().readEntry(
"Cv", Cv_);
325 this->coeffs().readEntry(
"ftVarMin", ftVarMin_);
326 reactionRateFlameArea_->read(this->coeffs());
Macros for easy insertion into run-time selection tables.
compressible::turbulenceModel & turb
static autoPtr< CombustionModel > New(ReactionThermo &thermo, const compressibleTurbulenceModel &turb, const word &combustionProperties=combustionPropertiesName)
Selector.
const dimensionSet & dimensions() const noexcept
Return dimensions.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const Time & time() const noexcept
Return Time associated with the objectRegistry.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
virtual ReactionThermo & thermo()
Return access to the thermo package.
const dictionary & coeffs() const
Return const dictionary of the model.
Switch active() const noexcept
Is combustion active?
const fvMesh & mesh() const
Return const access to the mesh database.
virtual void correct()
Correct combustion rate.
virtual ~FSD()
Destructor.
virtual bool read()
Update properties.
singleStepReactingMixture< ThermoType > * singleMixturePtr_
Pointer to singleStepReactingMixture mixture.
volScalarField wFuel_
Fuel consumption rate.
virtual bool read()
Update properties from given dictionary.
Abstract base class for turbulence models (RAS, LES and laminar).
scalar getScalar(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Same as get<scalar>(const word&, keyType::option).
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Abstract class for reaction rate per flame area unit.
static autoPtr< reactionRateFlameArea > New(const dictionary &dict, const fvMesh &mesh, const combustionModel &combModel)
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Calculate the divergence of the given field.
Calculate the gradient of the given field.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
LESModel< EddyDiffusivity< turbulenceModel > > LESModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
const dimensionedScalar c
Speed of light in a vacuum.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
static constexpr const zero Zero
Global zero (0).
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
psiReactionThermo & thermo
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.