EulerImplicit.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 -------------------------------------------------------------------------------
10 License
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 
28 #include "EulerImplicit.H"
30 #include "simpleMatrix.H"
31 #include "Reaction.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class ChemistryModel>
37 (
38  typename ChemistryModel::reactionThermo& thermo
39 )
40 :
42  coeffsDict_(this->subDict("EulerImplicitCoeffs")),
43  cTauChem_(coeffsDict_.get<scalar>("cTauChem")),
44  eqRateLimiter_(coeffsDict_.lookup("equilibriumRateLimiter")),
45  cTp_(this->nEqns())
46 {}
47 
48 
49 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
50 
51 template<class ChemistryModel>
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
57 
58 template<class ChemistryModel>
60 (
61  const label index,
62  const scalar pr,
63  const scalar pf,
64  const scalar corr,
65  const label lRef,
66  const label rRef,
67  const scalar p,
68  const scalar T,
70 ) const
71 {
73  this->reactions_[index];
74 
75  forAll(R.lhs(), s)
76  {
77  const label si = R.lhs()[s].index;
78  const scalar sl = R.lhs()[s].stoichCoeff;
79  RR[si][rRef] -= sl*pr*corr;
80  RR[si][lRef] += sl*pf*corr;
81  }
82 
83  forAll(R.rhs(), s)
84  {
85  const label si = R.rhs()[s].index;
86  const scalar sr = R.rhs()[s].stoichCoeff;
87  RR[si][lRef] -= sr*pf*corr;
88  RR[si][rRef] += sr*pr*corr;
89  }
90 }
91 
92 
93 template<class ChemistryModel>
95 (
96  scalarField& c,
97  scalar& T,
98  scalar& p,
99  scalar& deltaT,
100  scalar& subDeltaT
101 ) const
102 {
103  const label nSpecie = this->nSpecie();
105 
106  for (label i=0; i<nSpecie; i++)
107  {
108  c[i] = max(0, c[i]);
109  }
110 
111  // Calculate the absolute enthalpy
112  const scalar cTot = sum(c);
113  typename ChemistryModel::thermoType mixture
114  (
115  (this->specieThermo_[0].W()*c[0])*this->specieThermo_[0]
116  );
117  for (label i=1; i<nSpecie; i++)
118  {
119  mixture += (this->specieThermo_[i].W()*c[i])*this->specieThermo_[i];
120  }
121  const scalar ha = mixture.Ha(p, T);
122  const scalar deltaTEst = min(deltaT, subDeltaT);
123 
124  forAll(this->reactions(), i)
125  {
126  scalar pf, cf, pr, cr;
127  label lRef, rRef;
128 
129  const scalar omegai =
130  this->omegaI(i, c, T, p, pf, cf, lRef, pr, cr, rRef);
131 
132  scalar corr = 1;
133  if (eqRateLimiter_)
134  {
135  if (omegai < 0)
136  {
137  corr = 1/(1 + pr*deltaTEst);
138  }
139  else
140  {
141  corr = 1/(1 + pf*deltaTEst);
142  }
143  }
144 
145  updateRRInReactionI(i, pr, pf, corr, lRef, rRef, p, T, RR);
146  }
147 
148  // Calculate the stable/accurate time-step
149  scalar tMin = GREAT;
150 
151  for (label i=0; i<nSpecie; i++)
152  {
153  scalar d = 0;
154  for (label j=0; j<nSpecie; j++)
155  {
156  d -= RR(i, j)*c[j];
157  }
158 
159  if (d < -SMALL)
160  {
161  tMin = min(tMin, -(c[i] + SMALL)/d);
162  }
163  else
164  {
165  d = max(d, SMALL);
166  const scalar cm = max(cTot - c[i], 1e-5);
167  tMin = min(tMin, cm/d);
168  }
169  }
170 
171  subDeltaT = cTauChem_*tMin;
172  deltaT = min(deltaT, subDeltaT);
173 
174  // Add the diagonal and source contributions from the time-derivative
175  for (label i=0; i<nSpecie; i++)
176  {
177  RR(i, i) += 1/deltaT;
178  RR.source()[i] = c[i]/deltaT;
179  }
180 
181  // Solve for the new composition
182  c = RR.LUsolve();
183 
184  // Limit the composition
185  for (label i=0; i<nSpecie; i++)
186  {
187  c[i] = max(0, c[i]);
188  }
189 
190  // Update the temperature
191  mixture = (this->specieThermo_[0].W()*c[0])*this->specieThermo_[0];
192  for (label i=1; i<nSpecie; i++)
193  {
194  mixture += (this->specieThermo_[i].W()*c[i])*this->specieThermo_[i];
195  }
196  T = mixture.THa(ha, p, T);
197 }
198 
199 
200 // ************************************************************************* //
Foam::constant::thermodynamic::RR
const scalar RR
Universal gas constant: default in [J/(kmol K)].
Definition: thermodynamicConstants.C:46
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::simpleMatrix
A simple square matrix solver with scalar coefficients.
Definition: simpleMatrix.H:49
s
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))
Definition: gmvOutputSpray.H:25
Foam::chemistrySolver
An abstract base class for solving chemistry.
Definition: chemistrySolver.H:52
nSpecie
label nSpecie
Definition: readInitialConditions.H:17
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::EulerImplicit::EulerImplicit
EulerImplicit(typename ChemistryModel::reactionThermo &thermo)
Construct from thermo.
Definition: EulerImplicit.C:37
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::EulerImplicit::solve
virtual void solve(scalarField &c, scalar &T, scalar &p, scalar &deltaT, scalar &subDeltaT) const
Update the concentrations and return the chemical time.
Definition: EulerImplicit.C:95
R
#define R(A, B, C, D, E, F, K, M)
Foam::label
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:62
Foam::Field< scalar >
EulerImplicit.H
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Reaction.H
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::EulerImplicit::~EulerImplicit
virtual ~EulerImplicit()
Destructor.
Definition: EulerImplicit.C:52
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
simpleMatrix.H
Foam::Reaction
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:56
mixture
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
Foam::EulerImplicit::updateRRInReactionI
void updateRRInReactionI(const label index, const scalar pr, const scalar pf, const scalar corr, const label lRef, const label rRef, const scalar p, const scalar T, simpleMatrix< scalar > &RR) const
Definition: EulerImplicit.C:60
Foam::mixture
Definition: mixture.H:54