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-2024 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 
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 (
32  const speciesTable& st,
33  const dimensionSet& dims,
34  const dictionary& dict
35 )
36 :
37  species_(st),
38  Vmax_(dict.lookup<scalar>("Vmax", dims*dimMoles/dimVolume)),
39  Km_(dict.lookup<scalar>("Km", dimMoles/dimVolume)),
40  s_(st[dict.lookup("S")])
41 {}
42 
43 
44 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
45 
47 {}
48 
49 
51 {}
52 
53 
54 inline Foam::scalar Foam::MichaelisMentenReactionRate::operator()
55 (
56  const scalar p,
57  const scalar T,
58  const scalarField& c,
59  const label li
60 ) const
61 {
62  return Vmax_/(Km_ + c[s_]);
63 }
64 
65 
67 (
68  const scalar p,
69  const scalar T,
70  const scalarField& c,
71  const label li
72 ) const
73 {
74  return 0;
75 }
76 
77 
79 {
80  return true;
81 }
82 
83 
85 (
86  const scalar p,
87  const scalar T,
88  const scalarField& c,
89  const label li,
90  scalarField& ddc
91 ) const
92 {
93  ddc = 0;
94  ddc[s_] = - Vmax_/sqr(Km_ + c[s_]);
95 }
96 
97 
99 {
100  writeEntry(os, "Vmax", Vmax_);
101  writeEntry(os, "Km", Km_);
102  writeEntry(os, "S", species_[s_]);
103 }
104 
105 
106 inline Foam::Ostream& Foam::operator<<
107 (
108  Ostream& os,
109  const MichaelisMentenReactionRate& mmrr
110 )
111 {
112  mmrr.write(os);
113  return os;
114 }
115 
116 
117 // ************************************************************************* //
Michaelis-Menten reaction rate for enzymatic reactions.
void preEvaluate() const
Pre-evaluation hook.
void postEvaluate() const
Post-evaluation hook.
MichaelisMentenReactionRate(const speciesTable &species, const dimensionSet &dims, const dictionary &dict)
Construct from dictionary.
void write(Ostream &os) const
Write to stream.
void ddc(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &ddc) const
The derivative of the rate w.r.t. concentration.
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
The derivative of the rate w.r.t. temperature.
bool hasDdc() const
Is the rate a function of concentration?
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual Ostream & write(const char)=0
Write character.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
A wordList with hashed indices for faster lookup by name.
const dimensionedScalar c
Speed of light in a vacuum.
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimVolume
const dimensionSet dimMoles
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p