LifeV
IonicNoblePurkinje Class Reference

IonicModel - This class implements an ionic model. More...

#include <IonicNoblePurkinje.hpp>

+ Inheritance diagram for IonicNoblePurkinje:
+ Collaboration diagram for IonicNoblePurkinje:

Public Member Functions

void computeGatingRhs (const std::vector< Real > &v, std::vector< Real > &rhs)
 Methods. More...
 
void computeRhs (const std::vector< Real > &v, std::vector< Real > &rhs)
 This methods contains the actual evaluation of the rhs of all state variablesin the model (0D version) More...
 
Real computeLocalPotentialRhs (const std::vector< Real > &v)
 This methods contains the actual evaluation of the rhs of the voltage equation only (0D version) More...
 
void computeGatingVariablesWithRushLarsen (std::vector< Real > &v, const Real dt)
 This methods contains the actual evaluation of the rhs of the voltage equation only (0D version) More...
 
void showMe ()
 Display information about the model. More...
 
- Public Member Functions inherited from ElectroIonicModel
virtual matrix_Type getJac (const vector_Type &v, Real h=1.0e-8)
 This methods computes the Jacobian numerically. More...
 
virtual std::vector< std::vector< Real > > getJac (const std::vector< Real > &v, Real h=1.0e-8)
 This methods computes the Jacobian numerically. More...
 
virtual void computeGatingRhs (const std::vector< vectorPtr_Type > &v, std::vector< vectorPtr_Type > &rhs)
 This methods computes the right hand side of the gating variables in the 3D case. More...
 
virtual void computeNonGatingRhs (const std::vector< vectorPtr_Type > &v, std::vector< vectorPtr_Type > &rhs)
 This methods computes the right hand side of the state variables that are not gating variables in the 3D case. More...
 
virtual void computeGatingVariablesWithRushLarsen (std::vector< vectorPtr_Type > &v, const Real dt)
 Compute the new value of the gating variables in 3D with the Rush Larsen method specified in the 0D version of the ionic model. More...
 
virtual void computeRhs (const std::vector< vectorPtr_Type > &v, std::vector< vectorPtr_Type > &rhs)
 Compute the right hand side of the ionic model in 3D. More...
 
virtual void computePotentialRhsICI (const std::vector< vectorPtr_Type > &v, std::vector< vectorPtr_Type > &rhs, matrix_Type &massMatrix)
 Compute the right hand side of the voltage equation linearly interpolating the ionic currents. More...
 
virtual void computePotentialRhsSVI (const std::vector< vectorPtr_Type > &v, std::vector< vectorPtr_Type > &rhs, FESpace< mesh_Type, MapEpetra > &uFESpace)
 Compute the right hand side of the voltage equation using SVI. More...
 
virtual void computePotentialRhsSVI (const std::vector< vectorPtr_Type > &v, std::vector< vectorPtr_Type > &rhs, FESpace< mesh_Type, MapEpetra > &uFESpace, const QuadratureRule &qr)
 Compute the right hand side of the voltage equation using SVI specifying the quadrature rule. More...
 
virtual void initialize (std::vector< Real > &v)
 Initialize the ionic model with a given vector of state variable (0D version) More...
 
virtual void initialize (std::vector< vectorPtr_Type > &v)
 Initialize the ionic model in 3D with a given vector of state variable vector pointers. More...
 
virtual void computeNonGatingRhs (const std::vector< Real > &v, std::vector< Real > &rhs)
 In the case this method is improperly used, it should use this default implementation. More...
 
 ElectroIonicModel ()
 Empty Constructor. More...
 
 ElectroIonicModel (int n)
 Constructor. More...
 
 ElectroIonicModel (int n, int g)
 Constructor. More...
 
 ElectroIonicModel (const ElectroIonicModel &Ionic)
 Copy Constructor. More...
 
virtual ~ElectroIonicModel ()
 Destructor. More...
 
virtual void setup (Teuchos::ParameterList &parameterList)
 setup the parameters of the ionic model from an xml file More...
 
const short int Size () const
 returns the number of equations of the ionic model More...
 
const short int numberOfGatingVariables () const
 returns the number of gating variables in the ionic model More...
 
const Real membraneCapacitance () const
 returns the value of the membrane capacitance in the model More...
 
const Real appliedCurrent () const
 returns the value of the applied current in the model/point More...
 
vectorPtr_Type appliedCurrentPtr ()
 returns the pointer to the applied current FE vector in the 3D case More...
 
const std::vector< RealrestingConditions () const
 returns the vector with the resting values of the variables in the ionic model More...
 
const function_Type pacaingProtocol () const
 returns the function describing the pacing protocol for the ionic model More...
 
void setMembraneCapacitance (const Real p)
 set the membrane capacitance in the ionic model More...
 
void setAppliedCurrent (const Real p)
 set the applied current in the ionic model/point More...
 
void setAppliedCurrentPtr (const vectorPtr_Type p)
 set the pointer to the applied current in the 3D ionic model More...
 
