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-2019 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 SourceFiles
31  EulerImplicit.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef EulerImplicit_H
36 #define EulerImplicit_H
37 
38 #include "chemistrySolver.H"
39 #include "Switch.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 template <class Type>
47 class simpleMatrix;
48 
49 /*---------------------------------------------------------------------------*\
50  Class EulerImplicit Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 template<class ChemistryModel>
54 class EulerImplicit
55 :
56  public chemistrySolver<ChemistryModel>
57 {
58  // Private Data
59 
60  //- Coefficients dictionary
61  dictionary coeffsDict_;
62 
63 
64  // Model constants
65 
66  //- Chemistry timescale
67  scalar cTauChem_;
68 
69  //- Equilibrium rate limiter flag (on/off)
70  Switch eqRateLimiter_;
71 
72  // Solver data
73  mutable scalarField cTp_;
74 
75 
76 public:
77 
78  //- Runtime type information
79  TypeName("EulerImplicit");
80 
81 
82  // Constructors
83 
84  //- Construct from thermo
85  EulerImplicit(typename ChemistryModel::reactionThermo& thermo);
86 
87 
88  //- Destructor
89  virtual ~EulerImplicit();
90 
91 
92  // Member Functions
93 
95  (
96  const label index,
97  const scalar pr,
98  const scalar pf,
99  const scalar corr,
100  const label lRef,
101  const label rRef,
102  const scalar p,
103  const scalar T,
105  ) const;
106 
107  //- Update the concentrations and return the chemical time
108  virtual void solve
109  (
110  scalarField& c,
111  scalar& T,
112  scalar& p,
113  scalar& deltaT,
114  scalar& subDeltaT
115  ) const;
116 };
117 
118 
119 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
120 
121 } // End namespace Foam
122 
123 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
124 
125 #ifdef NoRepository
126  #include "EulerImplicit.C"
127 #endif
128 
129 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130 
131 #endif
132 
133 // ************************************************************************* //
An abstract base class for solving chemistry.
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
TypeName("EulerImplicit")
Runtime type information.
rhoReactionThermo & thermo
Definition: createFields.H:28
A simple square matrix solver with scalar coefficients.
Definition: simpleMatrix.H:47
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:93
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
EulerImplicit(typename ChemistryModel::reactionThermo &thermo)
Construct from thermo.
Definition: EulerImplicit.C:35
const dimensionedScalar c
Speed of light in a vacuum.
const scalar RR
Universal gas constant (default in [J/kmol/K])
An Euler implicit solver for chemistry.
Definition: EulerImplicit.H:53
volScalarField & p
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:58
virtual ~EulerImplicit()
Destructor.
Definition: EulerImplicit.C:50
Namespace for OpenFOAM.