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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 (
32  const scalar A,
33  const scalar beta,
34  const scalar Ta,
36 )
37 :
38  A_(A),
39  beta_(beta),
40  Ta_(Ta),
41  b_(b)
42 {}
43 
44 
46 (
47  const speciesTable&,
48  const dictionary& dict
49 )
50 :
51  A_(dict.get<scalar>("A")),
52  beta_(dict.get<scalar>("beta")),
53  Ta_(dict.get<scalar>("Ta")),
54  b_(dict.lookup("b"))
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
59 
60 inline Foam::scalar Foam::JanevReactionRate::operator()
61 (
62  const scalar p,
63  const scalar T,
64  const scalarField&
65 ) const
66 {
67  scalar lta = A_;
68 
69  if (mag(beta_) > VSMALL)
70  {
71  lta *= pow(T, beta_);
72  }
73 
74  scalar expArg = 0.0;
75 
76  if (mag(Ta_) > VSMALL)
77  {
78  expArg -= Ta_/T;
79  }
80 
81  scalar lnT = log(T);
82 
83  for (int n=0; n<nb_; n++)
84  {
85  expArg += b_[n]*pow(lnT, n);
86  }
87 
88  lta *= exp(expArg);
89 
90  return lta;
91 }
92 
93 
95 {
96  os.writeKeyword("A") << A_ << nl;
97  os.writeKeyword("beta") << beta_ << nl;
98  os.writeKeyword("Ta") << Ta_ << nl;
99  os.writeKeyword("b") << b_ << nl;
100 }
101 
102 
103 inline Foam::Ostream& Foam::operator<<
104 (
105  Ostream& os,
106  const JanevReactionRate& jrr
107 )
108 {
109  jrr.write(os);
110  return os;
111 }
112 
113 
114 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
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:27
Foam::Field< scalar >
Foam::hashedWordList
A wordList with hashed named lookup, which can be faster in some situations than using the normal lis...
Definition: hashedWordList.H:54
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::JanevReactionRate
Janev, Langer, Evans and Post reaction rate.
Definition: JanevReactionRate.H:59
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
os
OBJstream os(runTime.globalPath()/outputName)
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::JanevReactionRate::JanevReactionRate
JanevReactionRate(const scalar A, const scalar beta, const scalar Ta, const FixedList< scalar, nb_ > b)
Construct from components.
Definition: JanevReactionRateI.H:31
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::FixedList< scalar, nb_ >
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::JanevReactionRate::write
void write(Ostream &os) const
Write to stream.
Definition: JanevReactionRateI.H:94