KASKADE 7 development version
Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Kaskade::PhysiologicalAlievPanfilov Struct Reference

Phenomenological model based on Aliev and Panfilov, scaled to physiological values. More...

#include <membraneModels.hh>

Detailed Description

Phenomenological model based on Aliev and Panfilov, scaled to physiological values.

The AlievPanfilov model uses dimensionless quantities that need to be scaled to physiological ranges. Moreover, the model neglects the membrane capacity.

Here, we rescale the model to be used in the monodomain system

\[ \begin{aligned} \chi \left( C_m \dot U + I_{\rm ion}(U,V)\right) &= \mathrm{div}(\sigma\nabla U) \\ \dot V &= g(U,V) \end{aligned}, \]

while reproducing the original behaviour designed for the simplified system

\[ \begin{aligned} \frac{d}{d\tau} u &= \mathrm{div}(\sigma_{\rm AP}\nabla u) + I_{\rm AP} \\ \frac{d}{d\tau} \dot v &= g_{\rm AP}. \end{aligned} \]

By inserting \( U = u/10 - 0.08 \) and \( t = s \tau \) with \( s = 0.0129 \) as suggested by Aliev and Panfilov into the monodomain system, we obtain

\[ \begin{aligned} \frac{d}{d\tau} u &= \mathrm{div}\left(\frac{s\sigma}{\chi C_m}\nabla u\right) - \frac{10 s}{C_m} I_{\rm ion} \\ \frac{d}{d\tau} v &= s g \end{aligned}. \]

Comparing the right hand sides leads to

\[ \begin{aligned} I_{\rm ion}(U,V) &= -\frac{C_m}{10 s} I_{\rm AP}(u,v) \\ g(U,V) &= \frac{1}{s} g_{\rm AP}(u,v). \end{aligned} \]

The action potential duration is around 310 ms. With standard values for \( C_m, \chi \), and \( \sigma \) from [1], Tab. 10.1 a conduction velocity around 13 cm/s is obtained, about a factor 4 too slow for myocardial CV. We thus scale the membrane capacitance default down by a factor of 16 to obtain a CV of 50 cm/s.

References

[1] K. Horgmo-Jæger, A. Tveito. Differential Equations for Studies in Computational Electrophysiology. Springer 2023.

Definition at line 378 of file membraneModels.hh.

Inheritance diagram for Kaskade::PhysiologicalAlievPanfilov:
Kaskade::MembraneModelBase< PhysiologicalAlievPanfilov, 1 >

Public Types

typedef Dune::FieldVector< double, nGatingGating
 
typedef Dune::FieldMatrix< double, nGating, nGatingGatingJacobian
 

Public Member Functions

 PhysiologicalAlievPanfilov (double Cm=1e-2/16)
 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...
 

Member Typedef Documentation

◆ Gating

Definition at line 385 of file membraneModels.hh.

◆ GatingJacobian

Definition at line 386 of file membraneModels.hh.

Constructor & Destructor Documentation

◆ PhysiologicalAlievPanfilov()

Kaskade::PhysiologicalAlievPanfilov::PhysiologicalAlievPanfilov ( double  Cm = 1e-2/16)

Constructor.

Parameters
Cmthe membrane capacitance density used for scaling the ion current (in F/m^2). The default value from [1] is scaled down by a factor of 16 to obtain a conduction velocity of 50cm/s.

Member Function Documentation

◆ current()

double Kaskade::PhysiologicalAlievPanfilov::current ( double  u,
Gating const &  v 
) const

The transmembrane ion current.

Parameters
utransmembrane voltage
vgating variables

◆ current_du() [1/2]

double Kaskade::PhysiologicalAlievPanfilov::current_du ( double  u,
Gating const &  v 
) const

The transmembrane ion current derivative w.r.t. the transmembrane voltage.

Parameters
utransmembrane voltage
vgating variables

◆ current_du() [2/2]

double Kaskade::MembraneModelBase< PhysiologicalAlievPanfilov , n >::current_du ( double  u,
Gating const &  v,
double  h = 1e-5 
) const
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.

Parameters
utransmembrane voltage [V]
vgating variables
hstep size for numerical differentiation

Definition at line 109 of file membraneModels.hh.

◆ current_dv()

Gating Kaskade::PhysiologicalAlievPanfilov::current_dv ( double  u,
Gating const &  v 
) const

The transmembrane ion current derivative w.r.t. the gating variables.

Parameters
utransmembrane voltage
vgating variables

◆ gatingRhs()

Gating Kaskade::PhysiologicalAlievPanfilov::gatingRhs ( double  u,
Gating const &  v,
double  Istim = 0 
) const

The right hand side for the evolution of gating variables.

Parameters
utransmembrane voltage
vgating variables

◆ gatingRhs_du()

Gating Kaskade::PhysiologicalAlievPanfilov::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.

Parameters
utransmembrane voltage
vgating variables

◆ gatingRhs_dv()

GatingJacobian Kaskade::PhysiologicalAlievPanfilov::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.

Parameters
utransmembrane voltage
vgating variables

◆ name()

std::string const & Kaskade::MembraneModelBase< PhysiologicalAlievPanfilov , n >::name ( ) const
inlineinherited

Human-readable name of the membrane model.

Definition at line 93 of file membraneModels.hh.

◆ restState()

std::pair< double, Gating > const & Kaskade::MembraneModelBase< PhysiologicalAlievPanfilov , n >::restState ( ) const
inlineinherited

Value of the resting state fixed point.

Definition at line 98 of file membraneModels.hh.

◆ setRestState()

void Kaskade::MembraneModelBase< PhysiologicalAlievPanfilov , n >::setRestState ( double  uFix,
Gating const &  vFix 
)
inlineprotectedinherited

Accepts an (approximate) resting state.

Parameters
uFixapproximate resting state of the transmembrane voltage
vFixapproximate 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.

Member Data Documentation

◆ nGating

int const Kaskade::PhysiologicalAlievPanfilov::nGating = 1
static

number of gating variables

Definition at line 383 of file membraneModels.hh.


The documentation for this struct was generated from the following file: