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-2017 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 
79 {
80  os.writeKeyword("A") << A_ << token::END_STATEMENT << nl;
81  os.writeKeyword("beta") << beta_ << token::END_STATEMENT << nl;
82  os.writeKeyword("Ta") << Ta_ << token::END_STATEMENT << nl;
83 }
84 
85 
86 inline Foam::Ostream& Foam::operator<<
87 (
88  Ostream& os,
89  const ArrheniusReactionRate& arr
90 )
91 {
92  arr.write(os);
93  return os;
94 }
95 
96 
97 // ************************************************************************* //
#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
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:53
static const char nl
Definition: Ostream.H:262
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:
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.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576