LandauTellerReactionRate.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-2024 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::LandauTellerReactionRate
26 
27 Description
28  Landau-Teller reaction rate.
29 
30 SourceFiles
31  LandauTellerReactionRateI.H
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef LandauTellerReactionRate_H
36 #define LandauTellerReactionRate_H
37 
38 #include "speciesTable.H"
39 #include "scalarField.H"
40 #include "typeInfo.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 // Forward declaration of friend functions and operators
48 
49 class LandauTellerReactionRate;
50 
51 Ostream& operator<<(Ostream&, const LandauTellerReactionRate&);
52 
53 
54 /*---------------------------------------------------------------------------*\
55  Class LandauTellerReactionRate Declaration
56 \*---------------------------------------------------------------------------*/
57 
59 {
60  // Private Data
61 
62  scalar beta_;
63  scalar A_;
64  scalar Ta_;
65  scalar B_;
66  scalar C_;
67 
68 
69 public:
70 
71  // Constructors
72 
73  //- Construct from components
75  (
76  const scalar A,
77  const scalar beta,
78  const scalar Ta,
79  const scalar B,
80  const scalar C
81  );
82 
83  //- Construct from dictionary
85  (
86  const speciesTable& species,
87  const dimensionSet& dims,
88  const dictionary& dict
89  );
90 
91 
92  // Member Functions
93 
94  //- Return the type name
95  static word type()
96  {
97  return "LandauTeller";
98  }
99 
100  //- Pre-evaluation hook
101  inline void preEvaluate() const;
102 
103  //- Post-evaluation hook
104  inline void postEvaluate() const;
105 
106  //- Return the rate
107  inline scalar operator()
108  (
109  const scalar p,
110  const scalar T,
111  const scalarField& c,
112  const label li
113  ) const;
114 
115  //- The derivative of the rate w.r.t. temperature
116  inline scalar ddT
117  (
118  const scalar p,
119  const scalar T,
120  const scalarField& c,
121  const label li
122  ) const;
123 
124  //- Is the rate a function of concentration?
125  inline bool hasDdc() const;
126 
127  //- The derivative of the rate w.r.t. concentration
128  inline void ddc
129  (
130  const scalar p,
131  const scalar T,
132  const scalarField& c,
133  const label li,
135  ) const;
136 
137  //- Write to stream
138  inline void write(Ostream& os) const;
139 
140 
141  // Ostream Operator
142 
143  inline friend Ostream& operator<<
144  (
145  Ostream&,
147  );
148 };
149 
150 
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 
153 } // End namespace Foam
154 
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
161 #endif
162 
163 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("B", Foam::dimless, 18.678)
Graphite solid properties.
Definition: C.H:51
Landau-Teller reaction rate.
void preEvaluate() const
Pre-evaluation hook.
void postEvaluate() const
Post-evaluation hook.
static word type()
Return the type name.
void write(Ostream &os) const
Write to stream.
void ddc(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &ddc) const
The derivative of the rate w.r.t. concentration.
LandauTellerReactionRate(const scalar A, const scalar beta, const scalar Ta, const scalar B, const scalar C)
Construct from components.
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
The derivative of the rate w.r.t. temperature.
bool hasDdc() const
Is the rate a function of concentration?
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
A wordList with hashed indices for faster lookup by name.
A class for handling words, derived from string.
Definition: word.H:62
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
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
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p
Basic run-time type information using word as the type's name. Used to enhance the standard RTTI to c...