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-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 \*---------------------------------------------------------------------------*/
25 
26 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
29 (
30  const scalar A,
31  const scalar beta,
32  const scalar Ta,
33  const FixedList<scalar, nCoeff_> coeffs
34 )
35 :
36  A_(A),
37  beta_(beta),
38  Ta_(Ta),
39  coeffs_(coeffs)
40 {}
41 
42 
44 (
45  const speciesTable&,
46  const dictionary& dict
47 )
48 :
49  A_(dict.lookup<scalar>("A")),
50  beta_(dict.lookup<scalar>("beta")),
51  Ta_(dict.lookup<scalar>("Ta")),
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 
94 inline Foam::scalar Foam::powerSeriesReactionRate::ddT
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 
127 {
128  return NullObjectRef<List<Tuple2<label, scalar>>>();
129 }
130 
131 
133 (
134  const scalar p,
135  const scalar T,
136  const scalarField& c,
137  const label li,
138  scalarField& dcidc
139 ) const
140 {}
141 
142 
143 inline Foam::scalar Foam::powerSeriesReactionRate::dcidT
144 (
145  const scalar p,
146  const scalar T,
147  const scalarField& c,
148  const label li
149 ) const
150 {
151  return 0;
152 }
153 
154 
156 {
157  writeEntry(os, "A", A_);
158  writeEntry(os, "beta", beta_);
159  writeEntry(os, "Ta", Ta_);
160  writeEntry(os, "coeffs", coeffs_);
161 }
162 
163 
164 inline Foam::Ostream& Foam::operator<<
165 (
166  Ostream& os,
167  const powerSeriesReactionRate& psrr
168 )
169 {
170  psrr.write(os);
171  return os;
172 }
173 
174 
175 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
virtual Ostream & write(const char)=0
Write character.
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 write(Ostream &os) const
Write to stream.
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
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
dimensionedScalar exp(const dimensionedScalar &ds)
powerSeriesReactionRate(const scalar A, const scalar beta, const scalar Ta, const FixedList< scalar, nCoeff_ > coeffs)
Construct from components.
Power series reaction rate.
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)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
void preEvaluate() const
Pre-evaluation hook.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
A wordList with hashed indices for faster lookup by name.
void postEvaluate() const
Post-evaluation hook.
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
volScalarField & p
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.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844