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 | 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 
30 template<class ReactionRate, class FallOffFunction>
33 (
34  const ReactionRate& k0,
35  const ReactionRate& kInf,
36  const FallOffFunction& F,
37  const thirdBodyEfficiencies& tbes
38 )
39 :
40  k0_(k0),
41  kInf_(kInf),
42  F_(F),
43  thirdBodyEfficiencies_(tbes)
44 {}
45 
46 
47 template<class ReactionRate, class FallOffFunction>
50 (
51  const speciesTable& species,
52  const dictionary& dict
53 )
54 :
55  k0_(species, dict.subDict("k0")),
56  kInf_(species, dict.subDict("kInf")),
57  F_(dict.subDict("F")),
58  thirdBodyEfficiencies_(species, dict.subDict("thirdBodyEfficiencies"))
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
64 template<class ReactionRate, class FallOffFunction>
65 inline Foam::scalar
67 (
68  const scalar p,
69  const scalar T,
70  const scalarField& c
71 ) const
72 {
73  const scalar k0 = k0_(p, T, c);
74  const scalar kInf = kInf_(p, T, c);
75 
76  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
77 
78  return kInf*(Pr/(1 + Pr))*F_(T, Pr);
79 }
80 
81 
82 template<class ReactionRate, class FallOffFunction>
84 (
85  Ostream& os
86 ) const
87 {
88  os.beginBlock("k0");
89  k0_.write(os);
90  os.endBlock();
91 
92  os.beginBlock("kInf");
93  kInf_.write(os);
94  os.endBlock();
95 
96  os.beginBlock("F");
97  F_.write(os);
98  os.endBlock();
99 
100  os.beginBlock("thirdBodyEfficiencies");
101  thirdBodyEfficiencies_.write(os);
102  os.endBlock();
103 }
104 
105 
106 template<class ReactionRate, class FallOffFunction>
107 inline Foam::Ostream& Foam::operator<<
108 (
109  Ostream& os,
111 )
112 {
113  forr.write(os);
114  return os;
115 }
116 
117 
118 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::FallOffReactionRate::operator
friend Ostream & operator(Ostream &, const FallOffReactionRate< ReactionRate, FallOffFunction > &)
F
volVectorField F(fluid.F())
Foam::thirdBodyEfficiencies
Third body efficiencies.
Definition: thirdBodyEfficiencies.H:57
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
Foam::FallOffReactionRate
General class for handling unimolecular/recombination fall-off reactions.
Definition: FallOffReactionRate.H:49
Pr
dimensionedScalar Pr("Pr", dimless, laminarTransport)
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
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::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::FallOffReactionRate::FallOffReactionRate
FallOffReactionRate(const ReactionRate &k0, const ReactionRate &kInf, const FallOffFunction &F, const thirdBodyEfficiencies &tbes)
Construct from components.
Definition: FallOffReactionRateI.H:33
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::FallOffReactionRate::write
void write(Ostream &os) const
Write to stream.
Definition: FallOffReactionRateI.H:84