BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / thermophysicalModels / specie / reaction / Reactions / NonEquilibriumReversibleReaction / NonEquilibriumReversibleReaction.C
blob6c90e7c0823e4e690c89babcb3b30a54ec197a63
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "NonEquilibriumReversibleReaction.H"
28 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
30 template<class ReactionThermo, class ReactionRate>
31 Foam::NonEquilibriumReversibleReaction<ReactionThermo, ReactionRate>::
32 NonEquilibriumReversibleReaction
34     const Reaction<ReactionThermo>& reaction,
35     const ReactionRate& forwardReactionRate,
36     const ReactionRate& reverseReactionRate
39     Reaction<ReactionThermo>(reaction),
40     fk_(forwardReactionRate),
41     rk_(reverseReactionRate)
46 template<class ReactionThermo, class ReactionRate>
47 Foam::NonEquilibriumReversibleReaction<ReactionThermo, ReactionRate>::
48 NonEquilibriumReversibleReaction
50     const speciesTable& species,
51     const HashPtrTable<ReactionThermo>& thermoDatabase,
52     Istream& is
55     Reaction<ReactionThermo>(species, thermoDatabase, is),
56     fk_(species, is),
57     rk_(species, is)
61 template<class ReactionThermo, class ReactionRate>
62 Foam::NonEquilibriumReversibleReaction<ReactionThermo, ReactionRate>::
63 NonEquilibriumReversibleReaction
65     const speciesTable& species,
66     const HashPtrTable<ReactionThermo>& thermoDatabase,
67     const dictionary& dict
70     Reaction<ReactionThermo>(species, thermoDatabase, dict),
71     fk_(species, dict.subDict("forward")),
72     rk_(species, dict.subDict("reverse"))
76 template<class ReactionThermo, class ReactionRate>
77 Foam::NonEquilibriumReversibleReaction<ReactionThermo, ReactionRate>::
78 NonEquilibriumReversibleReaction
80     const NonEquilibriumReversibleReaction<ReactionThermo, ReactionRate>& nerr,
81     const speciesTable& species
84     Reaction<ReactionThermo>(nerr, species),
85     fk_(nerr.fk_),
86     rk_(nerr.rk_)
90 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
92 template<class ReactionThermo, class ReactionRate>
93 Foam::scalar
94 Foam::NonEquilibriumReversibleReaction<ReactionThermo, ReactionRate>::kf
96     const scalar T,
97     const scalar p,
98     const scalarField& c
99 ) const
101     return fk_(T, p, c);
105 template<class ReactionThermo, class ReactionRate>
106 Foam::scalar
107 Foam::NonEquilibriumReversibleReaction<ReactionThermo, ReactionRate>::kr
109     const scalar,
110     const scalar T,
111     const scalar p,
112     const scalarField& c
113 ) const
115     return rk_(T, p, c);
119 template<class ReactionThermo, class ReactionRate>
120 Foam::scalar
121 Foam::NonEquilibriumReversibleReaction<ReactionThermo, ReactionRate>::kr
123     const scalar T,
124     const scalar p,
125     const scalarField& c
126 ) const
128     return rk_(T, p, c);
132 template<class ReactionThermo, class ReactionRate>
133 void Foam::NonEquilibriumReversibleReaction<ReactionThermo, ReactionRate>::write
135     Ostream& os
136 ) const
138     Reaction<ReactionThermo>::write(os);
140     os  << indent << "forward" << nl;
141     os  << indent << token::BEGIN_BLOCK << nl;
142     os  << incrIndent;
143     fk_.write(os);
144     os  << decrIndent;
145     os  << indent << token::END_BLOCK << nl;
147     os  << indent << "reverse" << nl;
148     os  << indent << token::BEGIN_BLOCK << nl;
149     os  << incrIndent;
150     rk_.write(os);
151     os  << decrIndent;
152     os  << indent << token::END_BLOCK << nl;
156 // ************************************************************************* //