EulerImplicit.H
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-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::EulerImplicit
28
29Description
30 An Euler implicit solver for chemistry
31
32SourceFiles
33 EulerImplicit.C
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef EulerImplicit_H
38#define EulerImplicit_H
39
40#include "chemistrySolver.H"
41#include "Switch.H"
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45namespace Foam
46{
47
48template <class Type>
49class simpleMatrix;
50
51/*---------------------------------------------------------------------------*\
52 Class EulerImplicit Declaration
53\*---------------------------------------------------------------------------*/
54
55template<class ChemistryModel>
56class EulerImplicit
57:
58 public chemistrySolver<ChemistryModel>
59{
60 // Private data
61
62 //- Coefficients dictionary
63 dictionary coeffsDict_;
64
65
66 // Model constants
67
68 //- Chemistry timescale
69 scalar cTauChem_;
70
71 //- Equilibrium rate limiter flag (on/off)
72 Switch eqRateLimiter_;
73
74 // Solver data
75 mutable scalarField cTp_;
76
77
78public:
79
80 //- Runtime type information
81 TypeName("EulerImplicit");
82
83
84 // Constructors
85
86 //- Construct from thermo
87 EulerImplicit(typename ChemistryModel::reactionThermo& thermo);
88
89
90 //- Destructor
91 virtual ~EulerImplicit();
92
93
94 // Member Functions
95
97 (
98 const label index,
99 const scalar pr,
100 const scalar pf,
101 const scalar corr,
102 const label lRef,
103 const label rRef,
104 const scalar p,
105 const scalar T,
107 ) const;
108
109 //- Update the concentrations and return the chemical time
110 virtual void solve
111 (
112 scalarField& c,
113 scalar& T,
114 scalar& p,
115 scalar& deltaT,
116 scalar& subDeltaT
117 ) const;
118};
119
120
121// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
122
123} // End namespace Foam
124
125// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
126
127#ifdef NoRepository
128 #include "EulerImplicit.C"
129#endif
130
131// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
132
133#endif
134
135// ************************************************************************* //
An Euler implicit solver for chemistry.
Definition: EulerImplicit.H:58
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:61
TypeName("EulerImplicit")
Runtime type information.
virtual ~EulerImplicit()
Destructor.
Definition: EulerImplicit.C:53
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
An abstract base class for solving chemistry.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A simple square matrix solver with scalar coefficients.
Definition: simpleMatrix.H:67
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
volScalarField & p
const volScalarField & T
Namespace for OpenFOAM.
CEqn solve()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73