JanevReactionRate.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-2021 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::JanevReactionRate
26 
27 Description
28  Janev, Langer, Evans and Post reaction rate.
29 
30 SourceFiles
31  JanevReactionRateI.H
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef JanevReactionRate_H
36 #define JanevReactionRate_H
37 
38 #include "scalarField.H"
39 #include "typeInfo.H"
40 #include "FixedList.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 // Forward declaration of friend functions and operators
48 
49 class JanevReactionRate;
50 
51 Ostream& operator<<(Ostream&, const JanevReactionRate&);
52 
53 
54 /*---------------------------------------------------------------------------*\
55  Class JanevReactionRate Declaration
56 \*---------------------------------------------------------------------------*/
57 
59 {
60  // Private Data
61 
62  scalar A_;
63  scalar beta_;
64  scalar Ta_;
65 
66  static const label nb_ = 9;
68 
69 
70 public:
71 
72  // Constructors
73 
74  //- Construct from components
75  inline JanevReactionRate
76  (
77  const scalar A,
78  const scalar beta,
79  const scalar Ta,
81  );
82 
83  //- Construct from dictionary
84  inline JanevReactionRate
85  (
86  const speciesTable& species,
87  const dictionary& dict
88  );
89 
90 
91  // Member Functions
92 
93  //- Return the type name
94  static word type()
95  {
96  return "Janev";
97  }
98 
99  //- Pre-evaluation hook
100  inline void preEvaluate() const;
101 
102  //- Post-evaluation hook
103  inline void postEvaluate() const;
104 
105  inline scalar operator()
106  (
107  const scalar p,
108  const scalar T,
109  const scalarField& c,
110  const label li
111  ) const;
112 
113  inline scalar ddT
114  (
115  const scalar p,
116  const scalar T,
117  const scalarField& c,
118  const label li
119  ) const;
120 
121  //- Third-body efficiencies (beta = 1-alpha)
122  // non-empty only for third-body reactions
123  // with enhanced molecularity (alpha != 1)
124  inline const List<Tuple2<label, scalar>>& beta() const;
125 
126  //- Species concentration derivative of the pressure dependent term
127  inline void dcidc
128  (
129  const scalar p,
130  const scalar T,
131  const scalarField& c,
132  const label li,
134  ) const;
135 
136  //- Temperature derivative of the pressure dependent term
137  inline scalar dcidT
138  (
139  const scalar p,
140  const scalar T,
141  const scalarField& c,
142  const label li
143  ) const;
144 
145  //- Write to stream
146  inline void write(Ostream& os) const;
147 
148 
149  // Ostream Operator
150 
151  inline friend Ostream& operator<<
152  (
153  Ostream&,
154  const JanevReactionRate&
155  );
156 };
157 
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
161 } // End namespace Foam
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 #include "JanevReactionRateI.H"
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 #endif
170 
171 // ************************************************************************* //
dictionary dict
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 postEvaluate() const
Post-evaluation hook.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
void preEvaluate() const
Pre-evaluation hook.
static word type()
Return the type name.
const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
const dimensionedScalar c
Speed of light in a vacuum.
Janev, Langer, Evans and Post reaction rate.
A class for handling words, derived from string.
Definition: word.H:59
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar dcidT(const scalar p, const scalar T, const scalarField &c, const label li) const
Temperature derivative of the pressure dependent term.
void dcidc(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
A wordList with hashed indices for faster lookup by name.
Ostream & operator<<(Ostream &, const ensightPart &)
JanevReactionRate(const scalar A, const scalar beta, const scalar Ta, const FixedList< scalar, nb_ > b)
Construct from components.
volScalarField & p
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
void write(Ostream &os) const
Write to stream.
Namespace for OpenFOAM.