ArrheniusReactionRateI.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 )
34 :
35  A_(A),
36  beta_(beta),
37  Ta_(Ta)
38 {}
39 
40 
42 (
43  const speciesTable&,
44  const dictionary& dict
45 )
46 :
47  A_(dict.lookup<scalar>("A")),
48  beta_(dict.lookup<scalar>("beta")),
49  Ta_(dict.lookup<scalar>("Ta"))
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
54 
55 inline Foam::scalar Foam::ArrheniusReactionRate::operator()
56 (
57  const scalar p,
58  const scalar T,
59  const scalarField&,
60  const label
61 ) const
62 {
63  scalar ak = A_;
64 
65  if (mag(beta_) > vSmall)
66  {
67  ak *= pow(T, beta_);
68  }
69 
70  if (mag(Ta_) > vSmall)
71  {
72  ak *= exp(-Ta_/T);
73  }
74 
75  return ak;
76 }
77 
78 
79 inline Foam::scalar Foam::ArrheniusReactionRate::ddT
80 (
81  const scalar p,
82  const scalar T,
83  const scalarField&,
84  const label
85 ) const
86 {
87  scalar ak = A_;
88 
89  if (mag(beta_) > vSmall)
90  {
91  ak *= pow(T, beta_);
92  }
93 
94  if (mag(Ta_) > vSmall)
95  {
96  ak *= exp(-Ta_/T);
97  }
98 
99  return ak*(beta_+Ta_/T)/T;
100 }
101 
102 
105 {
106  return NullObjectRef<List<Tuple2<label, scalar>>>();
107 }
108 
109 
111 (
112  const scalar p,
113  const scalar T,
114  const scalarField& c,
115  const label li,
116  scalarField& dcidc
117 ) const
118 {}
119 
120 
121 inline Foam::scalar Foam::ArrheniusReactionRate::dcidT
122 (
123  const scalar p,
124  const scalar T,
125  const scalarField& c,
126  const label li
127 ) const
128 {
129  return 0;
130 }
131 
132 
134 {
135  writeEntry(os, "A", A_);
136  writeEntry(os, "beta", beta_);
137  writeEntry(os, "Ta", Ta_);
138 }
139 
140 
141 inline Foam::Ostream& Foam::operator<<
142 (
143  Ostream& os,
144  const ArrheniusReactionRate& arr
145 )
146 {
147  arr.write(os);
148  return os;
149 }
150 
151 
152 // ************************************************************************* //
dictionary dict
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
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
ArrheniusReactionRate(const scalar A, const scalar beta, const scalar Ta)
Construct from components.
dimensionedScalar exp(const dimensionedScalar &ds)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
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.
Arrhenius reaction rate given by:
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.
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar dcidT(const scalar p, const scalar T, const scalarField &c, const label li) const
Temperature derivative of the pressure dependent term.
volScalarField & p
void write(Ostream &os) const
Write to stream.
const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812