powerSeriesReactionRateI.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 \*---------------------------------------------------------------------------*/
25 
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 (
32  const scalar A,
33  const scalar beta,
34  const FixedList<scalar, nCoeff_> coeffs
35 )
36 :
37  beta_(beta),
38  A_(A),
39  coeffs_(coeffs)
40 {}
41 
42 
44 (
45  const speciesTable&,
46  const dimensionSet& dims,
47  const dictionary& dict
48 )
49 :
50  beta_(dict.lookup<scalar>("beta", dimless)),
51  A_(dict.lookup<scalar>("A", dims/pow(dimTemperature, beta_))),
52  coeffs_(dict.lookup("coeffs"))
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
57 
59 {}
60 
61 
63 {}
64 
65 
66 inline Foam::scalar Foam::powerSeriesReactionRate::operator()
67 (
68  const scalar p,
69  const scalar T,
70  const scalarField&,
71  const label
72 ) const
73 {
74  scalar lta = A_;
75 
76  if (mag(beta_) > vSmall)
77  {
78  lta *= pow(T, beta_);
79  }
80 
81  scalar expArg = 0;
82 
83  forAll(coeffs_, n)
84  {
85  expArg += coeffs_[n]/pow(T, n + 1);
86  }
87 
88  lta *= exp(expArg);
89 
90  return lta;
91 }
92 
93 
95 (
96  const scalar p,
97  const scalar T,
98  const scalarField&,
99  const label
100 ) const
101 {
102  scalar lta = A_;
103 
104  if (mag(beta_) > vSmall)
105  {
106  lta *= pow(T, beta_);
107  }
108 
109  scalar expArg = 0;
110  scalar deriv = 0;
111 
112  forAll(coeffs_, n)
113  {
114  scalar cT = coeffs_[n]/pow(T, n + 1);
115  expArg += cT;
116  deriv -= (n + 1)*cT;
117  }
118 
119  lta *= exp(expArg);
120 
121  return lta*(beta_+deriv)/T;
122 }
123 
124 
126 {
127  return false;
128 }
129 
130 
132 (
133  const scalar p,
134  const scalar T,
135  const scalarField& c,
136  const label li,
137  scalarField& ddc
138 ) const
139 {
140  ddc = 0;
141 }
142 
143 
145 {
146  writeEntry(os, "A", A_);
147  writeEntry(os, "beta", beta_);
148  writeEntry(os, "coeffs", coeffs_);
149 }
150 
151 
152 inline Foam::Ostream& Foam::operator<<
153 (
154  Ostream& os,
155  const powerSeriesReactionRate& psrr
156 )
157 {
158  psrr.write(os);
159  return os;
160 }
161 
162 
163 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual Ostream & write(const char)=0
Write character.
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.
Power series reaction rate.
void preEvaluate() const
Pre-evaluation hook.
void postEvaluate() const
Post-evaluation hook.
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.
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
The derivative of the rate w.r.t. temperature.
powerSeriesReactionRate(const scalar A, const scalar beta, const FixedList< scalar, nCoeff_ > coeffs)
Construct from components.
bool hasDdc() const
Is the rate a function of concentration?
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar exp(const dimensionedScalar &ds)
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
const dimensionSet dimless
const dimensionSet dimTemperature
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< scalar > mag(const dimensioned< Type > &)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p