EulerImplicit.H
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-2020 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 Class
25  Foam::EulerImplicit
26 
27 Description
28  An Euler implicit solver for chemistry
29 
30  Euler implicit integration based on the reaction rate Jacobian is used to
31  solver for the composition and Euler explicit integration for the
32  temperature. Separating the integration of temperature from composition
33  in this manner is significantly more stable for exothermic systems
34 
35 SourceFiles
36  EulerImplicit.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef EulerImplicit_H
41 #define EulerImplicit_H
42 
43 #include "chemistrySolver.H"
44 #include "simpleMatrix.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class EulerImplicit Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 template<class ChemistryModel>
56 class EulerImplicit
57 :
58  public chemistrySolver<ChemistryModel>
59 {
60  // Private Data
61 
62  //- Coefficients dictionary
63  dictionary coeffsDict_;
64 
65  //- Chemistry timescale coefficient
66  scalar cTauChem_;
67 
68  //- Field encapsulating the composition, temperature and pressure
69  mutable scalarField cTp_;
70 
71  //- Reaction rate field
72  mutable scalarField R_;
73 
74  //- Reaction Jacobian
75  mutable scalarSquareMatrix J_;
76 
77  //- Euler implicit integration matrix for composition
78  mutable simpleMatrix<scalar> E_;
79 
80 
81 public:
82 
83  //- Runtime type information
84  TypeName("EulerImplicit");
85 
86 
87  // Constructors
88 
89  //- Construct from thermo
91 
92 
93  //- Destructor
94  virtual ~EulerImplicit();
95 
96 
97  // Member Functions
98 
99  //- Update the concentrations and return the chemical time
100  virtual void solve
101  (
102  scalar& p,
103  scalar& T,
104  scalarField& c,
105  const label li,
106  scalar& deltaT,
107  scalar& subDeltaT
108  ) const;
109 };
110 
111 
112 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
113 
114 } // End namespace Foam
115 
116 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
117 
118 #ifdef NoRepository
119  #include "EulerImplicit.C"
120 #endif
121 
122 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
123 
124 #endif
125 
126 // ************************************************************************* //
An abstract base class for solving chemistry.
fluidReactionThermo & thermo
Definition: createFields.H:28
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Base-class for multi-component fluid thermodynamic properties.
TypeName("EulerImplicit")
Runtime type information.
const dimensionedScalar c
Speed of light in a vacuum.
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:59
EulerImplicit(const fluidReactionThermo &thermo)
Construct from thermo.
Definition: EulerImplicit.C:34
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
An Euler implicit solver for chemistry.
Definition: EulerImplicit.H:55
volScalarField & p
virtual ~EulerImplicit()
Destructor.
Definition: EulerImplicit.C:51
Namespace for OpenFOAM.