KASKADE 7 development version
Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Kaskade::MembraneModelBase< Derived, n > Struct Template Reference

Convenience base class for membrane models providing numerical differentiation. More...

#include <membraneModels.hh>

Detailed Description

template<class Derived, int n>
struct Kaskade::MembraneModelBase< Derived, n >

Convenience base class for membrane models providing numerical differentiation.

Definition at line 56 of file membraneModels.hh.

Public Types

typedef Dune::FieldVector< double, nGatingGating
 Vector type holding gating variables. More...
 
typedef Dune::FieldMatrix< double, nGating, nGatingGatingJacobian
 Matrix type for the Jacobian of the gating dynamics. More...
 

Public Member Functions

 MembraneModelBase ()
 Default constructor. More...
 
 MembraneModelBase (std::string const &name)
 Constructor specifying the model data. 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...
 
Gating current_dv (double u, Gating const &v) const
 The transmembrane ion current derivative w.r.t. the 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...
 

Static Public Attributes

static int const nGating = n
 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

template<class Derived , int n>
typedef Dune::FieldVector<double,nGating> Kaskade::MembraneModelBase< Derived, n >::Gating

Vector type holding gating variables.

Definition at line 66 of file membraneModels.hh.

◆ GatingJacobian

template<class Derived , int n>
typedef Dune::FieldMatrix<double,nGating,nGating> Kaskade::MembraneModelBase< Derived, n >::GatingJacobian

Matrix type for the Jacobian of the gating dynamics.

Definition at line 71 of file membraneModels.hh.

Constructor & Destructor Documentation

◆ MembraneModelBase() [1/2]

template<class Derived , int n>
Kaskade::MembraneModelBase< Derived, n >::MembraneModelBase ( )
inline

Default constructor.

Name and resting state are default initialized, i.e. set to "" and 0.0, respectively.

Definition at line 78 of file membraneModels.hh.

◆ MembraneModelBase() [2/2]

template<class Derived , int n>
Kaskade::MembraneModelBase< Derived, n >::MembraneModelBase ( std::string const &  name)
inline

Constructor specifying the model data.

Parameters
namethe human-readable model name
uFixthe resting state of the transmembrane voltage
vFixthe resting state of the gating variables

Definition at line 86 of file membraneModels.hh.

Member Function Documentation

◆ current_du()

template<class Derived , int n>
double Kaskade::MembraneModelBase< Derived, n >::current_du ( double  u,
Gating const &  v,
double  h = 1e-5 
) const
inline

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.

Referenced by Kaskade::MembraneModelBase< Derived, n >::setRestState().

◆ current_dv()

template<class Derived , int n>
Gating Kaskade::MembraneModelBase< Derived, n >::current_dv ( double  u,
Gating const &  v 
) const
inline

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

This convenience method performs numerical differentiation and should be overwritten in derived classes by a better implementation.

Parameters
utransmembrane voltage
vgating variables

Definition at line 126 of file membraneModels.hh.

◆ gatingRhs_du()

template<class Derived , int n>
Gating Kaskade::MembraneModelBase< Derived, n >::gatingRhs_du ( double  u,
Gating const &  v 
) const
inline

The derivative of the right hand side for the evolution of gating variables w.r.t. the transmembrane voltage.

This default implementation uses numerical differentiation for computing the derivative.

Parameters
utransmembrane voltage
vgating variables

Definition at line 146 of file membraneModels.hh.

◆ gatingRhs_dv()

template<class Derived , int n>
GatingJacobian Kaskade::MembraneModelBase< Derived, n >::gatingRhs_dv ( double  u,
Gating const &  v 
) const
inline

The derivative of the right hand side for the evolution of gating variables w.r.t. the gating variables.

Parameters
utransmembrane voltage
vgating variables

Definition at line 158 of file membraneModels.hh.

Referenced by Kaskade::MembraneModelBase< Derived, n >::setRestState().

◆ name()

template<class Derived , int n>
std::string const & Kaskade::MembraneModelBase< Derived, n >::name ( ) const
inline

Human-readable name of the membrane model.

Definition at line 93 of file membraneModels.hh.

◆ restState()

template<class Derived , int n>
std::pair< double, Gating > const & Kaskade::MembraneModelBase< Derived, n >::restState ( ) const
inline

Value of the resting state fixed point.

Definition at line 98 of file membraneModels.hh.

◆ setRestState()

template<class Derived , int n>
void Kaskade::MembraneModelBase< Derived, n >::setRestState ( double  uFix,
Gating const &  vFix 
)
inlineprotected

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

template<class Derived , int n>
int const Kaskade::MembraneModelBase< Derived, n >::nGating = n
static

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