Chaste Release::3.1
AbstractCardiacCellInterface Class Reference

#include <AbstractCardiacCellInterface.hpp>

Inheritance diagram for AbstractCardiacCellInterface:
Collaboration diagram for AbstractCardiacCellInterface:

List of all members.

Public Member Functions

 AbstractCardiacCellInterface (boost::shared_ptr< AbstractIvpOdeSolver > pOdeSolver, unsigned voltageIndex, boost::shared_ptr< AbstractStimulusFunction > pIntracellularStimulus)
virtual ~AbstractCardiacCellInterface ()
virtual void SetTimestep (double dt)=0
virtual unsigned GetNumberOfStateVariables () const =0
virtual unsigned GetNumberOfParameters () const =0
virtual std::vector< doubleGetStdVecStateVariables ()=0
virtual const std::vector
< std::string > & 
rGetStateVariableNames () const =0
virtual void SetStateVariables (const std::vector< double > &rVariables)=0
virtual void SetStateVariable (unsigned index, double newValue)=0
virtual void SetStateVariable (const std::string &rName, double newValue)=0
virtual double GetAnyVariable (const std::string &rName, double time)=0
virtual double GetParameter (const std::string &rParameterName)=0
virtual double GetParameter (unsigned parameterIndex)=0
virtual void SetParameter (const std::string &rParameterName, double value)=0
virtual void SetParameter (unsigned parameterIdx, double value)=0
virtual void SolveAndUpdateState (double tStart, double tEnd)=0
virtual OdeSolution Compute (double tStart, double tEnd, double tSamp=0.0)=0
virtual void ComputeExceptVoltage (double tStart, double tEnd)=0
virtual double GetIIonic (const std::vector< double > *pStateVariables=NULL)=0
virtual void SetVoltage (double voltage)=0
virtual double GetVoltage ()=0
unsigned GetVoltageIndex ()
void SetStimulusFunction (boost::shared_ptr< AbstractStimulusFunction > pStimulus)
double GetStimulus (double time)
void SetIntracellularStimulusFunction (boost::shared_ptr< AbstractStimulusFunction > pStimulus)
double GetIntracellularStimulus (double time)
double GetIntracellularAreaStimulus (double time)
void SetUsedInTissueSimulation (bool tissue=true)
virtual void UseCellMLDefaultStimulus ()
bool HasCellMLDefaultStimulus ()
virtual
AbstractLookupTableCollection
GetLookupTableCollection ()
boost::shared_ptr
< AbstractStimulusFunction
GetStimulusFunction ()
const boost::shared_ptr
< AbstractStimulusFunction
GetStimulusFunction () const
const boost::shared_ptr
< AbstractIvpOdeSolver
GetSolver () const
virtual void SetVoltageDerivativeToZero (bool clamp=true)
void SetFixedVoltage (double voltage)
virtual void SetStretch (double stretch)
virtual double GetIntracellularCalciumConcentration ()

Protected Attributes

unsigned mVoltageIndex
boost::shared_ptr
< AbstractIvpOdeSolver
mpOdeSolver
boost::shared_ptr
< AbstractStimulusFunction
mpIntracellularStimulus
bool mSetVoltageDerivativeToZero
bool mIsUsedInTissue
bool mHasDefaultStimulusFromCellML
double mFixedVoltage

Private Member Functions

template<class Archive >
void serialize (Archive &archive, const unsigned int version)

Friends

class boost::serialization::access

Detailed Description

This class defines a common interface to AbstractCardiacCell and AbstractCvodeCell, primarily for the benefit of single-cell cardiac simulations, so that calling code doesn't need to know which solver is being used.

Strictly speaking this isn't an interface, since some methods have implementations defined. But the name AbstractCardiacCell was already taken.

Note that serialization is not defined for this class. AbstractCvodeCell does not have (or need) serialization at present, so adding serialization here would break archive backwards compatibility for little gain.

Definition at line 62 of file AbstractCardiacCellInterface.hpp.


Constructor & Destructor Documentation

AbstractCardiacCellInterface::AbstractCardiacCellInterface ( boost::shared_ptr< AbstractIvpOdeSolver pOdeSolver,
unsigned  voltageIndex,
boost::shared_ptr< AbstractStimulusFunction pIntracellularStimulus 
)

