EulerImplicit.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-2024 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 "EulerImplicit.H"
27 #include "SubField.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class ChemistryModel>
34 (
36 )
37 :
38  chemistrySolver<ChemistryModel>(thermo),
39  cTauChem_
40  (
41  this->subDict("EulerImplicitCoeffs").template lookup<scalar>("cTauChem")
42  ),
43  cTp_(this->nEqns()),
44  R_(this->nEqns()),
45  J_(this->nEqns()),
46  E_(this->nEqns() - 2)
47 {}
48 
49 
50 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
51 
52 template<class ChemistryModel>
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
58 
59 template<class ChemistryModel>
61 (
62  scalar& p,
63  scalar& T,
64  scalarField& c,
65  const label li,
66  scalar& deltaT,
67  scalar& subDeltaT
68 ) const
69 {
70  const label nSpecie = this->nSpecie();
71 
72  // Map the composition, temperature and pressure into cTp
73  for (int i=0; i<nSpecie; i++)
74  {
75  cTp_[i] = max(0, c[i]);
76  }
77  cTp_[nSpecie] = T;
78  cTp_[nSpecie + 1] = p;
79 
80  // Calculate the reaction rate and Jacobian
81  this->jacobian(0, cTp_, li, R_, J_);
82 
83  // Calculate the stable/accurate time-step
84  scalar tMin = great;
85  const scalar cTot = sum(c);
86 
87  for (label i=0; i<nSpecie; i++)
88  {
89  if (R_[i] < -small)
90  {
91  tMin = min(tMin, -(cTp_[i] + small)/R_[i]);
92  }
93  else
94  {
95  tMin = min
96  (
97  tMin,
98  max(cTot - cTp_[i], 1e-5)/max(R_[i], small)
99  );
100  }
101  }
102 
103  subDeltaT = cTauChem_*tMin;
104  deltaT = min(deltaT, subDeltaT);
105 
106  // Assemble the Euler implicit matrix for the composition
107  scalarField& source = E_.source();
108  for (label i=0; i<nSpecie; i++)
109  {
110  E_(i, i) = 1/deltaT - J_(i, i);
111  source[i] = R_[i] + E_(i, i)*cTp_[i];
112 
113  for (label j=0; j<nSpecie; j++)
114  {
115  if (i != j)
116  {
117  E_(i, j) = -J_(i, j);
118  source[i] += E_(i, j)*cTp_[j];
119  }
120  }
121  }
122 
123  // Solve for the new composition
124  scalarField::subField(cTp_, nSpecie) = E_.LUsolve();
125 
126  // Limit the composition and transfer back into c
127  for (label i=0; i<nSpecie; i++)
128  {
129  c[i] = max(0, cTp_[i]);
130  }
131 
132  // Euler explicit integrate the temperature.
133  // Separating the integration of temperature from composition
134  // is significantly more stable for exothermic systems
135  T += deltaT*R_[nSpecie];
136 }
137 
138 
139 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
virtual void solve(scalar &p, scalar &T, scalarField &c, const label li, scalar &deltaT, scalar &subDeltaT) const
Update the concentrations and return the chemical time.
Definition: EulerImplicit.C:61
EulerImplicit(const fluidMulticomponentThermo &thermo)
Construct from thermo.
Definition: EulerImplicit.C:34
virtual ~EulerImplicit()
Destructor.
Definition: EulerImplicit.C:53
SubField< scalar > subField
Declare type of subField.
Definition: Field.H:101
An abstract base class for solving chemistry.
Base-class for multi-component fluid thermodynamic properties.
const dimensionedScalar c
Speed of light in a vacuum.
const doubleScalar e
Definition: doubleScalar.H:106
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
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const label nSpecie
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31