void setAppliedCurrent (const vector_Type &p)
 set the pointer to the applied current in the 3D ionic model More...
 
void setAppliedCurrentFromFunction (function_Type &f, feSpacePtr_Type feSpacePtr, Real time=0.0)
 Interpolate the function f on the FE space feSpacePtr at time time. More...
 
void setAppliedCurrentFromElectroStimulus (ElectroStimulus &stimulus, feSpacePtr_Type feSpacePtr, Real time=0.0)
 Interpolate the function of the electro stimulus. More...
 
void setPacingProtocol (function_Type pacingProtocol)
 Set the pacing protocol as boost function. More...
 
void setRestingCondtions (Real value, int j)
 Set component of the resting conditions. More...
 
void setRestingCondtions (std::vector< Real > &restingConditions)
 Set resting conditions. More...
 
void addAppliedCurrent (Real &rhs)
 Simple wrapper to add the applied current. More...
 
void addAppliedCurrent (std::vector< Real > &rhs)
 Simple wrapper to add the applied current. More...
 
ElectroIonicModeloperator= (const ElectroIonicModel &Ionic)
 Assignment operator. More...
 

Private Attributes

Real M_gi
 Solves the ionic model. More...
 
Real M_vNa
 
Real M_vK
 
Real M_Cm
 
Real M_Itot
 

Type definitions

typedef ElectroIonicModel super
 
typedef std::shared_ptr< VectorEpetravectorPtr_Type
 
typedef std::shared_ptr< VectorElementalelvecPtr_Type
 
typedef RegionMesh< LinearTetramesh_Type
 

Constructors & Destructor

 IonicNoblePurkinje ()
 Constructor. More...
 
 IonicNoblePurkinje (Teuchos::ParameterList &parameterList)
 
 IonicNoblePurkinje (const IonicNoblePurkinje &model)
 
virtual ~IonicNoblePurkinje ()
 Destructor. More...
 

Overloads

IonicNoblePurkinjeoperator= (const IonicNoblePurkinje &model)
 Operator. More...
 

Setters and getters

const Realgi () const
 
const RealCm () const
 
const RealvK () const
 
const RealvNa () const
 
const RealItotal () const
 
void setgi (const Real &p)
 
void setCm (const Real &p)
 
void setvNa (const Real &p)
 
void setvK (const Real &p)
 
static Real GeneralFunctionAlphaAndBeta (const Real &V, const Real &C1, const Real &C2, const Real &C3, const Real &C4, const Real &C5, const Real &V_0)
 
static Real mInf (const Real &V)
 
static Real hInf (const Real &V)
 
static Real nInf (const Real &V)
 

Additional Inherited Members

- Public Types inherited from ElectroIonicModel
typedef VectorEpetra vector_Type
 
typedef std::shared_ptr< VectorEpetravectorPtr_Type
 
typedef std::shared_ptr< VectorElementalelvecPtr_Type
 
typedef RegionMesh< LinearTetramesh_Type
 
typedef MatrixEpetra< Realmatrix_Type
 
typedef std::shared_ptr< matrix_TypematrixPtr_Type
 
typedef FESpace< mesh_Type, MapEpetrafeSpace_Type
 
typedef std::shared_ptr< feSpace_TypefeSpacePtr_Type
 
typedef std::function< Real(const Real &t, const Real &x, const Real &y, const Real &z, const ID &i) > function_Type
 
typedef FactorySingleton< Factory< ElectroIonicModel, std::string > > IonicModelFactory
 
- Protected Attributes inherited from ElectroIonicModel
short int M_numberOfEquations
 
short int M_numberOfGatingVariables
 
std::vector< RealM_restingConditions
 
Real M_membraneCapacitance
 
Real M_appliedCurrent
 
vectorPtr_Type M_appliedCurrentPtr
 
function_Type M_pacingProtocol
 

Detailed Description

IonicModel - This class implements an ionic model.

Definition at line 58 of file IonicNoblePurkinje.hpp.

Member Typedef Documentation

◆ super

Definition at line 64 of file IonicNoblePurkinje.hpp.

◆ vectorPtr_Type

typedef std::shared_ptr<VectorEpetra> vectorPtr_Type

Definition at line 65 of file IonicNoblePurkinje.hpp.

◆ elvecPtr_Type

typedef std::shared_ptr<VectorElemental> elvecPtr_Type

Definition at line 66 of file IonicNoblePurkinje.hpp.

◆ mesh_Type

Definition at line 67 of file IonicNoblePurkinje.hpp.

Constructor & Destructor Documentation

◆ IonicNoblePurkinje() [1/3]

Constructor.

Constructors.

Definition at line 50 of file IonicNoblePurkinje.cpp.

◆ IonicNoblePurkinje() [2/3]

IonicNoblePurkinje ( Teuchos::ParameterList &  parameterList)
Parameters
Epetracommunicator
listof parameters in an xml file

Definition at line 65 of file IonicNoblePurkinje.cpp.

◆ IonicNoblePurkinje() [3/3]

Parameters
IonicNoblePurkinjeobject