Create a new cardiac cell. The state variables of the cell will be set to AbstractOdeSystemInformation::GetInitialConditions(). Note that calls to SetDefaultInitialConditions() on a particular instance of this class will not modify its state variables. You can modify them directly with rGetStateVariables().

Parameters:
pOdeSolverthe ODE solver to use when simulating this cell (ignored for some subclasses; can be an empty pointer in these cases)
voltageIndexthe index of the transmembrane potential within the system; see mVoltageIndex
pIntracellularStimulusthe intracellular stimulus current

Definition at line 43 of file AbstractCardiacCellInterface.cpp.

AbstractCardiacCellInterface::~AbstractCardiacCellInterface ( ) [virtual]

Virtual destructor

Definition at line 58 of file AbstractCardiacCellInterface.cpp.


Member Function Documentation

virtual OdeSolution AbstractCardiacCellInterface::Compute ( double  tStart,
double  tEnd,
double  tSamp = 0.0 
) [pure virtual]

Simulates this cell's behaviour between the time interval [tStart, tEnd], and return state variable values. The timestep used will depend on the subclass implementation.

Parameters:
tStartbeginning of the time interval to simulate
tEndend of the time interval to simulate
tSampsampling interval for returned results (defaults to dt)

Implemented in AbstractBackwardEulerCardiacCell< SIZE >, AbstractCardiacCell, AbstractCvodeCell, and AbstractRushLarsenCardiacCell.

virtual void AbstractCardiacCellInterface::ComputeExceptVoltage ( double  tStart,
double  tEnd 
) [pure virtual]

Simulates this cell's behaviour between the time interval [tStart, tEnd], but does not update the voltage. The timestep used will depend on the subclass implementation.

Parameters:
tStartbeginning of the time interval to simulate
tEndend of the time interval to simulate

Implemented in AbstractBackwardEulerCardiacCell< SIZE >, AbstractCardiacCell, AbstractCvodeCell, AbstractRushLarsenCardiacCell, and FakeBathCell.

virtual double AbstractCardiacCellInterface::GetAnyVariable ( const std::string &  rName,
double  time 
) [pure virtual]

All subclasses must implement a method that returns a state variable value

This needs to be declared here as AbstractCardiacProblem uses it.

Parameters:
rNamevariable name
timethe time at which to get the variable (only needed when evaluating derived quantities).
Returns:
value of the variable at the current time

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

virtual double AbstractCardiacCellInterface::GetIIonic ( const std::vector< double > *  pStateVariables = NULL) [pure virtual]

Computes the total current flowing through the cell membrane, using the current values of the state variables.

Should return a value in units of microamps/cm^2. Note that many cell models do not use these dimensions (let alone these units) and so a complex conversion is required. There are 2 main cases:

  • Cell model uses amps per unit capacitance. Often in this case the units used for the cell capacitance don't make sense (e.g. uF/cm^2 is used, and dV/dt=I/C_m). Hence we suggest examining the equation for dV/dt given in the model to determine what the cell model really considers the value for C_m to be, and scaling by Chaste's C_m / cell model C_m (the latter implicitly being dimensionless).
  • Cell model uses amps. In this case you need to divide by an estimate of the cell surface area. Assuming the model represents a single cell, and gives C_m in farads, then scaling by Chaste's C_m / model C_m seems reasonable. If the 'cell model' doesn't actually represent a single whole cell, then you'll have to be more careful. In both cases additional scaling may be required to obtain correct units once the dimensions have been sorted out.

Chaste's value for C_m can be obtained from HeartConfig::Instance()->GetCapacitance() and is measured in uF/cm^2.

For state-variable interpolation (SVI) we need to interpolate state variables at nodes onto quadrature points and then pass these into this method (see optional argument). Otherwise for ionic-current interpolation (ICI) we use the cell's internal state variables.

Parameters:
pStateVariablesoptionally can be supplied to evaluate the ionic current at the given state; by default the cell's internal state will be used.

Implemented in FakeBathCell, FitzHughNagumo1961OdeSystem, CML_noble_varghese_kohl_noble_1998_basic_with_sac, CorriasBuistICCModified, and CorriasBuistSMCModified.

