KASKADE 7 development version
|
Phenomenological model by Aliev and Panfilov. More...
#include <membraneModels.hh>
Phenomenological model by Aliev and Panfilov.
A phenomenological model with just one gating variables. This is supposed to be used in the monodomain-like system
\[ \begin{aligned} \dot u &= \mathrm{div}(\sigma\nabla u) - k u(u-a)(u-1) - uv \\ \dot v &= \epsilon(u,v)(-v - ku(u-1-a)) \end{aligned} \]
with \( \epsilon(u,v) = \epsilon_0 + \mu_1 v / (u+\mu_2) \).
The (dimensionless) transmembrane voltage \( u \) covers a range [0,1], and can be translated back to physiological values \( U \) by \( U = u/10 - 0.08 \) [V]. Similarly, the time \( t \) in the model is dimensionless, and can be translated back to real time \( T \) by \( T = 0.0129 t \) [s].
[1] Aliev, Panfilov. Chaos, Solitions, and Fractals 7(3):293-301, 1996.
Definition at line 280 of file membraneModels.hh.
Public Types | |
typedef Dune::FieldVector< double, nGating > | Gating |
typedef Dune::FieldMatrix< double, nGating, nGating > | GatingJacobian |
Public Member Functions | |
AlievPanfilov () | |
Default 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 = 1 |
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::AlievPanfilov::Gating |
Definition at line 287 of file membraneModels.hh.
typedef Dune::FieldMatrix<double,nGating,nGating> Kaskade::AlievPanfilov::GatingJacobian |
Definition at line 288 of file membraneModels.hh.
Kaskade::AlievPanfilov::AlievPanfilov | ( | ) |
Default constructor.
double Kaskade::AlievPanfilov::current | ( | double | u, |
Gating const & | v | ||
) | const |
The transmembrane ion current.
u | transmembrane voltage |
v | gating variables |
double Kaskade::AlievPanfilov::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::AlievPanfilov::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 285 of file membraneModels.hh.