LifeV
IonicFitzHughNagumo Class Reference

#include <IonicFitzHughNagumo.hpp>

+ Inheritance diagram for IonicFitzHughNagumo:
+ Collaboration diagram for IonicFitzHughNagumo:

Constructors & Destructor

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

Overloads

IonicFitzHughNagumooperator= (const IonicFitzHughNagumo &model)
 Operator. More...
 

Setters and getters

Real M_G
 Model Parameters. More...
 
Real M_Vth
 
Real M_Vp
 
Real M_Eta1
 
Real M_Eta2
 
Real M_Eta3
 
Real M_Eta
 
Real M_Gamma
 
const RealG () const
 
const RealVth () const
 
const RealVp () const
 
const RealEta1 () const
 
const RealEta2 () const
 
const RealEta3 () const
 
const RealEta () const
 
const RealGamma () const
 
void setG (const Real &G)
 
void setVth (const Real &Vth)
 
void setVp (const Real &Vp)
 
void setEta1 (const Real &Eta1)
 
void setEta2 (const Real &Eta2)
 
void setEta3 (const Real &Eta3)
 
void setup (Teuchos::ParameterList &parameterList)
 Methods. More...
 
void computeGatingRhs (const std::vector< Real > &v, std::vector< Real > &rhs)
 A new ionic model must have the following 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...
 
void computeRhs (const VectorSmall< 2 > &v, VectorSmall< 2 > &rhs)
 
Real computeLocalPotentialRhs (const std::vector< Real > &v)
 This methods contains the actual evaluation of the rhs of the voltage equation only (0D version) More...
 
matrix_Type getJac (const vector_Type &v, Real h=1.0e-8)
 This methods computes the Jacobian numerically. More...
 
std::vector< std::vector< Real > > getJac (const std::vector< Real > &v, Real h=1.0e-8)
 This methods computes the Jacobian numerically. More...
 
MatrixSmall< 2, 2 > getJac (const VectorSmall< 2 > &v, Real h=1.0e-8)
 
void showMe ()
 Display information about the model. More...
 

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
 
- Public Member Functions inherited from ElectroIonicModel
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 computeGatingVariablesWithRushLarsen (std::vector< Real > &, const Real)
 This methods contains the actual evaluation of the rhs of the voltage equation only (0D version) 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...
 
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...
 
- 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

Definition at line 73 of file IonicFitzHughNagumo.hpp.

Constructor & Destructor Documentation

◆ IonicFitzHughNagumo() [1/3]

Constructor.

Constructors.

Definition at line 45 of file IonicFitzHughNagumo.cpp.

◆ IonicFitzHughNagumo() [2/3]

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

Definition at line 61 of file IonicFitzHughNagumo.cpp.

◆ IonicFitzHughNagumo() [3/3]

Parameters
IonicFitzHughNagumoobject

Definition at line 70 of file IonicFitzHughNagumo.cpp.

◆ ~IonicFitzHughNagumo()

virtual ~IonicFitzHughNagumo ( )
inlinevirtual

Destructor.

Definition at line 94 of file IonicFitzHughNagumo.hpp.

Member Function Documentation

◆ operator=()

IonicFitzHughNagumo & operator= ( const IonicFitzHughNagumo model)

Operator.

Definition at line 88 of file IonicFitzHughNagumo.cpp.

◆ G()

const Real& G ( ) const
inline

Definition at line 109 of file IonicFitzHughNagumo.hpp.

◆ Vth()

const Real& Vth ( ) const
inline

Definition at line 113 of file IonicFitzHughNagumo.hpp.

◆ Vp()

const Real& Vp ( ) const
inline

Definition at line 117 of file IonicFitzHughNagumo.hpp.

◆ Eta1()

const Real& Eta1 ( ) const
inline

Definition at line 121 of file IonicFitzHughNagumo.hpp.

◆ Eta2()

const Real& Eta2 ( ) const
inline

Definition at line 125 of file IonicFitzHughNagumo.hpp.

◆ Eta3()

const Real& Eta3 ( ) const
inline

Definition at line 129 of file IonicFitzHughNagumo.hpp.

◆ Eta()

