Loading...
Searching...
No Matches
ReversibleReaction.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 Copyright (C) 2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "ReversibleReaction.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33template
34<
35 template<class> class ReactionType,
36 class ReactionThermo,
37 class ReactionRate
38>
41(
42 const ReactionType<ReactionThermo>& reaction,
43 const ReactionRate& k
44)
45:
46 ReactionType<ReactionThermo>(reaction),
47 k_(k)
48{}
49
51template
52<
53 template<class> class ReactionType,
54 class ReactionThermo,
55 class ReactionRate
56>
59(
60 const speciesTable& species,
61 const ReactionTable<ReactionThermo>& thermoDatabase,
62 const dictionary& dict
63)
64:
65 ReactionType<ReactionThermo>(species, thermoDatabase, dict),
66 k_(species, dict)
67{}
68
70template
71<
72 template<class> class ReactionType,
73 class ReactionThermo,
74 class ReactionRate
75>
78(
81)
82:
83 ReactionType<ReactionThermo>(rr, species),
84 k_(rr.k_)
85{}
86
87
88// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89
90template
91<
92 template<class> class ReactionType,
93 class ReactionThermo,
94 class ReactionRate
95>
96Foam::scalar Foam::ReversibleReaction
97<
98 ReactionType,
99 ReactionThermo,
100 ReactionRate
101>::kf
102(
103 const scalar p,
104 const scalar T,
105 const scalarField& c
106) const
107{
108 return k_(p, T, c);
109}
110
111
112template
113<
114 template<class> class ReactionType,
115 class ReactionThermo,
116 class ReactionRate
117>
118Foam::scalar Foam::ReversibleReaction
119<
120 ReactionType,
121 ReactionThermo,
122 ReactionRate
123>::kr
124(
125 const scalar kfwd,
126 const scalar p,
127 const scalar T,
128 const scalarField& c
129) const
130{
131 return kfwd/max(this->Kc(p, T), VSMALL);
132}
133
134
135template
136<
137 template<class> class ReactionType,
138 class ReactionThermo,
139 class ReactionRate
140>
141Foam::scalar Foam::ReversibleReaction
142<
143 ReactionType,
144 ReactionThermo,
145 ReactionRate
146>::kr
147(
148 const scalar p,
149 const scalar T,
150 const scalarField& c
151) const
152{
153 return kr(kf(p, T, c), p, T, c);
154}
155
156
157template
158<
159 template<class> class ReactionType,
160 class ReactionThermo,
161 class ReactionRate
162>
164<
165 ReactionType,
166 ReactionThermo,
167 ReactionRate
168>::write
169(
170 Ostream& os
171) const
172{
174 k_.write(os);
175}
176
177
178// ************************************************************************* //
label k
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
virtual void write(Ostream &os) const
Write.
Definition Reaction.C:382
Simple extension of Reaction to handle reversible reactions using equilibrium thermodynamics.
virtual scalar kr(const scalar kfwd, const scalar p, const scalar T, const scalarField &c) const
Reverse rate constant from the given formard rate constant.
ReversibleReaction(const ReactionType< ReactionThermo > &reaction, const ReactionRate &k)
Construct from components.
virtual scalar kf(const scalar p, const scalar T, const scalarField &c) const
Forward rate constant.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
CombustionModel< rhoReactionThermo > & reaction
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
hashedWordList speciesTable
A table of species as a hashedWordList.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
HashPtrTable< ThermoType > ReactionTable
Definition Reaction.H:52
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
runTime write()
dictionary dict