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-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 
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  Istream& is
51 )
52 :
53  k0_(species, is.readBegin("FallOffReactionRate(Istream&)")),
54  kInf_(species, is),
55  F_(is),
56  thirdBodyEfficiencies_(species, is)
57 {
58  is.readEnd("FallOffReactionRate(Istream&)");
59 }
60 
61 
62 template<class ReactionRate, class FallOffFunction>
65 (
66  const speciesTable& species,
67  const dictionary& dict
68 )
69 :
70  k0_(species, dict.subDict("k0")),
71  kInf_(species, dict.subDict("kInf")),
72  F_(dict.subDict("F")),
73  thirdBodyEfficiencies_(species, dict.subDict("thirdBodyEfficiencies"))
74 {}
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
79 template<class ReactionRate, class FallOffFunction>
80 inline Foam::scalar
82 (
83  const scalar p,
84  const scalar T,
85  const scalarField& c
86 ) const
87 {
88  scalar k0 = k0_(p, T, c);
89  scalar kInf = kInf_(p, T, c);
90 
91  scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
92 
93  return kInf*(Pr/(1 + Pr))*F_(T, Pr);
94 }
95 
96 
97 template<class ReactionRate, class FallOffFunction>
99 (
100  Ostream& os
101 ) const
102 {
103  os << indent << "k0" << nl;
104  os << indent << token::BEGIN_BLOCK << nl;
105  os << incrIndent;
106  k0_.write(os);
107  os << decrIndent;
108  os << indent << token::END_BLOCK << nl;
109 
110  os << indent << "kInf" << nl;
111  os << indent << token::BEGIN_BLOCK << nl;
112  os << incrIndent;
113  kInf_.write(os);
114  os << decrIndent;
115  os << indent << token::END_BLOCK << nl;
116 
117  os << indent << "F" << nl;
118  os << indent << token::BEGIN_BLOCK << nl;
119  os << incrIndent;
120  F_.write(os);
121  os << decrIndent;
122  os << indent << token::END_BLOCK << nl;
123 
124  os << indent << "thirdBodyEfficiencies" << nl;
125  os << indent << token::BEGIN_BLOCK << nl;
126  os << incrIndent;
127  thirdBodyEfficiencies_.write(os);
128  os << decrIndent;
129  os << indent << token::END_BLOCK << nl;
130 }
131 
132 
133 template<class ReactionRate, class FallOffFunction>
134 inline Foam::Ostream& Foam::operator<<
135 (
136  Ostream& os,
138 )
139 {
140  os << token::BEGIN_LIST
141  << forr.k0_ << token::SPACE
142  << forr.kInf_ << token::SPACE
143  << forr.F_ << token::SPACE
144  << forr.thirdBodyEfficiencies_
145  << token::END_LIST;
146  return os;
147 }
148 
149 
150 // ************************************************************************* //
void write(Ostream &os) const
Write to stream.
Istream & readBegin(const char *funcName)
Definition: Istream.C:86
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
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
General class for handling unimolecular/recombination fall-off reactions.
Istream & readEnd(const char *funcName)
Definition: Istream.C:103
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