NonEquilibriumReversibleReaction.C
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 
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template
33 <
34  template<class> class ReactionType,
35  class ReactionThermo,
36  class ReactionRate
37 >
39 <
40  ReactionType,
41  ReactionThermo,
42  ReactionRate
43 >::
44 NonEquilibriumReversibleReaction
45 (
46  const ReactionType<ReactionThermo>& reaction,
47  const ReactionRate& forwardReactionRate,
48  const ReactionRate& reverseReactionRate
49 )
50 :
51  ReactionType<ReactionThermo>(reaction),
52  fk_(forwardReactionRate),
53  rk_(reverseReactionRate)
54 {}
55 
56 
57 template
58 <
59  template<class> class ReactionType,
60  class ReactionThermo,
61  class ReactionRate
62 >
64 <
65  ReactionType,
66  ReactionThermo,
67  ReactionRate
68 >::
70 (
71  const speciesTable& species,
72  const HashPtrTable<ReactionThermo>& thermoDatabase,
73  const dictionary& dict
74 )
75 :
76  ReactionType<ReactionThermo>(species, thermoDatabase, dict),
77  fk_(species, dict.subDict("forward")),
78  rk_(species, dict.subDict("reverse"))
79 {}
80 
81 
82 template
83 <
84  template<class> class ReactionType,
85  class ReactionThermo,
86  class ReactionRate
87 >
89 <
90  ReactionType,
91  ReactionThermo,
92  ReactionRate
93 >::
95 (
97  <
98  ReactionType,
99  ReactionThermo,
100  ReactionRate
101  >& nerr,
102  const speciesTable& species
103 )
104 :
105  ReactionType<ReactionThermo>(nerr, species),
106  fk_(nerr.fk_),
107  rk_(nerr.rk_)
108 {}
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
113 template
114 <
115  template<class> class ReactionType,
116  class ReactionThermo,
117  class ReactionRate
118 >
119 Foam::scalar
121 <
122  ReactionType,
123  ReactionThermo,
124  ReactionRate
125 >::kf
126 (
127  const scalar p,
128  const scalar T,
129  const scalarField& c
130 ) const
131 {
132  return fk_(p, T, c);
133 }
134 
135 
136 template
137 <
138  template<class> class ReactionType,
139  class ReactionThermo,
140  class ReactionRate
141 >
142 Foam::scalar
144 <
145  ReactionType,
146  ReactionThermo,
147  ReactionRate
148 >::kr
149 (
150  const scalar,
151  const scalar p,
152  const scalar T,
153  const scalarField& c
154 ) const
155 {
156  return rk_(p, T, c);
157 }
158 
159 
160 template
161 <
162  template<class> class ReactionType,
163  class ReactionThermo,
164  class ReactionRate
165 >
166 Foam::scalar
168 <
169  ReactionType,
170  ReactionThermo,
171  ReactionRate
172 >::kr
173 (
174  const scalar p,
175  const scalar T,
176  const scalarField& c
177 ) const
178 {
179  return rk_(p, T, c);
180 }
181 
182 
183 template
184 <
185  template<class> class ReactionType,
186  class ReactionThermo,
187  class ReactionRate
188 >
190 <
191  ReactionType,
192  ReactionThermo,
193  ReactionRate
194 >::write
195 (
196  Ostream& os
197 ) const
198 {
200 
201  os.beginBlock("forward");
202  fk_.write(os);
203  os.endBlock();
204 
205  os.beginBlock("reverse");
206  rk_.write(os);
207  os.endBlock();
208 }
209 
210 
211 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::NonEquilibriumReversibleReaction
Simple extension of Reaction to handle reversible reactions using equilibrium thermodynamics.
Definition: NonEquilibriumReversibleReaction.H:61
Foam::Ostream::beginBlock
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:91
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
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::Ostream::endBlock
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:109
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
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
reaction
CombustionModel< rhoReactionThermo > & reaction
Definition: setRegionFluidFields.H:3
NonEquilibriumReversibleReaction.H
Foam::HashPtrTable
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
Definition: HashPtrTable.H:54
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56