solidArrheniusReactionRate.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::solidArrheniusReactionRate
26 
27 Description
28  Arrhenius reaction rate for solids
29 
30 SourceFiles
31  solidArrheniusReactionRateI.H
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef solidArrheniusReactionRate_H
36 #define solidArrheniusReactionRate_H
37 
38 #include "scalarField.H"
39 #include "typeInfo.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 // Forward declaration of friend functions and operators
47 
48 class solidArrheniusReactionRate;
49 
50 Ostream& operator<<(Ostream&, const solidArrheniusReactionRate&);
51 
52 
53 /*---------------------------------------------------------------------------*\
54  Class solidArrheniusReactionRate Declaration
55 \*---------------------------------------------------------------------------*/
56 
58 {
59  // Private Data
60 
61  scalar A_;
62  scalar Ta_;
63  scalar Tcrit_;
64 
65 
66 public:
67 
68  // Constructors
69 
70  //- Construct from components
72  (
73  const scalar A,
74  const scalar Ta,
75  const scalar Tcrit
76  );
77 
78 
79  //- Construct from dictionary
81  (
82  const speciesTable& species,
83  const dictionary& dict
84  );
85 
86 
87  //- Destructor
89  {}
90 
91 
92  // Member Functions
93 
94  //- Return the type name
95  static word type()
96  {
97  return "Arrhenius";
98  }
99 
100  inline scalar operator()
101  (
102  const scalar p,
103  const scalar T,
104  const scalarField& c
105  ) const;
106 
107  inline scalar ddT
108  (
109  const scalar p,
110  const scalar T,
111  const scalarField& c
112  ) const;
113 
114  //- Third-body efficiencies (beta = 1-alpha)
115  // non-empty only for third-body reactions
116  // with enhanced molecularity (alpha != 1)
117  inline const List<Tuple2<label, scalar>>& beta() const;
118 
119  //- Species concentration derivative of the pressure dependent term
120  // By default this value is 1 as it multiplies the third-body term
121  inline void dcidc
122  (
123  const scalar p,
124  const scalar T,
125  const scalarField& c,
127  ) const;
128 
129  //- Temperature derivative of the pressure dependent term
130  // By default this value is 0 since ddT of molecularity is approx.0
131  inline scalar dcidT
132  (
133  const scalar p,
134  const scalar T,
135  const scalarField& c
136  ) const;
137 
138 
139  //- Write to stream
140  inline void write(Ostream& os) const;
141 
142 
143  // Ostream Operator
144 
145  inline friend Ostream& operator<<
146  (
147  Ostream&,
149  );
150 };
151 
152 
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 
155 } // End namespace Foam
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
160 
161 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162 
163 #endif
164 
165 // ************************************************************************* //
dictionary dict
const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
scalar dcidT(const scalar p, const scalar T, const scalarField &c) const
Temperature derivative of the pressure dependent term.
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
static word type()
Return the type name.
void write(Ostream &os) const
Write to stream.
A class for handling words, derived from string.
Definition: word.H:59
Arrhenius reaction rate for solids.
solidArrheniusReactionRate(const scalar A, const scalar Ta, const scalar Tcrit)
Construct from components.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
scalar ddT(const scalar p, const scalar T, const scalarField &c) const
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.
void dcidc(const scalar p, const scalar T, const scalarField &c, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
Ostream & operator<<(Ostream &, const ensightPart &)
volScalarField & p
Namespace for OpenFOAM.