Loading...
Searching...
No Matches
StandardChemistryModel.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2020-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "UniformField.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36template<class ReactionThermo, class ThermoType>
37Foam::StandardChemistryModel<ReactionThermo, ThermoType>::StandardChemistryModel
38(
39 ReactionThermo& thermo
40)
41:
42 BasicChemistryModel<ReactionThermo>(thermo),
43 ODESystem(),
44 Y_(this->thermo().composition().Y()),
46 (
47 dynamic_cast<const reactingMixture<ThermoType>&>(this->thermo())
48 ),
50 (
51 dynamic_cast<const reactingMixture<ThermoType>&>
52 (this->thermo()).speciesData()
53 ),
54
55 nSpecie_(Y_.size()),
58 (
59 BasicChemistryModel<ReactionThermo>::template getOrDefault<scalar>
60 (
61 "Treact",
62 0.0
63 )
64 ),
66 c_(nSpecie_),
68{
69 // Create the fields for the chemistry sources
70 forAll(RR_, fieldi)
71 {
72 RR_.set
73 (
74 fieldi,
76 (
78 (
79 "RR." + Y_[fieldi].name(),
80 this->mesh().time().timeName(),
81 this->mesh(),
84 ),
85 this->mesh(),
87 )
88 );
89 }
90
91 Info<< "StandardChemistryModel: Number of species = " << nSpecie_
92 << " and reactions = " << nReaction_ << endl;
93}
94
95
96// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97
98template<class ReactionThermo, class ThermoType>
101{}
102
103
104// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105
106template<class ReactionThermo, class ThermoType>
108(
109 const scalarField& c,
110 const scalar T,
111 const scalar p,
112 scalarField& dcdt
113) const
114{
115 scalar pf, cf, pr, cr;
116 label lRef, rRef;
117
118 dcdt = Zero;
119
120 forAll(reactions_, i)
121 {
122 const Reaction<ThermoType>& R = reactions_[i];
123
124 scalar omegai = omega
125 (
126 R, c, T, p, pf, cf, lRef, pr, cr, rRef
127 );
128
129 forAll(R.lhs(), s)
130 {
131 const label si = R.lhs()[s].index;
132 const scalar sl = R.lhs()[s].stoichCoeff;
133 dcdt[si] -= sl*omegai;
134 }
135
136 forAll(R.rhs(), s)
137 {
138 const label si = R.rhs()[s].index;
139 const scalar sr = R.rhs()[s].stoichCoeff;
140 dcdt[si] += sr*omegai;
141 }
142 }
143}
144
145
146template<class ReactionThermo, class ThermoType>
148(
149 const label index,
150 const scalarField& c,
151 const scalar T,
152 const scalar p,
153 scalar& pf,
154 scalar& cf,
155 label& lRef,
156 scalar& pr,
157 scalar& cr,
158 label& rRef
159) const
160{
162 scalar w = omega(R, c, T, p, pf, cf, lRef, pr, cr, rRef);
163 return(w);
164}
165
166
167template<class ReactionThermo, class ThermoType>
169(
170 const Reaction<ThermoType>& R,
171 const scalarField& c,
172 const scalar T,
173 const scalar p,
174 scalar& pf,
175 scalar& cf,
176 label& lRef,
177 scalar& pr,
178 scalar& cr,
179 label& rRef
180) const
181{
182 const scalar kf = R.kf(p, T, c);
183 const scalar kr = R.kr(kf, p, T, c);
184
185 pf = 1.0;
186 pr = 1.0;
187
188 const label Nl = R.lhs().size();
189 const label Nr = R.rhs().size();
190
191 label slRef = 0;
192 lRef = R.lhs()[slRef].index;
193
194 pf = kf;
195 for (label s = 1; s < Nl; s++)
196 {
197 const label si = R.lhs()[s].index;
198
199 if (c[si] < c[lRef])
200 {
201 const scalar exp = R.lhs()[slRef].exponent;
202 pf *= pow(max(c[lRef], 0.0), exp);
203 lRef = si;
204 slRef = s;
205 }
206 else
207 {
208 const scalar exp = R.lhs()[s].exponent;
209 pf *= pow(max(c[si], 0.0), exp);
210 }
211 }
212 cf = max(c[lRef], 0.0);
213
214 {
215 const scalar exp = R.lhs()[slRef].exponent;
216 if (exp < 1.0)
217 {
218 if (cf > SMALL)
219 {
220 pf *= pow(cf, exp - 1.0);
221 }
222 else
223 {
224 pf = 0.0;
225 }
226 }
227 else
228 {
229 pf *= pow(cf, exp - 1.0);
230 }
231 }
232
233 label srRef = 0;
234 rRef = R.rhs()[srRef].index;
235
236 // Find the matrix element and element position for the rhs
237 pr = kr;
238 for (label s = 1; s < Nr; s++)
239 {
240 const label si = R.rhs()[s].index;
241 if (c[si] < c[rRef])
243 const scalar exp = R.rhs()[srRef].exponent;
244 pr *= pow(max(c[rRef], 0.0), exp);
245 rRef = si;
246 srRef = s;
247 }
248 else
249 {
250 const scalar exp = R.rhs()[s].exponent;
251 pr *= pow(max(c[si], 0.0), exp);
252 }
253 }
254 cr = max(c[rRef], 0.0);
255
256 {
257 const scalar exp = R.rhs()[srRef].exponent;
258 if (exp < 1.0)
260 if (cr>SMALL)
261 {
262 pr *= pow(cr, exp - 1.0);
263 }
264 else
265 {
266 pr = 0.0;
267 }
268 }
269 else
270 {
271 pr *= pow(cr, exp - 1.0);
272 }
274
275 return pf*cf - pr*cr;
276}
277
278
279template<class ReactionThermo, class ThermoType>
281(
282 const scalar time,
283 const scalarField& c,
285) const
286{
287 const scalar T = c[nSpecie_];
288 const scalar p = c[nSpecie_ + 1];
289
290 forAll(c_, i)
291 {
292 c_[i] = max(c[i], 0.0);
293 }
294
295 omega(c_, T, p, dcdt);
296
297 // Constant pressure
298 // dT/dt = ...
299 scalar rho = 0.0;
300 for (label i = 0; i < nSpecie_; i++)
301 {
302 const scalar W = specieThermo_[i].W();
303 rho += W*c_[i];
304 }
305 scalar cp = 0.0;
306 for (label i=0; i<nSpecie_; i++)
308 cp += c_[i]*specieThermo_[i].cp(p, T);
309 }
310 cp /= rho;
311
312 scalar dT = 0.0;
313 for (label i = 0; i < nSpecie_; i++)
314 {
315 const scalar hi = specieThermo_[i].ha(p, T);
316 dT += hi*dcdt[i];
317 }
318 dT /= rho*cp;
319
320 dcdt[nSpecie_] = -dT;
322 // dp/dt = ...
323 dcdt[nSpecie_ + 1] = 0.0;
324}
325
326
327template<class ReactionThermo, class ThermoType>
330 const scalar t,
331 const scalarField& c,
332 scalarField& dcdt,
334) const
335{
336 const scalar T = c[nSpecie_];
337 const scalar p = c[nSpecie_ + 1];
338
339 forAll(c_, i)
340 {
341 c_[i] = max(c[i], 0.0);
342 }
343
344 dfdc = Zero;
345
346 // Length of the first argument must be nSpecie_
347 omega(c_, T, p, dcdt);
348
349 forAll(reactions_, ri)
350 {
351 const Reaction<ThermoType>& R = reactions_[ri];
352
353 const scalar kf0 = R.kf(p, T, c_);
354 const scalar kr0 = R.kr(kf0, p, T, c_);
355
356 forAll(R.lhs(), j)
357 {
358 const label sj = R.lhs()[j].index;
359 scalar kf = kf0;
360 forAll(R.lhs(), i)
361 {
362 const label si = R.lhs()[i].index;
363 const scalar el = R.lhs()[i].exponent;
364 if (i == j)
365 {
366 if (el < 1.0)
367 {
368 if (c_[si] > SMALL)
369 {
370 kf *= el*pow(c_[si], el - 1.0);
371 }
372 else
373 {
374 kf = 0.0;
375 }
376 }
377 else
378 {
379 kf *= el*pow(c_[si], el - 1.0);
380 }
381 }
382 else
383 {
384 kf *= pow(c_[si], el);
385 }
386 }
387
388 forAll(R.lhs(), i)
389 {
390 const label si = R.lhs()[i].index;
391 const scalar sl = R.lhs()[i].stoichCoeff;
392 dfdc(si, sj) -= sl*kf;
393 }
394 forAll(R.rhs(), i)
395 {
396 const label si = R.rhs()[i].index;
397 const scalar sr = R.rhs()[i].stoichCoeff;
398 dfdc(si, sj) += sr*kf;
399 }
400 }
401
402 forAll(R.rhs(), j)
403 {
404 const label sj = R.rhs()[j].index;
405 scalar kr = kr0;
406 forAll(R.rhs(), i)
407 {
408 const label si = R.rhs()[i].index;
409 const scalar er = R.rhs()[i].exponent;
410 if (i == j)
411 {
412 if (er < 1.0)
413 {
414 if (c_[si] > SMALL)
415 {
416 kr *= er*pow(c_[si], er - 1.0);
417 }
418 else
419 {
420 kr = 0.0;
421 }
422 }
423 else
424 {
425 kr *= er*pow(c_[si], er - 1.0);
426 }
427 }
428 else
429 {
430 kr *= pow(c_[si], er);
431 }
432 }
433
434 forAll(R.lhs(), i)
435 {
436 const label si = R.lhs()[i].index;
437 const scalar sl = R.lhs()[i].stoichCoeff;
438 dfdc(si, sj) += sl*kr;
439 }
440 forAll(R.rhs(), i)
441 {
442 const label si = R.rhs()[i].index;
443 const scalar sr = R.rhs()[i].stoichCoeff;
444 dfdc(si, sj) -= sr*kr;
445 }
446 }
447 }
448
449 // Calculate the dcdT elements numerically
450 const scalar delta = 1.0e-3;
451
452 omega(c_, T + delta, p, dcdt_);
453 for (label i=0; i<nSpecie_; i++)
454 {
455 dfdc(i, nSpecie_) = dcdt_[i];
456 }
457
458 omega(c_, T - delta, p, dcdt_);
459 for (label i=0; i<nSpecie_; i++)
460 {
461 dfdc(i, nSpecie_) = 0.5*(dfdc(i, nSpecie_) - dcdt_[i])/delta;
462 }
463
464 dfdc(nSpecie_, nSpecie_) = 0;
465 dfdc(nSpecie_ + 1, nSpecie_) = 0;
466}
467
468
469template<class ReactionThermo, class ThermoType>
472{
473 auto ttc = volScalarField::New
474 (
475 "tc",
477 this->mesh(),
480 );
481 scalarField& tc = ttc.ref();
482
484 const scalarField& rho = trho();
485
486 const scalarField& T = this->thermo().T();
487 const scalarField& p = this->thermo().p();
488
489 const label nReaction = reactions_.size();
490
491 scalar pf, cf, pr, cr;
492 label lRef, rRef;
493
494 if (this->chemistry_)
495 {
496 forAll(rho, celli)
497 {
498 const scalar rhoi = rho[celli];
499 const scalar Ti = T[celli];
500 const scalar pi = p[celli];
501
502 scalar cSum = 0.0;
503
504 for (label i=0; i<nSpecie_; i++)
505 {
506 c_[i] = rhoi*Y_[i][celli]/specieThermo_[i].W();
507 cSum += c_[i];
508 }
509
510 forAll(reactions_, i)
511 {
512 const Reaction<ThermoType>& R = reactions_[i];
513
514 omega(R, c_, Ti, pi, pf, cf, lRef, pr, cr, rRef);
515
516 forAll(R.rhs(), s)
517 {
518 tc[celli] += R.rhs()[s].stoichCoeff*pf*cf;
519 }
520 }
521
522 tc[celli] = nReaction*cSum/tc[celli];
523 }
524 }
525
526 ttc.ref().correctBoundaryConditions();
528 return ttc;
529}
530
531
532template<class ReactionThermo, class ThermoType>
535{
536 auto tQdot = volScalarField::New
537 (
538 "Qdot",
540 this->mesh_,
542 );
543
544 if (this->chemistry_)
545 {
546 scalarField& Qdot = tQdot.ref();
547
548 forAll(Y_, i)
549 {
550 forAll(Qdot, celli)
551 {
552 const scalar hi = specieThermo_[i].Hc();
553 Qdot[celli] -= hi*RR_[i][celli];
554 }
555 }
556 tQdot.ref().correctBoundaryConditions();
557 }
559 return tQdot;
560}
561
562
563template<class ReactionThermo, class ThermoType>
566(
567 const label ri,
568 const label si
569) const
570{
571 scalar pf, cf, pr, cr;
572 label lRef, rRef;
573
575 (
576 "RR",
578 this->mesh(),
580 );
581 auto& RR = tRR.ref();
582
584 const scalarField& rho = trho();
585
586 const scalarField& T = this->thermo().T();
587 const scalarField& p = this->thermo().p();
588
589 forAll(rho, celli)
590 {
591 const scalar rhoi = rho[celli];
592 const scalar Ti = T[celli];
593 const scalar pi = p[celli];
594
595 for (label i=0; i<nSpecie_; i++)
596 {
597 const scalar Yi = Y_[i][celli];
598 c_[i] = rhoi*Yi/specieThermo_[i].W();
599 }
600
601 const scalar w = omegaI
602 (
603 ri,
604 c_,
605 Ti,
606 pi,
607 pf,
608 cf,
609 lRef,
610 pr,
611 cr,
612 rRef
613 );
614
615 RR[celli] = w*specieThermo_[si].W();
617
618 return tRR;
619}
620
621
622template<class ReactionThermo, class ThermoType>
624{
625 if (!this->chemistry_)
626 {
627 return;
628 }
629
631 const scalarField& rho = trho();
632
633 const scalarField& T = this->thermo().T();
634 const scalarField& p = this->thermo().p();
635
636 forAll(rho, celli)
637 {
638 const scalar rhoi = rho[celli];
639 const scalar Ti = T[celli];
640 const scalar pi = p[celli];
641
642 for (label i=0; i<nSpecie_; i++)
643 {
644 const scalar Yi = Y_[i][celli];
645 c_[i] = rhoi*Yi/specieThermo_[i].W();
646 }
647
648 omega(c_, Ti, pi, dcdt_);
649
650 for (label i=0; i<nSpecie_; i++)
651 {
652 RR_[i][celli] = dcdt_[i]*specieThermo_[i].W();
654 }
655}
656
657
658template<class ReactionThermo, class ThermoType>
659template<class DeltaTType>
660Foam::scalar Foam::StandardChemistryModel<ReactionThermo, ThermoType>::solve
661(
662 const DeltaTType& deltaT
663)
664{
666
667 scalar deltaTMin = GREAT;
668
669 if (!this->chemistry_)
670 {
671 return deltaTMin;
672 }
673
675 const scalarField& rho = trho();
676
677 const scalarField& T = this->thermo().T();
678 const scalarField& p = this->thermo().p();
679
680 scalarField c0(nSpecie_);
681
682 forAll(rho, celli)
683 {
684 scalar Ti = T[celli];
685
686 if (Ti > Treact_)
687 {
688 const scalar rhoi = rho[celli];
689 scalar pi = p[celli];
690
691 for (label i=0; i<nSpecie_; i++)
692 {
693 c_[i] = rhoi*Y_[i][celli]/specieThermo_[i].W();
694 c0[i] = c_[i];
695 }
696
697 // Initialise time progress
698 scalar timeLeft = deltaT[celli];
699
700 // Calculate the chemical source terms
701 while (timeLeft > SMALL)
702 {
703 scalar dt = timeLeft;
704 this->solve(c_, Ti, pi, dt, this->deltaTChem_[celli]);
705 timeLeft -= dt;
706 }
707
708 deltaTMin = min(this->deltaTChem_[celli], deltaTMin);
709
710 this->deltaTChem_[celli] =
711 min(this->deltaTChem_[celli], this->deltaTChemMax_);
712
713 for (label i=0; i<nSpecie_; i++)
714 {
715 RR_[i][celli] =
716 (c_[i] - c0[i])*specieThermo_[i].W()/deltaT[celli];
717 }
718 }
719 else
720 {
721 for (label i=0; i<nSpecie_; i++)
722 {
723 RR_[i][celli] = 0;
724 }
725 }
727
728 return deltaTMin;
729}
730
731
732template<class ReactionThermo, class ThermoType>
733Foam::scalar Foam::StandardChemistryModel<ReactionThermo, ThermoType>::solve
734(
735 const scalar deltaT
736)
737{
738 // Don't allow the time-step to change more than a factor of 2
739 return min
740 (
742 2*deltaT
743 );
744}
745
746
747template<class ReactionThermo, class ThermoType>
748Foam::scalar Foam::StandardChemistryModel<ReactionThermo, ThermoType>::solve
749(
750 const scalarField& deltaT
751)
752{
753 return this->solve<scalarField>(deltaT);
754}
755
756
757// ************************************************************************* //
scalar delta
constexpr scalar pi(M_PI)
#define R(A, B, C, D, E, F, K, M)
ReactionThermo & thermo()
Return access to the thermo package.
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
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())
DimensionedField< scalar, volMesh > Internal
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
ODESystem()
Construct null.
Definition ODESystem.H:53
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition Reaction.H:69
scalarField c_
Temporary concentration field.
virtual tmp< volScalarField > tc() const
Return the chemical time scale.
scalarField dcdt_
Temporary rate-of-change of concentration field.
virtual tmp< volScalarField::Internal > calculateRR(const label reactionI, const label speciei) const
Return reaction rate of the speciei in reactionI.
PtrList< volScalarField::Internal > & RR()
Write access to chemical source terms.
const PtrList< ThermoType > & specieThermo_
Thermodynamic data of the species.
PtrList< volScalarField > & Y_
Reference to the field of specie mass fractions.
virtual void omega(const scalarField &c, const scalar T, const scalar p, scalarField &dcdt) const
dc/dt = omega, rate of change in concentration, for each species
virtual void derivatives(const scalar t, const scalarField &c, scalarField &dcdt) const
Calculate the derivatives in dydx.
virtual ~StandardChemistryModel()
Destructor.
const PtrList< Reaction< ThermoType > > & reactions_
Reactions.
virtual scalar omegaI(label iReaction, const scalarField &c, const scalar T, const scalar p, scalar &pf, scalar &cf, label &lRef, scalar &pr, scalar &cr, label &rRef) const
Return the reaction rate for iReaction and the reference.
virtual label nReaction() const
The number of reactions.
scalar Treact_
Temperature below which the reaction rates are assumed 0.
PtrList< volScalarField::Internal > RR_
List of reaction rate per specie [kg/m3/s].
label nReaction_
Number of reactions.
virtual tmp< volScalarField > Qdot() const
Return the heat release rate [kg/m/s3].
virtual void jacobian(const scalar t, const scalarField &c, scalarField &dcdt, scalarSquareMatrix &dfdc) const
Calculate the Jacobian of the system.
virtual void calculate()
Calculates the reaction rates.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A class representing the concept of a uniform field which stores only the single value and providing ...
const fvMesh & mesh_
Reference to the mesh database.
void correct()
Correct function - updates due to mesh changes.
Switch chemistry_
Chemistry activation switch.
const fvMesh & mesh() const
Return const access to the mesh database.
const scalar deltaTChemMax_
Maximum chemical time step.
volScalarField::Internal deltaTChem_
Latest estimation of integration step.
virtual const volScalarField & T() const
Temperature [K].
virtual volScalarField & p()
Pressure [Pa].
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Foam::reactingMixture.
A class for managing temporary objects.
Definition tmp.H:75
static const word null
An empty word.
Definition word.H:84
basicSpecieMixture & composition
volScalarField & p
PtrList< volScalarField > & Y
const volScalarField & T
dynamicFvMesh & mesh
auto & name
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))
word timeName
Definition getTimeIndex.H:3
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)
const dimensionSet dimEnergy
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
SquareMatrix< scalar > scalarSquareMatrix
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition POSIX.C:1065
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
tmp< volScalarField > trho
const volScalarField & cp
CEqn solve()
scalar Qdot
psiReactionThermo & thermo
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299