JanevReactionRateI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2012 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,
34 )
35 :
36  A_(A),
37  beta_(beta),
38  Ta_(Ta),
39  b_(b)
40 {}
41 
42 
44 (
45  const speciesTable&,
46  Istream& is
47 )
48 :
49  A_(readScalar(is.readBegin("JanevReactionRate(Istream&)"))),
50  beta_(readScalar(is)),
51  Ta_(readScalar(is)),
52  b_(is)
53 {
54  is.readEnd("JanevReactionRate(Istream&)");
55 }
56 
57 
59 (
60  const speciesTable&,
61  const dictionary& dict
62 )
63 :
64  A_(readScalar(dict.lookup("A"))),
65  beta_(readScalar(dict.lookup("beta"))),
66  Ta_(readScalar(dict.lookup("Ta"))),
67  b_(dict.lookup("b"))
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
73 inline Foam::scalar Foam::JanevReactionRate::operator()
74 (
75  const scalar p,
76  const scalar T,
77  const scalarField&
78 ) const
79 {
80  scalar lta = A_;
81 
82  if (mag(beta_) > VSMALL)
83  {
84  lta *= pow(T, beta_);
85  }
86 
87  scalar expArg = 0.0;
88 
89  if (mag(Ta_) > VSMALL)
90  {
91  expArg -= Ta_/T;
92  }
93 
94  scalar lnT = log(T);
95 
96  for (int n=0; n<nb_; n++)
97  {
98  expArg += b_[n]*pow(lnT, n);
99  }
100 
101  lta *= exp(expArg);
102 
103  return lta;
104 }
105 
106 
108 {
109  os.writeKeyword("A") << A_ << nl;
110  os.writeKeyword("beta") << beta_ << nl;
111  os.writeKeyword("Ta") << Ta_ << nl;
112  os.writeKeyword("b") << b_ << nl;
113 }
114 
115 
116 inline Foam::Ostream& Foam::operator<<
117 (
118  Ostream& os,
119  const JanevReactionRate& jrr
120 )
121 {
122  os << token::BEGIN_LIST
123  << jrr.A_ << token::SPACE << jrr.beta_ << token::SPACE << jrr.Ta_;
124 
125  for (int n=0; n<JanevReactionRate::nb_; n++)
126  {
127  os << token::SPACE << jrr.b_[n];
128  }
129 
130  os << token::END_LIST;
131 
132  return os;
133 }
134 
135 
136 // ************************************************************************* //
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::Istream::readEnd
Istream & readEnd(const char *funcName)
Definition: Istream.C:105
p
p
Definition: pEqn.H:62
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
A
simpleMatrix< scalar > A(Nc)
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::hashedWordList
A wordList with hashed indices for faster lookup by name.
Definition: hashedWordList.H:57
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::JanevReactionRate::b_
FixedList< scalar, nb_ > b_
Definition: JanevReactionRate.H:59
Foam::JanevReactionRate::nb_
static const label nb_
Definition: JanevReactionRate.H:58
Foam::JanevReactionRate::A_
scalar A_
Definition: JanevReactionRate.H:54
Foam::JanevReactionRate
Janev, Langer, Evans and Post reaction rate.
Definition: JanevReactionRate.H:50
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::JanevReactionRate::Ta_
scalar Ta_
Definition: JanevReactionRate.H:56
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:253
readScalar
#define readScalar
Definition: doubleScalar.C:38
T
const volScalarField & T
Definition: createFields.H:25
Foam::token::BEGIN_LIST
@ BEGIN_LIST
Definition: token.H:100
Foam::JanevReactionRate::JanevReactionRate
JanevReactionRate(const scalar A, const scalar beta, const scalar Ta, const FixedList< scalar, nb_ > b)
Construct from components.
Definition: JanevReactionRateI.H:29
Foam::FixedList< scalar, nb_ >
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::JanevReactionRate::beta_
scalar beta_
Definition: JanevReactionRate.H:55
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::Istream::readBegin
Istream & readBegin(const char *funcName)
Definition: Istream.C:88
Foam::JanevReactionRate::write
void write(Ostream &os) const
Write to stream.
Definition: JanevReactionRateI.H:107
Foam::token::END_LIST
@ END_LIST
Definition: token.H:101
Foam::token::SPACE
@ SPACE
Definition: token.H:95