ArrheniusReactionRate.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 | Copyright (C) 2011-2016 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::ArrheniusReactionRate
26 
27 Description
28  Arrhenius reaction rate given by:
29 
30  k = A * T^beta * exp(-Ta/T)
31 
32 SourceFiles
33  ArrheniusReactionRateI.H
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef ArrheniusReactionRate_H
38 #define ArrheniusReactionRate_H
39 
40 #include "scalarField.H"
41 #include "typeInfo.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 // Forward declaration of friend functions and operators
49 
50 class ArrheniusReactionRate;
51 
52 Ostream& operator<<(Ostream&, const ArrheniusReactionRate&);
53 
54 
55 /*---------------------------------------------------------------------------*\
56  Class ArrheniusReactionRate Declaration
57 \*---------------------------------------------------------------------------*/
58 
60 {
61  // Private data
62 
63  scalar A_;
64  scalar beta_;
65  scalar Ta_;
66 
67 
68 public:
69 
70  // Constructors
71 
72  //- Construct from components
74  (
75  const scalar A,
76  const scalar beta,
77  const scalar Ta
78  );
79 
80  //- Construct from Istream
82  (
83  const speciesTable& species,
84  Istream& is
85  );
86 
87  //- Construct from dictionary
89  (
90  const speciesTable& species,
91  const dictionary& dict
92  );
93 
94 
95  // Member Functions
96 
97  //- Return the type name
98  static word type()
99  {
100  return "Arrhenius";
101  }
102 
103  inline scalar operator()
104  (
105  const scalar p,
106  const scalar T,
107  const scalarField& c
108  ) const;
109 
110  //- Write to stream
111  inline void write(Ostream& os) const;
112 
113 
114  // Ostream Operator
115 
116  inline friend Ostream& operator<<
117  (
118  Ostream&,
119  const ArrheniusReactionRate&
120  );
121 };
122 
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 } // End namespace Foam
127 
128 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
129 
130 #include "ArrheniusReactionRateI.H"
131 
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 
134 #endif
135 
136 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
ArrheniusReactionRate(const scalar A, const scalar beta, const scalar Ta)
Construct from components.
A class for handling words, derived from string.
Definition: word.H:59
void write(Ostream &os) const
Write to stream.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static word type()
Return the type name.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
A wordList with hashed indices for faster lookup by name.
const dimensionedScalar c
Speed of light in a vacuum.
Ostream & operator<<(Ostream &, const ensightPart &)
Arrhenius reaction rate given by:
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
volScalarField & p
Namespace for OpenFOAM.