MichaelisMentenReactionRateI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2018 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 
29 (
30  const speciesTable& st,
31  const dictionary& dict
32 )
33 :
34  species_(st),
35  Vmax_(readScalar(dict.lookup("Vmax"))),
36  Km_(readScalar(dict.lookup("Km"))),
37  s_(st[dict.lookup("S")])
38 {
39  beta_.append(Tuple2<label, scalar>(s_, 1.0));
40 }
41 
42 
43 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
44 
45 inline Foam::scalar Foam::MichaelisMentenReactionRate::operator()
46 (
47  const scalar p,
48  const scalar T,
49  const scalarField& c
50 ) const
51 {
52  return Vmax_/(Km_ + c[s_]);
53 }
54 
55 
57 (
58  const scalar p,
59  const scalar T,
60  const scalarField& c
61 ) const
62 {
63  return 0;
64 }
65 
66 
69 {
70  return beta_;
71 }
72 
73 
75 (
76  const scalar p,
77  const scalar T,
78  const scalarField& c,
79  scalarField& dcidc
80 ) const
81 {
82  dcidc[0] = -1.0/(Km_ + c[s_]);
83 }
84 
85 
87 (
88  const scalar p,
89  const scalar T,
90  const scalarField& c
91 ) const
92 {
93  return 0;
94 }
95 
96 
98 {
99  os.writeKeyword("Vmax") << Vmax_ << token::END_STATEMENT << nl;
100  os.writeKeyword("Km") << Km_ << token::END_STATEMENT << nl;
101  os.writeKeyword("S") << species_[s_] << token::END_STATEMENT << nl;
102 }
103 
104 
105 inline Foam::Ostream& Foam::operator<<
106 (
107  Ostream& os,
108  const MichaelisMentenReactionRate& mmrr
109 )
110 {
111  mmrr.write(os);
112  return os;
113 }
114 
115 
116 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A 2-tuple for storing two objects of different types.
Definition: HashTable.H:66
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
Michaelis-Menten reaction rate for enzymatic reactions.
void dcidc(const scalar p, const scalar T, const scalarField &c, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
scalar ddT(const scalar p, const scalar T, const scalarField &c) const
const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
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:265
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
const volScalarField & T
scalar dcidT(const scalar p, const scalar T, const scalarField &c) const
Temperature derivative of the pressure dependent term.
void write(Ostream &os) const
Write to stream.
A wordList with hashed indices for faster lookup by name.
const dimensionedScalar c
Speed of light in a vacuum.
MichaelisMentenReactionRate(const speciesTable &species, const dictionary &dict)
Construct from dictionary.
virtual Ostream & write(const token &)=0
Write next token to stream.
volScalarField & p
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576