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-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 \*---------------------------------------------------------------------------*/
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_(readScalar(dict.lookup("A"))),
50  beta_(readScalar(dict.lookup("beta"))),
51  Ta_(readScalar(dict.lookup("Ta"))),
52  coeffs_(dict.lookup("coeffs"))
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
57 
58 inline Foam::scalar Foam::powerSeriesReactionRate::operator()
59 (
60  const scalar p,
61  const scalar T,
62  const scalarField&
63 ) const
64 {
65  scalar lta = A_;
66 
67  if (mag(beta_) > vSmall)
68  {
69  lta *= pow(T, beta_);
70  }
71 
72  scalar expArg = 0;
73 
74  forAll(coeffs_, n)
75  {
76  expArg += coeffs_[n]/pow(T, n + 1);
77  }
78 
79  lta *= exp(expArg);
80 
81  return lta;
82 }
83 
84 
85 inline Foam::scalar Foam::powerSeriesReactionRate::ddT
86 (
87  const scalar p,
88  const scalar T,
89  const scalarField&
90 ) const
91 {
92  scalar lta = A_;
93 
94  if (mag(beta_) > vSmall)
95  {
96  lta *= pow(T, beta_);
97  }
98 
99  scalar expArg = 0;
100  scalar deriv = 0;
101 
102  forAll(coeffs_, n)
103  {
104  scalar cT = coeffs_[n]/pow(T, n + 1);
105  expArg += cT;
106  deriv -= (n + 1)*cT;
107  }
108 
109  lta *= exp(expArg);
110 
111  return lta*(beta_+deriv)/T;
112 }
113 
114 
117 {
118  return NullObjectRef<List<Tuple2<label, scalar>>>();
119 }
120 
121 
123 (
124  const scalar p,
125  const scalar T,
126  const scalarField& c,
127  scalarField& dcidc
128 ) const
129 {}
130 
131 
132 inline Foam::scalar Foam::powerSeriesReactionRate::dcidT
133 (
134  const scalar p,
135  const scalar T,
136  const scalarField& c
137 ) const
138 {
139  return 0;
140 }
141 
142 
144 {
145  writeEntry(os, "A", A_);
146  writeEntry(os, "beta", beta_);
147  writeEntry(os, "Ta", Ta_);
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 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
#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)
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:158
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
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:53
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const volScalarField & T
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
A wordList with hashed indices for faster lookup by name.
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
scalar dcidT(const scalar p, const scalar T, const scalarField &c) const
Temperature derivative of the pressure dependent term.
virtual Ostream & write(const token &)=0
Write next token to stream.
volScalarField & p
scalar ddT(const scalar p, const scalar T, const scalarField &c) const
void dcidc(const scalar p, const scalar T, const scalarField &c, 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:583