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-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,
34 )
35 :
36  A_(A),
37  beta_(beta),
38  Ta_(Ta),
39  b_(b)
40 {}
41 
42 
44 (
45  const speciesTable&,
46  const dictionary& dict
47 )
48 :
49  A_(readScalar(dict.lookup("A"))),
50  beta_(readScalar(dict.lookup("beta"))),
51  Ta_(readScalar(dict.lookup("Ta"))),
52  b_(dict.lookup("b"))
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
57 
58 inline Foam::scalar Foam::JanevReactionRate::operator()
59 (
60  const scalar p,
61  const scalar T,
62  const scalarField&
63 ) const
64 {
65  scalar lta = A_;
66 
67  if (mag(beta_) > VSMALL)
68  {
69  lta *= pow(T, beta_);
70  }
71 
72  scalar expArg = 0.0;
73 
74  if (mag(Ta_) > VSMALL)
75  {
76  expArg -= Ta_/T;
77  }
78 
79  scalar lnT = log(T);
80 
81  for (int n=0; n<nb_; n++)
82  {
83  expArg += b_[n]*pow(lnT, n);
84  }
85 
86  lta *= exp(expArg);
87 
88  return lta;
89 }
90 
91 
93 {
94  os.writeKeyword("A") << A_ << nl;
95  os.writeKeyword("beta") << beta_ << nl;
96  os.writeKeyword("Ta") << Ta_ << nl;
97  os.writeKeyword("b") << b_ << nl;
98 }
99 
100 
101 inline Foam::Ostream& Foam::operator<<
102 (
103  Ostream& os,
104  const JanevReactionRate& jrr
105 )
106 {
107  jrr.write(os);
108  return os;
109 }
110 
111 
112 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Janev, Langer, Evans and Post reaction rate.
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.
JanevReactionRate(const scalar A, const scalar beta, const scalar Ta, const FixedList< scalar, nb_ > b)
Construct from components.
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
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