const Real& Eta ( ) const
inline

Definition at line 133 of file IonicFitzHughNagumo.hpp.

◆ Gamma()

const Real& Gamma ( ) const
inline

Definition at line 137 of file IonicFitzHughNagumo.hpp.

◆ setG()

void setG ( const Real G)
inline

Definition at line 142 of file IonicFitzHughNagumo.hpp.

◆ setVth()

void setVth ( const Real Vth)
inline

Definition at line 146 of file IonicFitzHughNagumo.hpp.

◆ setVp()

void setVp ( const Real Vp)
inline

Definition at line 150 of file IonicFitzHughNagumo.hpp.

◆ setEta1()

void setEta1 ( const Real Eta1)
inline

Definition at line 155 of file IonicFitzHughNagumo.hpp.

◆ setEta2()

void setEta2 ( const Real Eta2)
inline

Definition at line 159 of file IonicFitzHughNagumo.hpp.

◆ setEta3()

void setEta3 ( const Real Eta3)
inline

Definition at line 165 of file IonicFitzHughNagumo.hpp.

◆ setup()

void setup ( Teuchos::ParameterList &  parameterList)
virtual

Methods.

Reimplemented from ElectroIonicModel.

Definition at line 109 of file IonicFitzHughNagumo.cpp.

◆ computeGatingRhs()

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

A new ionic model must have the following methods ///.

This methods contains the actual evaluation of the rhs of all state variables except for the voltage equation (0D version)

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

Implements ElectroIonicModel.

Definition at line 123 of file IonicFitzHughNagumo.cpp.

◆ computeRhs() [1/2]

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 IonicFitzHughNagumo.cpp.

◆ computeRhs() [2/2]

void computeRhs ( const VectorSmall< 2 > &  v,
VectorSmall< 2 > &  rhs 
)

Definition at line 142 of file IonicFitzHughNagumo.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 152 of file IonicFitzHughNagumo.cpp.

◆ getJac() [1/3]

matrix_Type getJac ( const vector_Type v,
Real  h = 1.0e-8 
)
virtual

This methods computes the Jacobian numerically.

Parameters
vvector of pointers to the state variables vectors
hdifferentiation step

Reimplemented from ElectroIonicModel.

Definition at line 179 of file IonicFitzHughNagumo.cpp.

◆ getJac() [2/3]

std::vector< std::vector< Real > > getJac ( const std::vector< Real > &  v,
Real  h = 1.0e-8 
)
virtual

This methods computes the Jacobian numerically.

Parameters
vvector of the state variables
hdifferentiation step

Reimplemented from ElectroIonicModel.

Definition at line 157 of file IonicFitzHughNagumo.cpp.

◆ getJac() [3/3]

MatrixSmall< 2, 2 > getJac ( const VectorSmall< 2 > &  v,
Real  h = 1.0e-8 
)

Definition at line 168 of file IonicFitzHughNagumo.cpp.

◆ showMe()

void showMe ( )
virtual

Display information about the model.

Implements ElectroIonicModel.

Definition at line 211 of file IonicFitzHughNagumo.cpp.

+ Here is the caller graph for this function:

Field Documentation

◆ M_G

Real M_G
private

Model Parameters.

Chemical kinetics parameters

Definition at line 195 of file IonicFitzHughNagumo.hpp.

◆ M_Vth

Real M_Vth
private

Definition at line 196 of file IonicFitzHughNagumo.hpp.

◆ M_Vp

Real M_Vp
private

Definition at line 197 of file IonicFitzHughNagumo.hpp.

◆ M_Eta1

Real M_Eta1
private

Definition at line 198 of file IonicFitzHughNagumo.hpp.

◆ M_Eta2

Real M_Eta2
private

Definition at line 199 of file IonicFitzHughNagumo.hpp.

◆ M_Eta3

Real M_Eta3
private

Definition at line 200 of file IonicFitzHughNagumo.hpp.

◆ M_Eta

Real M_Eta
private

Definition at line 201 of file IonicFitzHughNagumo.hpp.

◆ M_Gamma

Real M_Gamma
private

Definition at line 202 of file IonicFitzHughNagumo.hpp.


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