ArrheniusReactionRateI.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
33 )
34 :
35  A_(A),
36  beta_(beta),
37  Ta_(Ta)
38 {}
39 
40 
42 (
43  const speciesTable&,
44  Istream& is
45 )
46 :
47  A_(readScalar(is.readBegin("ArrheniusReactionRate(Istream&)"))),
48  beta_(readScalar(is)),
49  Ta_(readScalar(is))
50 {
51  is.readEnd("ArrheniusReactionRate(Istream&)");
52 }
53 
54 
56 (
57  const speciesTable&,
58  const dictionary& dict
59 )
60 :
61  A_(readScalar(dict.lookup("A"))),
62  beta_(readScalar(dict.lookup("beta"))),
63  Ta_(readScalar(dict.lookup("Ta")))
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
69 inline Foam::scalar Foam::ArrheniusReactionRate::operator()
70 (
71  const scalar p,
72  const scalar T,
73  const scalarField&
74 ) const
75 {
76  scalar ak = A_;
77 
78  if (mag(beta_) > VSMALL)
79  {
80  ak *= pow(T, beta_);
81  }
82 
83  if (mag(Ta_) > VSMALL)
84  {
85  ak *= exp(-Ta_/T);
86  }
87 
88  return ak;
89 }
90 
91 
93 {
94  os.writeKeyword("A") << A_ << token::END_STATEMENT << nl;
95  os.writeKeyword("beta") << beta_ << token::END_STATEMENT << nl;
96  os.writeKeyword("Ta") << Ta_ << token::END_STATEMENT << nl;
97 }
98 
99 
100 inline Foam::Ostream& Foam::operator<<
101 (
102  Ostream& os,
103  const ArrheniusReactionRate& arr
104 )
105 {
106  os << token::BEGIN_LIST
107  << arr.A_ << token::SPACE << arr.beta_ << token::SPACE << arr.Ta_
108  << token::END_LIST;
109  return os;
110 }
111 
112 
113 // ************************************************************************* //
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::Istream::readEnd
Istream & readEnd(const char *funcName)
Definition: Istream.C:105
p
p
Definition: pEqn.H:62
Foam::ArrheniusReactionRate::beta_
scalar beta_
Definition: ArrheniusReactionRate.H:56
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
A
simpleMatrix< scalar > A(Nc)
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::ArrheniusReactionRate::A_
scalar A_
Definition: ArrheniusReactionRate.H:55
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::ArrheniusReactionRate::Ta_
scalar Ta_
Definition: ArrheniusReactionRate.H:57
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::ArrheniusReactionRate::ArrheniusReactionRate
ArrheniusReactionRate(const scalar A, const scalar beta, const scalar Ta)
Construct from components.
Definition: ArrheniusReactionRateI.H:29
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::ArrheniusReactionRate::write
void write(Ostream &os) const
Write to stream.
Definition: ArrheniusReactionRateI.H:92
readScalar
#define readScalar
Definition: doubleScalar.C:38
Foam::ArrheniusReactionRate
Arrhenius reaction rate given by:
Definition: ArrheniusReactionRate.H:51
T
const volScalarField & T
Definition: createFields.H:25
Foam::token::BEGIN_LIST
@ BEGIN_LIST
Definition: token.H:100
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
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::token::END_LIST
@ END_LIST
Definition: token.H:101
Foam::token::SPACE
@ SPACE
Definition: token.H:95