Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
AbstractCardiacCellInterface Class Referenceabstract

#include <AbstractCardiacCellInterface.hpp>

+ Inheritance diagram for AbstractCardiacCellInterface:
+ Collaboration diagram for AbstractCardiacCellInterface:

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 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 boost::shared_ptr< RegularStimulusUseCellMLDefaultStimulus ()
 
bool HasCellMLDefaultStimulus ()
 
virtual AbstractLookupTableCollectionGetLookupTableCollection ()
 
boost::shared_ptr< AbstractStimulusFunctionGetStimulusFunction ()
 
const boost::shared_ptr< AbstractStimulusFunctionGetStimulusFunction () const
 
const boost::shared_ptr< AbstractIvpOdeSolverGetSolver () const
 
void SetSolver (boost::shared_ptr< AbstractIvpOdeSolver > pSolver)
 
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< AbstractIvpOdeSolvermpOdeSolver
 
boost::shared_ptr< AbstractStimulusFunctionmpIntracellularStimulus
 
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.

Definition at line 58 of file AbstractCardiacCellInterface.hpp.

Constructor & Destructor Documentation

◆ AbstractCardiacCellInterface()

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 53 of file AbstractCardiacCellInterface.cpp.

References Citations::Register().

◆ ~AbstractCardiacCellInterface()

AbstractCardiacCellInterface::~AbstractCardiacCellInterface ( )
virtual

Virtual destructor

Definition at line 71 of file AbstractCardiacCellInterface.cpp.

Member Function Documentation

◆ Compute()

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.

Returns
state variable values.
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 >, AbstractBackwardEulerCardiacCell< 0u >, AbstractCardiacCell, AbstractCvodeCell, AbstractGeneralizedRushLarsenCardiacCell, and AbstractRushLarsenCardiacCell.

◆ ComputeExceptVoltage()

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 >, AbstractBackwardEulerCardiacCell< 0u >, AbstractCardiacCell, AbstractCvodeCell, AbstractGeneralizedRushLarsenCardiacCell, AbstractRushLarsenCardiacCell, and FakeBathCell.

◆ GetAnyVariable()

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.

◆ GetIIonic()

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.

Returns
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().

◆ GetIntracellularAreaStimulus()

double AbstractCardiacCellInterface::GetIntracellularAreaStimulus ( double  time)
Returns
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 106 of file AbstractCardiacCellInterface.cpp.

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

Referenced by FitzHughNagumo1961OdeSystem::EvaluateYDerivatives().

◆ GetIntracellularCalciumConcentration()

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.

Returns the intracellular calcium concentraion *in milliMolar*.

Returns
intracellular calcium concentration

Reimplemented in FakeBathCell, and CML_noble_varghese_kohl_noble_1998_basic_with_sac.

Definition at line 173 of file AbstractCardiacCellInterface.cpp.

References EXCEPTION.

◆ GetIntracellularStimulus()

double AbstractCardiacCellInterface::GetIntracellularStimulus ( double  time)
Returns
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 100 of file AbstractCardiacCellInterface.cpp.

References mpIntracellularStimulus.

Referenced by GetIntracellularAreaStimulus(), and GetStimulus().

◆ GetLookupTableCollection()

virtual AbstractLookupTableCollection * AbstractCardiacCellInterface::GetLookupTableCollection ( )
inlinevirtual
Returns
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 346 of file AbstractCardiacCellInterface.hpp.

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

◆ GetNumberOfParameters()

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.

◆ GetNumberOfStateVariables()

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.

◆ GetParameter() [1/2]

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().

◆ GetParameter() [2/2]

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.

◆ GetSolver()

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 149 of file AbstractCardiacCellInterface.cpp.

References mpOdeSolver.

◆ GetStdVecStateVariables()

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.

◆ GetStimulus()

double AbstractCardiacCellInterface::GetStimulus ( double  time)
Returns
the value of the intracellular stimulus. Shorthand for GetIntracellularStimulus.
Parameters
timethe time at which to evaluate the stimulus

Definition at line 88 of file AbstractCardiacCellInterface.cpp.

References GetIntracellularStimulus().

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

◆ GetStimulusFunction() [1/2]

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

Definition at line 138 of file AbstractCardiacCellInterface.cpp.

References mpIntracellularStimulus.

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

◆ GetStimulusFunction() [2/2]

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 144 of file AbstractCardiacCellInterface.cpp.

References mpIntracellularStimulus.

◆ GetVoltage()

virtual double AbstractCardiacCellInterface::GetVoltage ( )
pure virtual
Returns
the current value of the cellular transmembrane potential.

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

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

◆ GetVoltageIndex()

unsigned AbstractCardiacCellInterface::GetVoltageIndex ( )
Returns
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 76 of file AbstractCardiacCellInterface.cpp.

References mVoltageIndex.

Referenced by AbstractGeneralizedRushLarsenCardiacCell::Compute(), AbstractGeneralizedRushLarsenCardiacCell::SolveAndUpdateState(), and AbstractRushLarsenCardiacCell::UpdateTransmembranePotential().

