Go to the documentation of this file.
37 namespace combustionModels
42 template<
class ReactionThermo,
class ThermoType>
45 const word& modelType,
48 const word& combustionProperties
58 reactionRateFlameArea_
71 this->
thermo().phasePropertyName(
"ft"),
82 Cv_(this->coeffs().getScalar(
"Cv")),
87 ftVarMin_(this->coeffs().getScalar(
"ftVarMin"))
93 template<
class ReactionThermo,
class ThermoType>
100 template<
class ReactionThermo,
class ThermoType>
101 void FSD<ReactionThermo, ThermoType>::calculateSourceNorm()
103 this->singleMixturePtr_->fresCorrect();
105 const label fuelI = this->singleMixturePtr_->fuelIndex();
114 (
s*YFuel - (YO2 - YO2OxiStream_))/(
s*YFuelFuelStream_ + YO2OxiStream_);
126 (ft_*cAux*mgft)().weightedAverage(this->
mesh().V())
127 /((ft_*cAux)().weightedAverage(this->
mesh().V()) + SMALL)
141 reactionRateFlameArea_->correct(
sigma);
143 const volScalarField& omegaFuel = reactionRateFlameArea_->omega();
146 const scalar ftStoich =
147 YO2OxiStream_.value()
149 s.value()*YFuelFuelStream_.value() + YO2OxiStream_.value()
158 this->
thermo().phasePropertyName(
"Pc"),
177 this->
thermo().phasePropertyName(
"omegaFuelBar"),
211 scalar deltaFt = 1.0/ftDim_;
215 if (ft_[celli] > ftMin_ && ft_[celli] < ftMax_)
217 scalar ftCell = ft_[celli];
219 if (ftVar[celli] > ftVarMin_)
221 scalar ftVarc = ftVar[celli];
223 max(ftCell*(ftCell*(1.0 - ftCell)/ftVarc - 1.0), 0.0);
224 scalar
b =
max(a/ftCell - a, 0.0);
226 for (
int i=1; i<ftDim_; i++)
228 scalar ft = i*deltaFt;
229 pc[celli] +=
pow(ft, a-1.0)*
pow(1.0 - ft,
b - 1.0)*deltaFt;
232 for (
int i=1; i<ftDim_; i++)
234 scalar ft = i*deltaFt;
235 omegaFuelBar[celli] +=
236 omegaFuel[celli]/omegaF[celli]
240 /(2.0*
sqr(0.01*omegaF[celli]))
243 *
pow(1.0 - ft,
b - 1.0)
246 omegaFuelBar[celli] /=
max(pc[celli], 1
e-4);
250 omegaFuelBar[celli] =
251 omegaFuel[celli]/omegaF[celli]
252 *
exp(-
sqr(ftCell - ftStoich)/(2.0*
sqr(0.01*omegaF[celli])));
257 omegaFuelBar[celli] = 0.0;
264 List<label> productsIndex(2, label(-1));
267 forAll(this->singleMixturePtr_->specieProd(), specieI)
269 if (this->singleMixturePtr_->specieProd()[specieI] < 0)
271 productsIndex[i] = specieI;
279 scalar YprodTotal = 0;
282 YprodTotal += this->singleMixturePtr_->Yprod0()[productsIndex[j]];
287 if (ft_[celli] < ftStoich)
289 pc[celli] = ft_[celli]*(YprodTotal/ftStoich);
293 pc[celli] = (1.0 - ft_[celli])*(YprodTotal/(1.0 - ftStoich));
297 tmp<volScalarField> tproducts
303 this->
thermo().phasePropertyName(
"products"),
318 label specieI = productsIndex[j];
325 max(scalar(1) - products/
max(pc, scalar(1
e-5)), scalar(0))
328 pc =
min(C_*
c, scalar(1));
332 this->wFuel_ == mgft*pc*omegaFuelBar;
336 template<
class ReactionThermo,
class ThermoType>
343 calculateSourceNorm();
348 template<
class ReactionThermo,
class ThermoType>
353 this->coeffs().readEntry(
"Cv", Cv_);
354 this->coeffs().readEntry(
"ftVarMin", ftVarMin_);
355 reactionRateFlameArea_->read(this->coeffs());
Defines the attributes of an object for which implicit objectRegistry management is supported,...
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from Foam::string.
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))
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
static constexpr const zero Zero
psiReactionThermo & thermo
static const word propertiesName
Calculate the divergence of the given field.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
dimensionedScalar exp(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Base class for combustion models using singleStepReactingMixture.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
const dimensionedScalar b
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label max(const labelHashSet &set, label maxValue=labelMin)
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
LESModel< EddyDiffusivity< turbulenceModel > > LESModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
static autoPtr< reactionRateFlameArea > New(const dictionary &dict, const fvMesh &mesh, const combustionModel &combModel)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Internal & ref(const bool updateAccessTime=true)
Calculate the gradient of the given field.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionedScalar e
const dimensionedScalar c
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
const dimensionSet dimVolume(pow3(dimLength))
Generic GeometricField class.
compressible::turbulenceModel & turb
Abstract base class for turbulence models (RAS, LES and laminar).
const dimensionSet dimless
Flame Surface Dennsity (FDS) combustion model.