MichaelisMentenReactionRate.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-2021 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 Class
25  Foam::MichaelisMentenReactionRate
26 
27 Description
28  Michaelis-Menten reaction rate for enzymatic reactions.
29 
30  Reference:
31  \verbatim
32  Michaelis, L., & Menten, M. L. (1913).
33  Die Kinetik der InwertinWirkung.
34  Biochem., (49), 333-369.
35  \endverbatim
36 
37 SourceFiles
38  MichaelisMentenReactionRateI.H
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef MichaelisMentenReactionRate_H
43 #define MichaelisMentenReactionRate_H
44 
45 #include "Reaction.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // Forward declaration of friend functions and operators
53 
54 class MichaelisMentenReactionRate;
55 
56 Ostream& operator<<(Ostream&, const MichaelisMentenReactionRate&);
57 
58 
59 /*---------------------------------------------------------------------------*\
60  Class MichaelisMentenReactionRate Declaration
61 \*---------------------------------------------------------------------------*/
62 
64 {
65  // Private Data
66 
67  //- List of species present in reaction system
68  const speciesTable& species_;
69 
70  //- The maximum reaction rate at saturating substrate concentration
71  scalar Vmax_;
72 
73  //- The Michaelis constant
74  // the substrate concentration at which the reaction rate is half Vmax_
75  scalar Km_;
76 
77  //- The substrate specie index
78  label s_;
79 
81 
82 
83 public:
84 
85  // Constructors
86 
87  //- Construct from dictionary
89  (
90  const speciesTable& species,
91  const dictionary& dict
92  );
93 
94 
95  // Member Functions
96 
97  //- Return the type name
98  static word type()
99  {
100  return "MichaelisMenten";
101  }
102 
103  //- Pre-evaluation hook
104  inline void preEvaluate() const;
105 
106  //- Post-evaluation hook
107  inline void postEvaluate() const;
108 
109  inline scalar operator()
110  (
111  const scalar p,
112  const scalar T,
113  const scalarField& c,
114  const label li
115  ) const;
116 
117  inline scalar ddT
118  (
119  const scalar p,
120  const scalar T,
121  const scalarField& c,
122  const label li
123  ) const;
124 
125  //- Third-body efficiencies (beta = 1-alpha)
126  // non-empty only for third-body reactions
127  // with enhanced molecularity (alpha != 1)
128  inline const List<Tuple2<label, scalar>>& beta() const;
129 
130  //- Species concentration derivative of the pressure dependent term
131  inline void dcidc
132  (
133  const scalar p,
134  const scalar T,
135  const scalarField& c,
136  const label li,
138  ) const;
139 
140  //- Temperature derivative of the pressure dependent term
141  inline scalar dcidT
142  (
143  const scalar p,
144  const scalar T,
145  const scalarField& c,
146  const label li
147  ) const;
148 
149  //- Write to stream
150  inline void write(Ostream& os) const;
151 
152 
153  // Ostream Operator
154 
155  inline friend Ostream& operator<<
156  (
157  Ostream&,
159  );
160 };
161 
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 } // End namespace Foam
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 #endif
174 
175 // ************************************************************************* //
dictionary dict
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
scalar dcidT(const scalar p, const scalar T, const scalarField &c, const label li) const
Temperature derivative of the pressure dependent term.
Michaelis-Menten reaction rate for enzymatic reactions.
void preEvaluate() const
Pre-evaluation hook.
static word type()
Return the type name.
const dimensionedScalar c
Speed of light in a vacuum.
void postEvaluate() const
Post-evaluation hook.
A class for handling words, derived from string.
Definition: word.H:59
void dcidc(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
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:54
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void write(Ostream &os) const
Write to stream.
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
A wordList with hashed indices for faster lookup by name.
Ostream & operator<<(Ostream &, const ensightPart &)
MichaelisMentenReactionRate(const speciesTable &species, const dictionary &dict)
Construct from dictionary.
volScalarField & p
Namespace for OpenFOAM.