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-2018 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_(readScalar(dict.lookup("A"))),
48  beta_(readScalar(dict.lookup("beta"))),
49  Ta_(readScalar(dict.lookup("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
61 {
62  scalar ak = A_;
63 
64  if (mag(beta_) > vSmall)
65  {
66  ak *= pow(T, beta_);
67  }
68 
69  if (mag(Ta_) > vSmall)
70  {
71  ak *= exp(-Ta_/T);
72  }
73 
74  return ak;
75 }
76 
77 
78 inline Foam::scalar Foam::ArrheniusReactionRate::ddT
79 (
80  const scalar p,
81  const scalar T,
82  const scalarField&
83 ) const
84 {
85  scalar ak = A_;
86 
87  if (mag(beta_) > vSmall)
88  {
89  ak *= pow(T, beta_);
90  }
91 
92  if (mag(Ta_) > vSmall)
93  {
94  ak *= exp(-Ta_/T);
95  }
96 
97  return ak*(beta_+Ta_/T)/T;
98 }
99 
100 
103 {
104  return NullObjectRef<List<Tuple2<label, scalar>>>();
105 }
106 
107 
109 (
110  const scalar p,
111  const scalar T,
112  const scalarField& c,
113  scalarField& dcidc
114 ) const
115 {}
116 
117 
118 inline Foam::scalar Foam::ArrheniusReactionRate::dcidT
119 (
120  const scalar p,
121  const scalar T,
122  const scalarField& c
123 ) const
124 {
125  return 0;
126 }
127 
128 
130 {
131  os.writeKeyword("A") << A_ << token::END_STATEMENT << nl;
132  os.writeKeyword("beta") << beta_ << token::END_STATEMENT << nl;
133  os.writeKeyword("Ta") << Ta_ << token::END_STATEMENT << nl;
134 }
135 
136 
137 inline Foam::Ostream& Foam::operator<<
138 (
139  Ostream& os,
140  const ArrheniusReactionRate& arr
141 )
142 {
143  arr.write(os);
144  return os;
145 }
146 
147 
148 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
ArrheniusReactionRate(const scalar A, const scalar beta, const scalar Ta)
Construct from components.
dimensionedScalar exp(const dimensionedScalar &ds)
scalar ddT(const scalar p, const scalar T, const scalarField &c) const
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:265
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
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, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual Ostream & write(const token &)=0
Write next token to stream.
volScalarField & p
void write(Ostream &os) const
Write to stream.
scalar dcidT(const scalar p, const scalar T, const scalarField &c) const
Temperature derivative of the pressure dependent term.
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:576