◆ HasCellMLDefaultStimulus()

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

Definition at line 133 of file AbstractCardiacCellInterface.cpp.

References mHasDefaultStimulusFromCellML.

Referenced by CellMLLoader::LoadCellMLFile().

◆ rGetStateVariableNames()

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.

◆ serialize()

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

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.

Definition at line 461 of file AbstractCardiacCellInterface.hpp.

References mHasDefaultStimulusFromCellML, mIsUsedInTissue, and mSetVoltageDerivativeToZero.

◆ SetFixedVoltage()

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 168 of file AbstractCardiacCellInterface.cpp.

References mFixedVoltage.

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

◆ SetIntracellularStimulusFunction()

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 94 of file AbstractCardiacCellInterface.cpp.

References mpIntracellularStimulus.

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

◆ SetParameter()

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().

◆ SetSolver()

void AbstractCardiacCellInterface::SetSolver ( boost::shared_ptr< AbstractIvpOdeSolver pSolver)

Change the ODE solver used by this cell. Note that this will only work for cell models which do not have a solver built-in (e.g. it will NOT work for AbstractCvodeCell derivatives).

Parameters
pSolverthe new solver to use

Definition at line 154 of file AbstractCardiacCellInterface.cpp.

References mpOdeSolver.

◆ SetStateVariable() [1/2]

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.

◆ SetStateVariable() [2/2]

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.

◆ SetStateVariables()

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.

◆ SetStimulusFunction()

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

Set the intracellular stimulus. Shorthand for SetIntracellularStimulusFunction.

Parameters
pStimulusnew stimulus function

Definition at line 82 of file AbstractCardiacCellInterface.cpp.

References SetIntracellularStimulusFunction().

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

◆ SetStretch()

virtual void AbstractCardiacCellInterface::SetStretch ( double  stretch)
inlinevirtual

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 403 of file AbstractCardiacCellInterface.hpp.

◆ SetTimestep()

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

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

Parameters
dtthe timestep

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

◆ SetUsedInTissueSimulation()

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 121 of file AbstractCardiacCellInterface.cpp.

References mIsUsedInTissue.

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

◆ SetVoltage()

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

Set the cellular transmembrane potential.

Parameters
voltagenew value

Implemented in AbstractCardiacCell, and AbstractCvodeCell.

◆ SetVoltageDerivativeToZero()

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

◆ SolveAndUpdateState()

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 >, AbstractBackwardEulerCardiacCell< 0u >, AbstractCardiacCell, AbstractCvodeCell, AbstractGeneralizedRushLarsenCardiacCell, and AbstractRushLarsenCardiacCell.

◆ UseCellMLDefaultStimulus()

boost::shared_ptr< RegularStimulus > 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, chaste_codegen will override this to provide a RegularStimulus as defined in the CellML.

Returns
a regular stimulus as defined in the CellML

Definition at line 126 of file AbstractCardiacCellInterface.cpp.

References EXCEPTION, and mHasDefaultStimulusFromCellML.

Referenced by CellMLLoader::LoadCellMLFile().

Friends And Related Symbol Documentation

◆ boost::serialization::access

friend class boost::serialization::access
friend

Needed for serialization.

Definition at line 453 of file AbstractCardiacCellInterface.hpp.

Member Data Documentation

◆ mFixedVoltage

double AbstractCardiacCellInterface::mFixedVoltage
protected

The value of the fixed voltage if mSetVoltageDerivativeToZero is set.

Definition at line 449 of file AbstractCardiacCellInterface.hpp.

Referenced by SetFixedVoltage(), and SetVoltageDerivativeToZero().

◆ mHasDefaultStimulusFromCellML

bool AbstractCardiacCellInterface::mHasDefaultStimulusFromCellML
protected

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

Definition at line 446 of file AbstractCardiacCellInterface.hpp.

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

◆ mIsUsedInTissue

bool AbstractCardiacCellInterface::mIsUsedInTissue
protected

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

Definition at line 443 of file AbstractCardiacCellInterface.hpp.

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

◆ mpIntracellularStimulus

boost::shared_ptr<AbstractStimulusFunction> AbstractCardiacCellInterface::mpIntracellularStimulus
protected

The intracellular stimulus current.

Definition at line 433 of file AbstractCardiacCellInterface.hpp.

Referenced by GetIntracellularStimulus(), GetStimulusFunction(), GetStimulusFunction(), and SetIntracellularStimulusFunction().

◆ mpOdeSolver

boost::shared_ptr<AbstractIvpOdeSolver> AbstractCardiacCellInterface::mpOdeSolver
protected

◆ mSetVoltageDerivativeToZero

bool AbstractCardiacCellInterface::mSetVoltageDerivativeToZero
protected

◆ mVoltageIndex

unsigned AbstractCardiacCellInterface::mVoltageIndex
protected

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(), AbstractCardiacCell::GetVoltage(), AbstractCvodeCell::GetVoltage(), GetVoltageIndex(), AbstractCardiacCell::SetVoltage(), and AbstractCvodeCell::SetVoltage().


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