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-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 \*---------------------------------------------------------------------------*/
25 
26 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
29 (
30  const speciesTable& st,
31  const dictionary& dict
32 )
33 :
34  species_(st),
35  Vmax_(dict.lookup<scalar>("Vmax")),
36  Km_(dict.lookup<scalar>("Km")),
37  s_(st[dict.lookup("S")])
38 {
39  beta_.append(Tuple2<label, scalar>(s_, 1.0));
40 }
41 
42 
43 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
44 
46 {}
47 
48 
50 {}
51 
52 
53 inline Foam::scalar Foam::MichaelisMentenReactionRate::operator()
54 (
55  const scalar p,
56  const scalar T,
57  const scalarField& c,
58  const label li
59 ) const
60 {
61  return Vmax_/(Km_ + c[s_]);
62 }
63 
64 
66 (
67  const scalar p,
68  const scalar T,
69  const scalarField& c,
70  const label li
71 ) const
72 {
73  return 0;
74 }
75 
76 
79 {
80  return beta_;
81 }
82 
83 
85 (
86  const scalar p,
87  const scalar T,
88  const scalarField& c,
89  const label li,
90  scalarField& dcidc
91 ) const
92 {
93  dcidc[0] = -1.0/(Km_ + c[s_]);
94 }
95 
96 
98 (
99  const scalar p,
100  const scalar T,
101  const scalarField& c,
102  const label li
103 ) const
104 {
105  return 0;
106 }
107 
108 
110 {
111  writeEntry(os, "Vmax", Vmax_);
112  writeEntry(os, "Km", Km_);
113  writeEntry(os, "S", species_[s_]);
114 }
115 
116 
117 inline Foam::Ostream& Foam::operator<<
118 (
119  Ostream& os,
120  const MichaelisMentenReactionRate& mmrr
121 )
122 {
123  mmrr.write(os);
124  return os;
125 }
126 
127 
128 // ************************************************************************* //
virtual Ostream & write(const char)=0
Write character.
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 2-tuple for storing two objects of different types.
Definition: HashTable.H:65
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.
const dimensionedScalar c
Speed of light in a vacuum.
void postEvaluate() const
Post-evaluation hook.
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 writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
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.
MichaelisMentenReactionRate(const speciesTable &species, const dictionary &dict)
Construct from dictionary.
volScalarField & p
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844