Definition at line 76 of file IonicNoblePurkinje.cpp.

◆ ~IonicNoblePurkinje()

virtual ~IonicNoblePurkinje ( )
inlinevirtual

Destructor.

Definition at line 89 of file IonicNoblePurkinje.hpp.

Member Function Documentation

◆ operator=()

IonicNoblePurkinje & operator= ( const IonicNoblePurkinje model)

Operator.

Definition at line 93 of file IonicNoblePurkinje.cpp.

◆ gi()

const Real& gi ( ) const
inline

Definition at line 104 of file IonicNoblePurkinje.hpp.

◆ Cm()

const Real& Cm ( ) const
inline

Definition at line 110 of file IonicNoblePurkinje.hpp.

◆ vK()

const Real& vK ( ) const
inline

Definition at line 115 of file IonicNoblePurkinje.hpp.

◆ vNa()

const Real& vNa ( ) const
inline

Definition at line 120 of file IonicNoblePurkinje.hpp.

◆ Itotal()

const Real& Itotal ( ) const
inline

Definition at line 125 of file IonicNoblePurkinje.hpp.

◆ setgi()

void setgi ( const Real p)
inline

Definition at line 131 of file IonicNoblePurkinje.hpp.

◆ setCm()

void setCm ( const Real p)
inline

Definition at line 137 of file IonicNoblePurkinje.hpp.

◆ setvNa()

void setvNa ( const Real p)
inline

Definition at line 142 of file IonicNoblePurkinje.hpp.

◆ setvK()

void setvK ( const Real p)
inline

Definition at line 147 of file IonicNoblePurkinje.hpp.

◆ GeneralFunctionAlphaAndBeta()

static Real GeneralFunctionAlphaAndBeta ( const Real V,
const Real C1,
const Real C2,
const Real C3,
const Real C4,
const Real C5,
const Real V_0 
)
inlinestatic

Definition at line 153 of file IonicNoblePurkinje.hpp.

◆ mInf()

static Real mInf ( const Real V)
inlinestatic

Definition at line 166 of file IonicNoblePurkinje.hpp.

◆ hInf()

static Real hInf ( const Real V)
inlinestatic

Definition at line 174 of file IonicNoblePurkinje.hpp.

◆ nInf()

static Real nInf ( const Real V)
inlinestatic

Definition at line 182 of file IonicNoblePurkinje.hpp.

◆ computeGatingRhs()

void computeGatingRhs ( const std::vector< Real > &  v,
std::vector< Real > &  rhs 
)
virtual

Methods.

Implements ElectroIonicModel.

Definition at line 111 of file IonicNoblePurkinje.cpp.

◆ computeRhs()

void computeRhs ( const std::vector< Real > &  v,
std::vector< Real > &  rhs 
)
virtual

This methods contains the actual evaluation of the rhs of all state variablesin the model (0D version)

Although this method can just wrap the computeGatingRhs ( const std::vector<Real>& v, std::vector<Real>& rhs ) and the computeLocalPotentialRhs ( const std::vector<Real>& v ) methods, for efficiency it may be better to duplicate the code.

Parameters
vvector of state variables including the voltage (with n elements)
rhsvector of right hand side state variables including the voltage (with n elements)

Implements ElectroIonicModel.

Definition at line 132 of file IonicNoblePurkinje.cpp.

◆ computeLocalPotentialRhs()

Real computeLocalPotentialRhs ( const std::vector< Real > &  v)
virtual

This methods contains the actual evaluation of the rhs of the voltage equation only (0D version)

Parameters
vvector of state variables including the voltage (with n elements)

Implements ElectroIonicModel.

Definition at line 160 of file IonicNoblePurkinje.cpp.

◆ computeGatingVariablesWithRushLarsen()

void computeGatingVariablesWithRushLarsen ( std::vector< Real > &  ,
const Real   
)
virtual

This methods contains the actual evaluation of the rhs of the voltage equation only (0D version)

Overload this method in order to solve the ionic model with the Rush-Larsen method in the monodomain solver

Parameters
vvector of state variables
dttime step of the simulation

Reimplemented from ElectroIonicModel.

Definition at line 188 of file IonicNoblePurkinje.cpp.

◆ showMe()

void showMe ( )
virtual

Display information about the model.

Implements ElectroIonicModel.

Definition at line 214 of file IonicNoblePurkinje.cpp.

Field Documentation

◆ M_gi

Real M_gi
private

Solves the ionic model.

Model Parameters Chemical kinetics parameters

Definition at line 216 of file IonicNoblePurkinje.hpp.

◆ M_vNa

Real M_vNa
private

Definition at line 217 of file IonicNoblePurkinje.hpp.

◆ M_vK

Real M_vK
private

Definition at line 218 of file IonicNoblePurkinje.hpp.

◆ M_Cm

Real M_Cm
private

Definition at line 219 of file IonicNoblePurkinje.hpp.

◆ M_Itot

Real M_Itot
private

Definition at line 220 of file IonicNoblePurkinje.hpp.


The documentation for this class was generated from the following files: