Loading...
Searching...
No Matches
diffusionMulticomponent.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2015-2018 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#include "fvcGrad.H"
30#include "reactingMixture.H"
31#include "fvCFD.H"
33
34// * * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * //
35
36template<class ReactionThermo, class ThermoType>
37void Foam::combustionModels::
38diffusionMulticomponent<ReactionThermo, ThermoType>::init()
39{
40 // Load default values
41 this->coeffs().readIfPresent("Ci", Ci_);
42 this->coeffs().readIfPresent("YoxStream", YoxStream_);
43 this->coeffs().readIfPresent("YfStream", YfStream_);
44 this->coeffs().readIfPresent("sigma", sigma_);
45 this->coeffs().readIfPresent("ftCorr", ftCorr_);
46 this->coeffs().readIfPresent("alpha", alpha_);
47 this->coeffs().readIfPresent("laminarIgn", laminarIgn_);
48
50
51 const speciesTable& species = this->thermo().composition().species();
52
53 scalarList specieStoichCoeffs(species.size());
54 const label nReactions = reactions_.size();
55
56 for (label k=0; k < nReactions; k++)
57 {
58 RijPtr_.set
59 (
60 k,
62 (
63 this->mesh_.newIOobject("Rijk" + Foam::name(k)),
64 this->mesh_,
67 )
68 );
69
70 RijPtr_[k].storePrevIter();
71
72 const List<specieCoeffs>& lhs = reactions_[k].lhs();
73 const List<specieCoeffs>& rhs = reactions_[k].rhs();
74
75 const label fuelIndex = species.find(fuelNames_[k]);
76 const label oxidantIndex = species.find(oxidantNames_[k]);
77
78 const scalar Wu = specieThermo_[fuelIndex].W();
79 const scalar Wox = specieThermo_[oxidantIndex].W();
80
81 forAll(lhs, i)
82 {
83 const label specieI = lhs[i].index;
84 specieStoichCoeffs[specieI] = -lhs[i].stoichCoeff;
85 qFuel_[k] +=
86 specieThermo_[specieI].hc()*lhs[i].stoichCoeff/Wu;
87 }
88
89 forAll(rhs, i)
90 {
91 const label specieI = rhs[i].index;
92 specieStoichCoeffs[specieI] = rhs[i].stoichCoeff;
93 qFuel_[k] -=
94 specieThermo_[specieI].hc()*rhs[i].stoichCoeff/Wu;
95 }
96
97 Info << "Fuel heat of combustion : " << qFuel_[k] << endl;
98
99 s_[k] =
100 (Wox*mag(specieStoichCoeffs[oxidantIndex]))
101 / (Wu*mag(specieStoichCoeffs[fuelIndex]));
102
103 Info << "stoichiometric oxygen-fuel ratio : " << s_[k] << endl;
104
105 stoicRatio_[k] = s_[k]*YfStream_[k]/YoxStream_[k];
106
107 Info << "stoichiometric air-fuel ratio : " << stoicRatio_[k] << endl;
108
109 const scalar fStoich = 1.0/(1.0 + stoicRatio_[k]);
110
111 Info << "stoichiometric mixture fraction : " << fStoich << endl;
113}
114
115
116// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
117
118template<class ReactionThermo, class ThermoType>
121(
122 const word& modelType,
123 ReactionThermo& thermo,
125 const word& combustionProperties
126)
127:
128 ChemistryCombustion<ReactionThermo>
129 (
130 modelType,
131 thermo,
132 turb,
133 combustionProperties
134 ),
135 reactions_
136 (
137 dynamic_cast<const reactingMixture<ThermoType>&>(thermo)
138 ),
139 specieThermo_
140 (
141 dynamic_cast<const reactingMixture<ThermoType>&>(thermo).
142 speciesData()
143 ),
144 RijPtr_(reactions_.size()),
145 Ci_(reactions_.size(), 1.0),
146 fuelNames_(this->coeffs().lookup("fuels")),
147 oxidantNames_(this->coeffs().lookup("oxidants")),
148 qFuel_(reactions_.size()),
149 stoicRatio_(reactions_.size()),
150 s_(reactions_.size()),
151 YoxStream_(reactions_.size(), 0.23),
152 YfStream_(reactions_.size(), 1.0),
153 sigma_(reactions_.size(), 0.02),
154 oxidantRes_(this->coeffs().lookup("oxidantRes")),
155 ftCorr_(reactions_.size(), Zero),
156 alpha_(1),
157 laminarIgn_(false)
158{
159 init();
160}
161
162
163// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
164
165template<class ReactionThermo, class ThermoType>
166Foam::tmp<Foam::volScalarField> Foam::combustionModels::
167diffusionMulticomponent<ReactionThermo, ThermoType>::tc() const
169 return this->chemistryPtr_->tc();
170}
171
172
173template<class ReactionThermo, class ThermoType>
176{
177 if (this->active())
178 {
180 const speciesTable& species = this->thermo().composition().species();
181
182 const label nReactions = reactions_.size();
183
184 PtrList<volScalarField> RijlPtr(nReactions);
185
186 for (label k=0; k < nReactions; k++)
187 {
188 RijlPtr.set
189 (
190 k,
192 (
193 this->mesh_.newIOobject("Rijl" + Foam::name(k)),
194 this->mesh_,
197 )
198 );
199
200 volScalarField& Rijl = RijlPtr[k];
201
202 // Obtain individual reaction rates for each reaction
203 const label fuelIndex = species.find(fuelNames_[k]);
204
205 if (laminarIgn_)
206 {
207 Rijl.ref() = -this->chemistryPtr_->calculateRR(k, fuelIndex);
208 }
209
210
211 // Look for the fuelStoic
212 const List<specieCoeffs>& rhs = reactions_[k].rhs();
213 const List<specieCoeffs>& lhs = reactions_[k].lhs();
214
215 // Set to zero RR's
216 forAll(lhs, l)
217 {
218 const label lIndex = lhs[l].index;
219 this->chemistryPtr_->RR(lIndex) = Zero;
220 }
221
222 forAll(rhs, l)
223 {
224 const label rIndex = rhs[l].index;
225 this->chemistryPtr_->RR(rIndex) = Zero;
226 }
227 }
228
229
230 for (label k=0; k < nReactions; k++)
231 {
232 const label fuelIndex = species.find(fuelNames_[k]);
233 const label oxidantIndex = species.find(oxidantNames_[k]);
234
235 const volScalarField& Yfuel =
236 this->thermo().composition().Y(fuelIndex);
237
238 const volScalarField& Yox =
239 this->thermo().composition().Y(oxidantIndex);
240
241 const volScalarField ft
242 (
243 "ft" + Foam::name(k),
244 (
245 s_[k]*Yfuel - (Yox - YoxStream_[k])
246 )
247 /(
248 s_[k]*YfStream_[k] + YoxStream_[k]
249 )
250 );
251
252 const scalar sigma = sigma_[k];
253
254 const scalar fStoich = 1.0/(1.0 + stoicRatio_[k]) + ftCorr_[k];
255
256 const volScalarField OAvailScaled
257 (
258 "OAvailScaled",
259 Yox/max(oxidantRes_[k], 1e-3)
260 );
261
262 const volScalarField preExp
263 (
264 "preExp" + Foam::name(k),
265 1.0 + sqr(OAvailScaled)
266 );
267
268 const volScalarField filter
269 (
271 * exp(-sqr(ft - fStoich)/(2*sqr(sigma)))
272 );
273
274 const volScalarField topHatFilter(pos(filter - 1e-3));
275
276 const volScalarField prob("prob" + Foam::name(k), preExp*filter);
277
278 const volScalarField RijkDiff
279 (
280 "RijkDiff",
281 Ci_[k]*this->turbulence().muEff()*prob*
282 (
283 mag(fvc::grad(Yfuel) & fvc::grad(Yox))
284 )
285 *pos(Yox)*pos(Yfuel)
286 );
287
288 volScalarField& Rijk = RijPtr_[k];
289
290 if (laminarIgn_)
291 {
292 Rijk =
293 min(RijkDiff, topHatFilter*RijlPtr[k]*pos(Yox)*pos(Yfuel));
294 }
295 else
296 {
297 Rijk = RijkDiff;
298 }
299
300 Rijk.relax(alpha_);
301
302 if (debug && this->mesh_.time().writeTime())
303 {
304 Rijk.write();
305 ft.write();
306 }
307
308 // Look for the fuelStoic
309 const List<specieCoeffs>& rhs = reactions_[k].rhs();
310 const List<specieCoeffs>& lhs = reactions_[k].lhs();
311
312 scalar fuelStoic = 1.0;
313 forAll(lhs, l)
314 {
315 const label lIndex = lhs[l].index;
316 if (lIndex == fuelIndex)
317 {
318 fuelStoic = lhs[l].stoichCoeff;
319 break;
320 }
321 }
322
323 const scalar MwFuel = specieThermo_[fuelIndex].W();
324
325 // Update left hand side species
326 forAll(lhs, l)
327 {
328 const label lIndex = lhs[l].index;
329
330 const scalar stoichCoeff = lhs[l].stoichCoeff;
331
332 this->chemistryPtr_->RR(lIndex) +=
333 -Rijk*stoichCoeff*specieThermo_[lIndex].W()/fuelStoic/MwFuel;
334
335 }
336
337 // Update right hand side species
338 forAll(rhs, r)
339 {
340 const label rIndex = rhs[r].index;
341
342 const scalar stoichCoeff = rhs[r].stoichCoeff;
343
344 this->chemistryPtr_->RR(rIndex) +=
345 Rijk*stoichCoeff*specieThermo_[rIndex].W()/fuelStoic/MwFuel;
346 }
348 }
349}
350
351
352template<class ReactionThermo, class ThermoType>
355(
357) const
358{
360 auto& Su = tSu.ref();
361
362 if (this->active())
363 {
364 const label specieI =
365 this->thermo().composition().species().find(Y.member());
366
367 Su += this->chemistryPtr_->RR(specieI);
368 }
369
370 return tSu;
371}
372
373
374template<class ReactionThermo, class ThermoType>
378{
379 auto tQdot = volScalarField::New
380 (
381 "Qdot",
383 this->mesh(),
386 );
387
388 if (this->active())
389 {
390 volScalarField& Qdot = tQdot.ref();
391 Qdot = this->chemistryPtr_->Qdot();
392 }
394 return tQdot;
395}
396
397
398template<class ReactionThermo, class ThermoType>
401{
403 {
404 this->coeffs().readIfPresent("Ci", Ci_);
405 this->coeffs().readIfPresent("sigma", sigma_);
406 this->coeffs().readIfPresent("oxidantRes", oxidantRes_);
407 this->coeffs().readIfPresent("ftCorr", ftCorr_);
408 this->coeffs().readIfPresent("alpha", alpha_);
409 this->coeffs().readIfPresent("laminarIgn", laminarIgn_);
410 return true;
411 }
412
413 return false;
414}
415
416
417// ************************************************************************* //
label k
compressible::turbulenceModel & turb
Chemistry model wrapper for combustion models.
autoPtr< BasicChemistryModel< ReactionThermo > > chemistryPtr_
Pointer to chemistry model.
virtual ReactionThermo & thermo()
Return access to the thermo package.
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
void relax(const scalar alpha)
Relax field (for steady-state solution).
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())
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
@ NO_REGISTER
Do not request registration (bool: false).
static word member(const word &name)
Return member (name without the extension).
Definition IOobject.C:313
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
const speciesTable & species() const
Return the table of species.
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
const fvMesh & mesh_
Reference to the mesh database.
const compressibleTurbulenceModel & turbulence() const
Return access to turbulence.
const dictionary & coeffs() const
Return const dictionary of the model.
Switch active() const noexcept
Is combustion active?
virtual bool read()
Update properties from given dictionary.
Diffusion based turbulent combustion model for multicomponent species.
virtual void correct()
Correct combustion rate.
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
virtual tmp< volScalarField > Qdot() const
Heat release rate calculated from fuel consumption rate matrix.
virtual bool read()
Update properties from given dictionary.
Abstract base class for turbulence models (RAS, LES and laminar).
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
label find(const word &val) const
Find index of the value (searches the hash).
virtual basicSpecieMixture & composition()=0
Return the composition of the multi-component mixture.
Lookup type of boundary radiation properties.
Definition lookup.H:60
Foam::reactingMixture.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
PtrList< volScalarField > & Y
dynamicFvMesh & mesh
compressible::turbulenceModel & turbulence
Calculate the gradient of the given field.
zeroField Su
Definition alphaSuSp.H:1
constexpr scalar pi(M_PI)
Namespace for handling debugging switches.
Definition debug.C:45
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
hashedWordList speciesTable
A table of species as a hashedWordList.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
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.
Definition hashSets.C:26
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionSet dimVolume(pow3(dimLength))
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
scalar Qdot
psiReactionThermo & thermo
volScalarField & e
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Hold specie index and its coefficients in the reaction rate expression.
Definition Reaction.H:87