KASKADE 7 development version
|
The 3-variable Fenton-Karma membrane model. More...
#include <membraneModels.hh>
The 3-variable Fenton-Karma membrane model.
[1] F. Fenton, A. Karma. Vortex dynamics in three-dimensional continuous myocardium with fiber rotation: Filament instability and fibrillation. Chaos 8: 20-47, 1998.
Definition at line 453 of file membraneModels.hh.
Public Types | |
enum | ParameterSet { BR , MBR , MLRI , GP } |
Different parameter sets provided. More... | |
typedef Dune::FieldVector< double, nGating > | Gating |
typedef Dune::FieldMatrix< double, nGating, nGating > | GatingJacobian |
Public Member Functions | |
FentonKarma (ParameterSet p, double Cm=1e-2) | |
Constructor. More... | |
double | current (double u, Gating const &v) const |
The transmembrane ion current. More... | |
double | current_du (double u, Gating const &v) const |
The transmembrane ion current derivative w.r.t. the transmembrane voltage. More... | |
Gating | current_dv (double u, Gating const &v) const |
The transmembrane ion current derivative w.r.t. the gating variables. More... | |
Gating | gatingRhs (double u, Gating const &v, double Istim=0) const |
The right hand side for the evolution of gating variables. More... | |
Gating | gatingRhs_du (double u, Gating const &v) const |
The derivative of the right hand side for the evolution of gating variables w.r.t. the transmembrane voltage. More... | |
GatingJacobian | gatingRhs_dv (double u, Gating const &v) const |
The derivative of the right hand side for the evolution of gating variables w.r.t. the gating variables. More... | |
std::string const & | name () const |
Human-readable name of the membrane model. More... | |
std::pair< double, Gating > const & | restState () const |
Value of the resting state fixed point. More... | |
double | current_du (double u, Gating const &v, double h=1e-5) const |
The transmembrane ion current derivative w.r.t. the transmembrane voltage. More... | |
Static Public Attributes | |
static int const | nGating = 2 |
number of gating variables More... | |
Protected Member Functions | |
void | setRestState (double uFix, Gating const &vFix) |
Accepts an (approximate) resting state. More... | |
typedef Dune::FieldVector<double,nGating> Kaskade::FentonKarma::Gating |
Definition at line 474 of file membraneModels.hh.
typedef Dune::FieldMatrix<double,nGating,nGating> Kaskade::FentonKarma::GatingJacobian |
Definition at line 475 of file membraneModels.hh.
Different parameter sets provided.
The original article [1] provides four different parameter sets for approximating different (more complex) membrane models. We allow to choose between them based on the constructor argument.
Enumerator | |
---|---|
BR | Beeler-Reuter model. |
MBR | modified Beeler-Reuter |
MLRI | modified Luo-Rudy I |
GP | guinea pig model |
Definition at line 468 of file membraneModels.hh.
Kaskade::FentonKarma::FentonKarma | ( | ParameterSet | p, |
double | Cm = 1e-2 |
||
) |
Constructor.
Cm | the membrane capacitance density used for scaling the ion current (in F/m^2). |
double Kaskade::FentonKarma::current | ( | double | u, |
Gating const & | v | ||
) | const |
The transmembrane ion current.
u | transmembrane voltage |
v | gating variables |
double Kaskade::FentonKarma::current_du | ( | double | u, |
Gating const & | v | ||
) | const |
The transmembrane ion current derivative w.r.t. the transmembrane voltage.
u | transmembrane voltage |
v | gating variables |
|
inlineinherited |
The transmembrane ion current derivative w.r.t. the transmembrane voltage.
This is a default implementation based on numerical differentiation. Overload this if a better implementation is available.
u | transmembrane voltage [V] |
v | gating variables |
h | step size for numerical differentiation |
Definition at line 109 of file membraneModels.hh.
The transmembrane ion current derivative w.r.t. the gating variables.
u | transmembrane voltage |
v | gating variables |
The right hand side for the evolution of gating variables.
u | transmembrane voltage |
v | gating variables |
The derivative of the right hand side for the evolution of gating variables w.r.t. the transmembrane voltage.
u | transmembrane voltage |
v | gating variables |
GatingJacobian Kaskade::FentonKarma::gatingRhs_dv | ( | double | u, |
Gating const & | v | ||
) | const |
The derivative of the right hand side for the evolution of gating variables w.r.t. the gating variables.
u | transmembrane voltage |
v | gating variables |
|
inlineinherited |
Human-readable name of the membrane model.
Definition at line 93 of file membraneModels.hh.
|
inlineinherited |
Value of the resting state fixed point.
Definition at line 98 of file membraneModels.hh.
|
inlineprotectedinherited |
Accepts an (approximate) resting state.
uFix | approximate resting state of the transmembrane voltage |
vFix | approximate resting state of the gating variables |
The method shall be used by derived classes to provide an approximate value of the resting state (which needs to be a stable fixed point). A few iterations of Newton's method are performed in order to obtain an accurate value.
Definition at line 181 of file membraneModels.hh.
|
static |
number of gating variables
Definition at line 459 of file membraneModels.hh.