FallOffReactionRateI.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 
28 template<class ReactionRate, class FallOffFunction>
31 (
32  const ReactionRate& k0,
33  const ReactionRate& kInf,
34  const FallOffFunction& F,
35  const thirdBodyEfficiencies& tbes
36 )
37 :
38  k0_(k0),
39  kInf_(kInf),
40  F_(F),
41  thirdBodyEfficiencies_(tbes)
42 {}
43 
44 
45 template<class ReactionRate, class FallOffFunction>
48 (
49  const speciesTable& species,
50  const dictionary& dict
51 )
52 :
53  k0_(species, dict.subDict("k0")),
54  kInf_(species, dict.subDict("kInf")),
55  F_(dict.subDict("F")),
56  thirdBodyEfficiencies_(species, dict.subDict("thirdBodyEfficiencies"))
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
62 template<class ReactionRate, class FallOffFunction>
63 inline Foam::scalar
65 (
66  const scalar p,
67  const scalar T,
68  const scalarField& c
69 ) const
70 {
71  scalar k0 = k0_(p, T, c);
72  scalar kInf = kInf_(p, T, c);
73 
74  scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
75 
76  return kInf*(Pr/(1 + Pr))*F_(T, Pr);
77 }
78 
79 
80 template<class ReactionRate, class FallOffFunction>
82 (
83  Ostream& os
84 ) const
85 {
86  os << indent << "k0" << nl;
87  os << indent << token::BEGIN_BLOCK << nl;
88  os << incrIndent;
89  k0_.write(os);
90  os << decrIndent;
91  os << indent << token::END_BLOCK << nl;
92 
93  os << indent << "kInf" << nl;
94  os << indent << token::BEGIN_BLOCK << nl;
95  os << incrIndent;
96  kInf_.write(os);
97  os << decrIndent;
98  os << indent << token::END_BLOCK << nl;
99 
100  os << indent << "F" << nl;
101  os << indent << token::BEGIN_BLOCK << nl;
102  os << incrIndent;
103  F_.write(os);
104  os << decrIndent;
105  os << indent << token::END_BLOCK << nl;
106 
107  os << indent << "thirdBodyEfficiencies" << nl;
108  os << indent << token::BEGIN_BLOCK << nl;
109  os << incrIndent;
110  thirdBodyEfficiencies_.write(os);
111  os << decrIndent;
112  os << indent << token::END_BLOCK << nl;
113 }
114 
115 
116 template<class ReactionRate, class FallOffFunction>
117 inline Foam::Ostream& Foam::operator<<
118 (
119  Ostream& os,
121 )
122 {
123  forr.write(os);
124  return os;
125 }
126 
127 
128 // ************************************************************************* //
dimensionedScalar Pr("Pr", dimless, laminarTransport)
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:223
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
General class for handling unimolecular/recombination fall-off reactions.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
void write(Ostream &os) const
Write to stream.
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 & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:237
const volScalarField & T
FallOffReactionRate(const ReactionRate &k0, const ReactionRate &kInf, const FallOffFunction &F, const thirdBodyEfficiencies &tbes)
Construct from components.
friend Ostream & operator(Ostream &, const FallOffReactionRate< ReactionRate, FallOffFunction > &)
A wordList with hashed indices for faster lookup by name.
const dimensionedScalar c
Speed of light in a vacuum.
virtual Ostream & write(const token &)=0
Write next token to stream.
volScalarField & p
Third body efficiencies.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:230