standardChemistryModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "standardChemistryModel.H"
27 #include "UniformField.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class ThermoType>
34 (
35  const fluidReactionThermo& thermo
36 )
37 :
38  basicChemistryModel(thermo),
39  ODESystem(),
40  Y_(this->thermo().composition().Y()),
41  mixture_(refCast<const multiComponentMixture<ThermoType>>(this->thermo())),
42  specieThermos_(mixture_.specieThermos()),
43  reactions_(mixture_.species(), specieThermos_, this->mesh(), *this),
44  nSpecie_(Y_.size()),
45  nReaction_(reactions_.size()),
46  Treact_(basicChemistryModel::template lookupOrDefault<scalar>("Treact", 0)),
47  RR_(nSpecie_),
48  c_(nSpecie_),
49  dcdt_(nSpecie_)
50 {
51  // Create the fields for the chemistry sources
52  forAll(RR_, fieldi)
53  {
54  RR_.set
55  (
56  fieldi,
58  (
59  IOobject
60  (
61  "RR." + Y_[fieldi].name(),
62  this->mesh().time().timeName(),
63  this->mesh(),
64  IOobject::NO_READ,
65  IOobject::NO_WRITE
66  ),
67  thermo.T().mesh(),
69  )
70  );
71  }
72 
73  Info<< "standardChemistryModel: Number of species = " << nSpecie_
74  << " and reactions = " << nReaction_ << endl;
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
79 
80 template<class ThermoType>
82 {}
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
87 template<class ThermoType>
89 (
90  const scalar p,
91  const scalar T,
92  const scalarField& c,
93  const label li,
94  scalarField& dcdt
95 ) const
96 {
97  dcdt = Zero;
98 
99  forAll(reactions_, i)
100  {
101  const Reaction<ThermoType>& R = reactions_[i];
102 
103  R.omega(p, T, c, li, dcdt);
104  }
105 }
106 
107 
108 template<class ThermoType>
110 (
111  const scalar t,
112  const scalarField& c,
113  const label li,
114  scalarField& dcdt
115 ) const
116 {
117  const scalar T = c[nSpecie_];
118  const scalar p = c[nSpecie_ + 1];
119 
120  forAll(c_, i)
121  {
122  c_[i] = max(c[i], 0);
123  }
124 
125  dcdt = Zero;
126 
127  // Evaluate contributions from reactions
128  omega(p, T, c_, li, dcdt);
129 
130  // Evaluate the effect on the thermodynamic system ...
131 
132  // c*Cp
133  scalar ccp = 0;
134  for (label i = 0; i < nSpecie_; i++)
135  {
136  ccp += c_[i]*specieThermos_[i].cp(p, T);
137  }
138 
139  // dT/dt
140  scalar& dTdt = dcdt[nSpecie_];
141  for (label i = 0; i < nSpecie_; i++)
142  {
143  dTdt -= dcdt[i]*specieThermos_[i].ha(p, T);
144  }
145  dTdt /= ccp;
146 
147  // dp/dt = 0 (pressure is assumed constant)
148 }
149 
150 
151 template<class ThermoType>
153 (
154  const scalar t,
155  const scalarField& c,
156  const label li,
157  scalarField& dcdt,
159 ) const
160 {
161  const scalar T = c[nSpecie_];
162  const scalar p = c[nSpecie_ + 1];
163 
164  forAll(c_, i)
165  {
166  c_[i] = max(c[i], 0);
167  }
168 
169  dcdt = Zero;
170  J = Zero;
171 
172  // Evaluate contributions from reactions
173  forAll(reactions_, ri)
174  {
175  const Reaction<ThermoType>& R = reactions_[ri];
176  scalar omegaI, kfwd, kbwd;
177  const labelList null;
178  R.dwdc(p, T, c_, li, J, dcdt, omegaI, kfwd, kbwd, false, null);
179  R.dwdT(p, T, c_, li, omegaI, kfwd, kbwd, J, false, null, nSpecie_);
180  }
181 
182  // Evaluate the effect on the thermodynamic system ...
183 
184  // c*Cp
185  scalar ccp = 0, dccpdT = 0;
186  for (label i = 0; i < nSpecie_; i++)
187  {
188  ccp += c_[i]*specieThermos_[i].cp(p, T);
189  dccpdT += c_[i]*specieThermos_[i].dcpdT(p, T);
190  }
191 
192  // dT/dt
193  scalar& dTdt = dcdt[nSpecie_];
194  for (label i = 0; i < nSpecie_; i++)
195  {
196  dTdt -= dcdt[i]*specieThermos_[i].ha(p, T);
197  }
198  dTdt /= ccp;
199 
200  // dp/dt = 0 (pressure is assumed constant)
201 
202  // d(dTdt)/dc
203  for (label i = 0; i < nSpecie_; i++)
204  {
205  scalar& d2Tdtdci = J(nSpecie_, i);
206  for (label j = 0; j < nSpecie_; j++)
207  {
208  const scalar d2cjdtdci = J(j, i);
209  d2Tdtdci -= d2cjdtdci*specieThermos_[j].ha(p, T);
210  }
211  d2Tdtdci -= specieThermos_[i].cp(p, T)*dTdt;
212  d2Tdtdci /= ccp;
213  }
214 
215  // d(dTdt)/dT
216  scalar& d2TdtdT = J(nSpecie_, nSpecie_);
217  for (label i = 0; i < nSpecie_; i++)
218  {
219  const scalar d2cidtdT = J(i, nSpecie_);
220  d2TdtdT -=
221  dcdt[i]*specieThermos_[i].cp(p, T)
222  + d2cidtdT*specieThermos_[i].ha(p, T);
223  }
224  d2TdtdT -= dTdt*dccpdT;
225  d2TdtdT /= ccp;
226 
227  // d(dpdt)/dc = 0 (pressure is assumed constant)
228 
229  // d(dpdt)/dT = 0 (pressure is assumed constant)
230 }
231 
232 
233 template<class ThermoType>
236 {
238  (
240  (
241  "tc",
242  this->mesh(),
243  dimensionedScalar(dimTime, small),
244  extrapolatedCalculatedFvPatchScalarField::typeName
245  )
246  );
247  scalarField& tc = ttc.ref();
248 
249  tmp<volScalarField> trho(this->thermo().rho());
250  const scalarField& rho = trho();
251 
252  const scalarField& T = this->thermo().T();
253  const scalarField& p = this->thermo().p();
254 
255  if (this->chemistry_)
256  {
257  reactionEvaluationScope scope(*this);
258 
259  forAll(rho, celli)
260  {
261  const scalar rhoi = rho[celli];
262  const scalar Ti = T[celli];
263  const scalar pi = p[celli];
264 
265  for (label i=0; i<nSpecie_; i++)
266  {
267  c_[i] = rhoi*Y_[i][celli]/specieThermos_[i].W();
268  }
269 
270  // A reaction's rate scale is calculated as it's molar
271  // production rate divided by the total number of moles in the
272  // system.
273  //
274  // The system rate scale is the average of the reactions' rate
275  // scales weighted by the reactions' molar production rates. This
276  // weighting ensures that dominant reactions provide the largest
277  // contribution to the system rate scale.
278  //
279  // The system time scale is then the reciprocal of the system rate
280  // scale.
281  //
282  // Contributions from forward and reverse reaction rates are
283  // handled independently and identically so that reversible
284  // reactions produce the same result as the equivalent pair of
285  // irreversible reactions.
286 
287  scalar sumW = 0, sumWRateByCTot = 0;
288  forAll(reactions_, i)
289  {
290  const Reaction<ThermoType>& R = reactions_[i];
291  scalar omegaf, omegar;
292  R.omega(pi, Ti, c_, celli, omegaf, omegar);
293 
294  scalar wf = 0;
295  forAll(R.rhs(), s)
296  {
297  wf += R.rhs()[s].stoichCoeff*omegaf;
298  }
299  sumW += wf;
300  sumWRateByCTot += sqr(wf);
301 
302  scalar wr = 0;
303  forAll(R.lhs(), s)
304  {
305  wr += R.lhs()[s].stoichCoeff*omegar;
306  }
307  sumW += wr;
308  sumWRateByCTot += sqr(wr);
309  }
310 
311  tc[celli] =
312  sumWRateByCTot == 0 ? vGreat : sumW/sumWRateByCTot*sum(c_);
313  }
314  }
315 
316  ttc.ref().correctBoundaryConditions();
317 
318  return ttc;
319 }
320 
321 
322 template<class ThermoType>
325 {
326  tmp<volScalarField> tQdot
327  (
329  (
330  "Qdot",
331  this->mesh_,
333  )
334  );
335 
336  if (this->chemistry_)
337  {
338  reactionEvaluationScope scope(*this);
339 
340  scalarField& Qdot = tQdot.ref();
341 
342  forAll(Y_, i)
343  {
344  forAll(Qdot, celli)
345  {
346  const scalar hi = specieThermos_[i].Hf();
347  Qdot[celli] -= hi*RR_[i][celli];
348  }
349  }
350  }
351 
352  return tQdot;
353 }
354 
355 
356 template<class ThermoType>
359 (
360  const label ri,
361  const label si
362 ) const
363 {
365  (
367  (
368  "RR",
369  this->mesh(),
371  )
372  );
373 
375 
376  tmp<volScalarField> trho(this->thermo().rho());
377  const scalarField& rho = trho();
378 
379  const scalarField& T = this->thermo().T();
380  const scalarField& p = this->thermo().p();
381 
382  reactionEvaluationScope scope(*this);
383 
384  scalar omegaf, omegar;
385 
386  forAll(rho, celli)
387  {
388  const scalar rhoi = rho[celli];
389  const scalar Ti = T[celli];
390  const scalar pi = p[celli];
391 
392  for (label i=0; i<nSpecie_; i++)
393  {
394  const scalar Yi = Y_[i][celli];
395  c_[i] = rhoi*Yi/specieThermos_[i].W();
396  }
397 
398  const Reaction<ThermoType>& R = reactions_[ri];
399  const scalar omegaI = R.omega(pi, Ti, c_, celli, omegaf, omegar);
400 
401  forAll(R.lhs(), s)
402  {
403  if (si == R.lhs()[s].index)
404  {
405  RR[celli] -= R.lhs()[s].stoichCoeff*omegaI;
406  }
407  }
408 
409  forAll(R.rhs(), s)
410  {
411  if (si == R.rhs()[s].index)
412  {
413  RR[celli] += R.rhs()[s].stoichCoeff*omegaI;
414  }
415  }
416 
417  RR[celli] *= specieThermos_[si].W();
418  }
419 
420  return tRR;
421 }
422 
423 
424 template<class ThermoType>
426 {
427  if (!this->chemistry_)
428  {
429  return;
430  }
431 
432  tmp<volScalarField> trho(this->thermo().rho());
433  const scalarField& rho = trho();
434 
435  const scalarField& T = this->thermo().T();
436  const scalarField& p = this->thermo().p();
437 
438  reactionEvaluationScope scope(*this);
439 
440  forAll(rho, celli)
441  {
442  const scalar rhoi = rho[celli];
443  const scalar Ti = T[celli];
444  const scalar pi = p[celli];
445 
446  for (label i=0; i<nSpecie_; i++)
447  {
448  const scalar Yi = Y_[i][celli];
449  c_[i] = rhoi*Yi/specieThermos_[i].W();
450  }
451 
452  omega(pi, Ti, c_, celli, dcdt_);
453 
454  for (label i=0; i<nSpecie_; i++)
455  {
456  RR_[i][celli] = dcdt_[i]*specieThermos_[i].W();
457  }
458  }
459 }
460 
461 
462 template<class ThermoType>
463 template<class DeltaTType>
465 (
466  const DeltaTType& deltaT
467 )
468 {
470 
471  scalar deltaTMin = great;
472 
473  if (!this->chemistry_)
474  {
475  return deltaTMin;
476  }
477 
478  tmp<volScalarField> trho(this->thermo().rho0());
479  const scalarField& rho = trho();
480 
481  const scalarField& T = this->thermo().T().oldTime();
482  const scalarField& p = this->thermo().p().oldTime();
483 
484  reactionEvaluationScope scope(*this);
485 
486  scalarField c0(nSpecie_);
487 
488  forAll(rho, celli)
489  {
490  scalar Ti = T[celli];
491 
492  if (Ti > Treact_)
493  {
494  const scalar rhoi = rho[celli];
495  scalar pi = p[celli];
496 
497  for (label i=0; i<nSpecie_; i++)
498  {
499  c_[i] = rhoi*Y_[i].oldTime()[celli]/specieThermos_[i].W();
500  c0[i] = c_[i];
501  }
502 
503  // Initialise time progress
504  scalar timeLeft = deltaT[celli];
505 
506  // Calculate the chemical source terms
507  while (timeLeft > small)
508  {
509  scalar dt = timeLeft;
510  this->solve(pi, Ti, c_, celli, dt, this->deltaTChem_[celli]);
511  timeLeft -= dt;
512  }
513 
514  deltaTMin = min(this->deltaTChem_[celli], deltaTMin);
515 
516  this->deltaTChem_[celli] =
517  min(this->deltaTChem_[celli], this->deltaTChemMax_);
518 
519  for (label i=0; i<nSpecie_; i++)
520  {
521  RR_[i][celli] =
522  (c_[i] - c0[i])*specieThermos_[i].W()/deltaT[celli];
523  }
524  }
525  else
526  {
527  for (label i=0; i<nSpecie_; i++)
528  {
529  RR_[i][celli] = 0;
530  }
531  }
532  }
533 
534  return deltaTMin;
535 }
536 
537 
538 template<class ThermoType>
540 (
541  const scalar deltaT
542 )
543 {
544  // Don't allow the time-step to change more than a factor of 2
545  return min
546  (
548  2*deltaT
549  );
550 }
551 
552 
553 template<class ThermoType>
555 (
556  const scalarField& deltaT
557 )
558 {
559  return this->solve<scalarField>(deltaT);
560 }
561 
562 
563 // ************************************************************************* //
A class representing the concept of a uniform field which stores only the single value and providing ...
Definition: UniformField.H:47
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
fluidReactionThermo & thermo
Definition: createFields.H:28
virtual tmp< volScalarField > Qdot() const
Return the heat release rate [kg/m/s^3].
Abstract base class for the systems of ordinary differential equations.
Definition: ODESystem.H:46
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
basicSpecieMixture & composition
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
void dwdc(const scalar p, const scalar T, const scalarField &c, const label li, scalarSquareMatrix &J, scalarField &dcdt, scalar &omegaI, scalar &kfwd, scalar &kbwd, const bool reduced, const List< label > &c2s) const
Derivative of the net reaction rate for each species involved.
Definition: Reaction.C:352
dimensionedSymmTensor sqr(const dimensionedVector &dv)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
virtual void jacobian(const scalar t, const scalarField &c, const label li, scalarField &dcdt, scalarSquareMatrix &J) const
Calculate the Jacobian of the system.
virtual void omega(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &dcdt) const
dc/dt = omega, rate of change in concentration, for each species
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< volScalarField > trho
Base-class for multi-component fluid thermodynamic properties.
scalar rho0
const List< specieCoeffs > & lhs() const
Return the components of the left hand side.
Definition: reactionI.H:36
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:55
scalar Qdot
Definition: solveChemistry.H:2
void omega(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &dcdt) const
Net reaction rate for individual species.
Definition: Reaction.C:290
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const dimensionSet dimTime
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
virtual void calculate()
Calculates the reaction rates.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
dynamicFvMesh & mesh
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-9/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
Foam::multiComponentMixture.
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
standardChemistryModel(const fluidReactionThermo &thermo)
Construct from thermo.
Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equa...
const Mesh & mesh() const
Return mesh.
virtual void derivatives(const scalar t, const scalarField &c, const label li, scalarField &dcdt) const
Calculate the derivatives in dydx.
const dimensionSet dimEnergy
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual tmp< volScalarField::Internal > calculateRR(const label reactionI, const label speciei) const
Return reaction rate of the speciei in reactionI.
virtual tmp< volScalarField > tc() const
Return the chemical time scale.
const List< specieCoeffs > & rhs() const
Return the components of the right hand side.
Definition: reactionI.H:42
rhoEqn solve()
#define R(A, B, C, D, E, F, K, M)
virtual const volScalarField & T() const =0
Temperature [K].
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const scalarList W(::W(thermo))
Class to define scope of reaction evaluation. Runs pre-evaluate.
PtrList< volScalarField > & Y
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimVolume
messageStream Info
void dwdT(const scalar p, const scalar T, const scalarField &c, const label li, const scalar omegaI, const scalar kfwd, const scalar kbwd, scalarSquareMatrix &J, const bool reduced, const List< label > &c2s, const label indexT) const
Derivative of the net reaction rate for each species involved.
Definition: Reaction.C:499
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual ~standardChemistryModel()
Destructor.
Base class for chemistry models.