Referenced by MonodomainCorrectionTermAssembler< ELEM_DIM, SPACE_DIM >::ComputeVectorTerm(), and BidomainCorrectionTermAssembler< ELEM_DIM, SPACE_DIM >::ComputeVectorTerm().

double AbstractCardiacCellInterface::GetIntracellularAreaStimulus ( double  time)

Get the value of the intracellular stimulus. This will always be in units of uA/cm^2.

Parameters:
timethe time at which to evaluate the stimulus

Definition at line 93 of file AbstractCardiacCellInterface.cpp.

References GetIntracellularStimulus(), HeartConfig::GetSurfaceAreaToVolumeRatio(), HeartConfig::Instance(), and mIsUsedInTissue.

Referenced by FitzHughNagumo1961OdeSystem::EvaluateYDerivatives().

double AbstractCardiacCellInterface::GetIntracellularCalciumConcentration ( ) [virtual]

[Ca_i] is needed for mechanics, so we explicitly have a Get method (rather than use a get by name type method, to avoid inefficiency when using different types of cells). This method by default throws an exception, so should be implemented in the concrete class if intracellular (cytosolic) calcium concentration is one of the state variables.

Reimplemented in FakeBathCell, and CML_noble_varghese_kohl_noble_1998_basic_with_sac.

Definition at line 154 of file AbstractCardiacCellInterface.cpp.

References EXCEPTION.

double AbstractCardiacCellInterface::GetIntracellularStimulus ( double  time)

Get the value of the intracellular stimulus. This will have units of uA/cm^2 for single-cell problems, or uA/cm^3 in a tissue simulation.

Parameters:
timethe time at which to evaluate the stimulus

Definition at line 87 of file AbstractCardiacCellInterface.cpp.

References mpIntracellularStimulus.

Referenced by GetIntracellularAreaStimulus(), and GetStimulus().

virtual AbstractLookupTableCollection* AbstractCardiacCellInterface::GetLookupTableCollection ( ) [inline, virtual]

If this cell uses lookup tables, returns the singleton object containing the tables for this cell's concrete class. Otherwise, returns NULL.

Must be implemented by subclasses iff they support lookup tables.

Definition at line 359 of file AbstractCardiacCellInterface.hpp.

Referenced by HeartConfigRelatedCellFactory< SPACE_DIM >::CreateCellWithIntracellularStimulus().

virtual unsigned AbstractCardiacCellInterface::GetNumberOfParameters ( ) const [pure virtual]

All subclasses must implement this method to get the number of parameters

Returns:
the number of parameters in the cell model ODE system

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

virtual unsigned AbstractCardiacCellInterface::GetNumberOfStateVariables ( ) const [pure virtual]

All subclasses must implement this method to get the number of state variables.

This needs to be declared here as AbstractCorrectionTermAssembler uses it.

Returns:
the number of state variables.

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

Referenced by AbstractCardiacTissue< ELEMENT_DIM, SPACE_DIM >::SolveCellSystems().

virtual double AbstractCardiacCellInterface::GetParameter ( const std::string &  rParameterName) [pure virtual]

All subclasses must implement a method that returns a parameter value.

This needs to be declared here as HeartConfigCellFactory uses it.

Parameters:
rParameterNamethe name of a parameter to get the value of,
Returns:
the parameter's value.

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

Referenced by HeartConfigRelatedCellFactory< SPACE_DIM >::SetCellParameters().

virtual double AbstractCardiacCellInterface::GetParameter ( unsigned  parameterIndex) [pure virtual]

All subclasses must implement a method that returns a parameter value.

This needs to be declared here as HeartConfigCellFactory uses it.

Parameters:
parameterIndexthe index of a parameter to get the value of,
Returns:
the parameter's value.

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

const boost::shared_ptr< AbstractIvpOdeSolver > AbstractCardiacCellInterface::GetSolver ( ) const

For boost archiving use only (that's why the 'consts' are required)

Returns:
pointer to the ODE solver being used

Definition at line 135 of file AbstractCardiacCellInterface.cpp.

References mpOdeSolver.

virtual std::vector<double> AbstractCardiacCellInterface::GetStdVecStateVariables ( ) [pure virtual]

All subclasses must implement a method that returns state variables as a std::vector (even if they work with another format), for compatibility with the PDE solvers.

This needs to be declared here as AbstractCorrectionTermAssembler uses it.

Returns:
the cell models' internal state variables

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

Referenced by AbstractCardiacTissue< ELEMENT_DIM, SPACE_DIM >::SolveCellSystems().

double AbstractCardiacCellInterface::GetStimulus ( double  time)

Get the value of the intracellular stimulus. Shorthand for GetIntracellularStimulus.

Parameters:
timethe time at which to evaluate the stimulus

Definition at line 75 of file AbstractCardiacCellInterface.cpp.

References GetIntracellularStimulus().

Referenced by CorriasBuistSMCModified::EvaluateYDerivatives(), CorriasBuistICCModified::EvaluateYDerivatives(), and CML_noble_varghese_kohl_noble_1998_basic_with_sac::EvaluateYDerivatives().

const boost::shared_ptr< AbstractStimulusFunction > AbstractCardiacCellInterface::GetStimulusFunction ( ) const

For boost archiving use only (that's why the 'consts' are required)

Returns:
The Intracellular stimulus function pointer

Definition at line 130 of file AbstractCardiacCellInterface.cpp.

References mpIntracellularStimulus.

boost::shared_ptr< AbstractStimulusFunction > AbstractCardiacCellInterface::GetStimulusFunction ( )
Returns:
The Intracellular stimulus function pointer

Definition at line 124 of file AbstractCardiacCellInterface.cpp.

References mpIntracellularStimulus.

Referenced by AbstractPurkinjeCellFactory< ELEMENT_DIM, SPACE_DIM >::CreateJunction().

virtual double AbstractCardiacCellInterface::GetVoltage ( ) [pure virtual]

Get the current value of the cellular transmembrane potential.

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

Referenced by PurkinjeVentricularJunctionStimulus::GetStimulus(), and SetVoltageDerivativeToZero().

unsigned AbstractCardiacCellInterface::GetVoltageIndex ( )

Get the index of the cellular transmembrane potential within this system. Usually it will be an index into the state variable vector, but this is not guaranteed. It will however always be suitable for use with AbstractParameterisedSystem::GetAnyVariable and OdeSolution::GetVariableAtIndex.

Definition at line 63 of file AbstractCardiacCellInterface.cpp.

References mVoltageIndex.

Referenced by AdaptiveBidomainProblem::AddCurrentSolutionToAdaptiveMesh(), AdaptiveBidomainProblem::InitializeSolutionOnAdaptedMesh(), and AbstractRushLarsenCardiacCell::UpdateTransmembranePotential().

bool AbstractCardiacCellInterface::HasCellMLDefaultStimulus ( )
Returns:
Whether the cell was generated from a CellML file with stimulus metadata.

Definition at line 119 of file AbstractCardiacCellInterface.cpp.

References mHasDefaultStimulusFromCellML.

Referenced by CellMLLoader::LoadCellMLFile().

virtual const std::vector<std::string>& AbstractCardiacCellInterface::rGetStateVariableNames ( ) const [pure virtual]

All subclasses must implement this method to get variable names.

Needs to be declared here as AdaptiveBidomainProblem uses it.

Returns:
the state variable names in the cell's ODE system.

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

Referenced by AdaptiveBidomainProblem::AddCurrentSolutionToAdaptiveMesh(), and AdaptiveBidomainProblem::InitializeSolutionOnAdaptedMesh().

template<class Archive >
void AbstractCardiacCellInterface::serialize ( Archive &  archive,
const unsigned int  version 
) [inline, private]

Main boost serialization method, some member variables are handled by Save/Load constructs.

Parameters:
archivethe archive file.
versionthe version of archiving, defined in the macro at the bottom of the class.

Reimplemented in AbstractBackwardEulerCardiacCell< SIZE >, AbstractCardiacCell, AbstractCvodeCell, AbstractRushLarsenCardiacCell, FakeBathCell, CML_noble_varghese_kohl_noble_1998_basic_with_sac, CorriasBuistICCModified, and CorriasBuistSMCModified.

Definition at line 461 of file AbstractCardiacCellInterface.hpp.

References mHasDefaultStimulusFromCellML, mIsUsedInTissue, and mSetVoltageDerivativeToZero.

void AbstractCardiacCellInterface::SetFixedVoltage ( double  voltage)

When the voltage derivative has been set to zero by SetVoltageDerivativeToZero, this method sets the transmembrane potential to use when computing the other derivatives and total ionic current.

Parameters:
voltagethe value of the transmembrane potential

Definition at line 149 of file AbstractCardiacCellInterface.cpp.

References mFixedVoltage.

Referenced by AbstractCvodeCell::SetVoltage(), and AbstractCardiacCell::SetVoltage().

void AbstractCardiacCellInterface::SetIntracellularStimulusFunction ( boost::shared_ptr< AbstractStimulusFunction pStimulus)

Set the intracellular stimulus. This should have units of uA/cm^2 for single-cell problems, or uA/cm^3 in a tissue simulation.

Parameters:
pStimulusnew stimulus function

Definition at line 81 of file AbstractCardiacCellInterface.cpp.

References mpIntracellularStimulus.

Referenced by HeartConfigRelatedCellFactory< SPACE_DIM >::SetCellIntracellularStimulus(), and SetStimulusFunction().

virtual void AbstractCardiacCellInterface::SetParameter ( const std::string &  rParameterName,
double  value 
) [pure virtual]

All subclasses must implement a method that sets a parameter value.

This needs to be declared here as HeartConfigCellFactory uses it.

Parameters:
rParameterNamethe parameter name to set the value of,
valuevalue to set it to.

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

Referenced by HeartConfigRelatedCellFactory< SPACE_DIM >::SetCellParameters().

virtual void AbstractCardiacCellInterface::SetParameter ( unsigned  parameterIdx,
double  value 
) [pure virtual]

All subclasses must implement a method that sets a parameter value.

This needs to be declared here to avoid confusion with the above method.ns

Parameters:
parameterIdxthe index of the parameter to set the value of,
valuevalue to set it to.

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

virtual void AbstractCardiacCellInterface::SetStateVariable ( const std::string &  rName,
double  newValue 
) [pure virtual]

Set the value of a single state variable in the ODE system.

This needs to be declared here as AdaptiveBidomainProblem uses the above, to avoid compiler confusion...!

Parameters:
rNamename of the state variable to be set
newValuenew value of the state variable

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

virtual void AbstractCardiacCellInterface::SetStateVariable ( unsigned  index,
double  newValue 
) [pure virtual]

Set the value of a single state variable in the ODE system.

This needs to be declared here as AdaptiveBidomainProblem uses it.

Parameters:
indexindex of the state variable to be set
newValuenew value of the state variable

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

Referenced by AdaptiveBidomainProblem::InitializeSolutionOnAdaptedMesh().

virtual void AbstractCardiacCellInterface::SetStateVariables ( const std::vector< double > &  rVariables) [pure virtual]

All subclasses must implement this method to set the state variables.

This needs to be declared here as AbstractCardiacTissue uses it.

Parameters:
rVariablesthe state variables to take a copy of.

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

Referenced by AbstractCardiacTissue< ELEMENT_DIM, SPACE_DIM >::SolveCellSystems().

void AbstractCardiacCellInterface::SetStimulusFunction ( boost::shared_ptr< AbstractStimulusFunction pStimulus)

Set the intracellular stimulus. Shorthand for SetIntracellularStimulusFunction.

Parameters:
pStimulusnew stimulus function

Definition at line 69 of file AbstractCardiacCellInterface.cpp.

References SetIntracellularStimulusFunction().

Referenced by AbstractPurkinjeCellFactory< ELEMENT_DIM, SPACE_DIM >::CreateJunction().

virtual void AbstractCardiacCellInterface::SetStretch ( double  stretch) [inline, virtual]

In electromechanics problems, the stretch is passed back to cell-model in case mechano-electric feedback has been modelled. We define an empty method here. Stretch-dependent cell models should overload this method to use the input stretch accordingly.

Parameters:
stretchthe stretch of the cell in the axial direction

Reimplemented in CML_noble_varghese_kohl_noble_1998_basic_with_sac.

Definition at line 407 of file AbstractCardiacCellInterface.hpp.

virtual void AbstractCardiacCellInterface::SetTimestep ( double  dt) [pure virtual]

Set the timestep (or maximum timestep) to use for simulating this cell.

Parameters:
dtthe timestep

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

void AbstractCardiacCellInterface::SetUsedInTissueSimulation ( bool  tissue = true)

Set whether this cell object exists in the context of a tissue simulation, or can be used for single cell simulations. This affects the units of the intracellular stimulus (see GetIntracellularStimulus) and so is used by GetIntracellularAreaStimulus to perform a units conversion if necessary.

Parameters:
tissuetrue if cell is in a tissue

Definition at line 108 of file AbstractCardiacCellInterface.cpp.

References mIsUsedInTissue.

Referenced by AbstractCardiacTissue< ELEMENT_DIM, SPACE_DIM >::SetUpHaloCells().

virtual void AbstractCardiacCellInterface::SetVoltage ( double  voltage) [pure virtual]

Set the cellular transmembrane potential.

Parameters:
voltagenew value

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

void AbstractCardiacCellInterface::SetVoltageDerivativeToZero ( bool  clamp = true) [virtual]

Set whether to clamp the voltage by setting its derivative to zero.

Parameters:
clampwhether to clamp

Reimplemented in AbstractCvodeCell.

Definition at line 140 of file AbstractCardiacCellInterface.cpp.

References DOUBLE_UNSET, GetVoltage(), mFixedVoltage, and mSetVoltageDerivativeToZero.

Referenced by AbstractCardiacCell::ComputeExceptVoltage().

virtual void AbstractCardiacCellInterface::SolveAndUpdateState ( double  tStart,
double  tEnd 
) [pure virtual]

Simulate this cell's behaviour between the time interval [tStart, tEnd], updating the internal state variable values. The timestep used will depend on the subclass implementation.

Parameters:
tStartbeginning of the time interval to simulate
tEndend of the time interval to simulate

Implemented in AbstractBackwardEulerCardiacCell< SIZE >, AbstractCardiacCell, AbstractCvodeCell, and AbstractRushLarsenCardiacCell.

void AbstractCardiacCellInterface::UseCellMLDefaultStimulus ( ) [virtual]

Use CellML metadata to set up the default stimulus for this cell. By default this method will always throw an exception. For suitably annotated models, PyCml will override this to provide a RegularStimulus as defined in the CellML.

Definition at line 113 of file AbstractCardiacCellInterface.cpp.

References EXCEPTION, and mHasDefaultStimulusFromCellML.

Referenced by CellMLLoader::LoadCellMLFile().


Friends And Related Function Documentation


Member Data Documentation

The value of the fixed voltage if mSetVoltageDerivativeToZero is set.

Definition at line 449 of file AbstractCardiacCellInterface.hpp.

Referenced by SetFixedVoltage(), and SetVoltageDerivativeToZero().

Whether this cell has a default stimulus specified by CellML metadata.

Definition at line 446 of file AbstractCardiacCellInterface.hpp.

Referenced by HasCellMLDefaultStimulus(), serialize(), AbstractCardiacCell::serialize(), and UseCellMLDefaultStimulus().

Whether this cell exists in a tissue, or is an isolated cell.

Definition at line 443 of file AbstractCardiacCellInterface.hpp.

Referenced by GetIntracellularAreaStimulus(), serialize(), AbstractCardiacCell::serialize(), and SetUsedInTissueSimulation().

The index of the voltage within this system. Usually it will be an index into the state variable vector, but this is not guaranteed. It will however always be suitable for use with AbstractParameterisedSystem::GetAnyVariable and OdeSolution::GetVariableAtIndex.

Definition at line 427 of file AbstractCardiacCellInterface.hpp.

Referenced by FitzHughNagumo1961OdeSystem::GetIIonic(), AbstractCvodeCell::GetVoltage(), AbstractCardiacCell::GetVoltage(), GetVoltageIndex(), AbstractCvodeCell::SetVoltage(), and AbstractCardiacCell::SetVoltage().


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