Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
AbstractCvodeCell Class Reference

#include <AbstractCvodeCell.hpp>

+ Inheritance diagram for AbstractCvodeCell:
+ Collaboration diagram for AbstractCvodeCell:

Public Member Functions

 AbstractCvodeCell (boost::shared_ptr< AbstractIvpOdeSolver > pSolver, unsigned numberOfStateVariables, unsigned voltageIndex, boost::shared_ptr< AbstractStimulusFunction > pIntracellularStimulus)
 
virtual ~AbstractCvodeCell ()
 
double GetVoltage ()
 
void SetVoltage (double voltage)
 
void SetMaxTimestep (double maxDt)
 
void SetTimestep (double maxDt)
 
double GetTimestep ()
 
virtual void SolveAndUpdateState (double tStart, double tEnd)
 
OdeSolution Compute (double tStart, double tEnd, double tSamp=0.0)
 
void ComputeExceptVoltage (double tStart, double tEnd)
 
void SetVoltageDerivativeToZero (bool clamp=true)
 
unsigned GetNumberOfStateVariables () const
 
unsigned GetNumberOfParameters () const
 
std::vector< doubleGetStdVecStateVariables ()
 
const std::vector< std::string > & rGetStateVariableNames () const
 
void SetStateVariables (const std::vector< double > &rVariables)
 
void SetStateVariables (const N_Vector &rVariables)
 
void SetStateVariable (unsigned index, double newValue)
 
void SetStateVariable (const std::string &rName, double newValue)
 
double GetAnyVariable (const std::string &rName, double time=0.0)
 
double GetParameter (const std::string &rParameterName)
 
double GetParameter (unsigned parameterIndex)
 
void SetParameter (const std::string &rParameterName, double value)
 
- Public Member Functions inherited from AbstractCardiacCellInterface
 AbstractCardiacCellInterface (boost::shared_ptr< AbstractIvpOdeSolver > pOdeSolver, unsigned voltageIndex, boost::shared_ptr< AbstractStimulusFunction > pIntracellularStimulus)
 
virtual ~AbstractCardiacCellInterface ()
 
virtual double GetIIonic (const std::vector< double > *pStateVariables=NULL)=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)
 
void SetFixedVoltage (double voltage)
 
virtual void SetStretch (double stretch)
 
virtual double GetIntracellularCalciumConcentration ()
 
- Public Member Functions inherited from AbstractCvodeSystem
 AbstractCvodeSystem (unsigned numberOfStateVariables)
 
virtual ~AbstractCvodeSystem ()
 
virtual void EvaluateYDerivatives (realtype time, const N_Vector y, N_Vector ydot)=0
 
virtual void EvaluateAnalyticJacobian (realtype time, N_Vector y, N_Vector ydot, CHASTE_CVODE_DENSE_MATRIX jacobian, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
 
void SetForceReset (bool autoReset)
 
void SetMinimalReset (bool minimalReset)
 
bool GetMinimalReset ()
 Get whether we want to run with minimal reset or not (no reinitialisation of the solver if variables change)
 
bool GetForceReset ()
 Get whether we will force a solver reset on every call to Solve()
 
void ResetSolver ()
 
OdeSolution Solve (realtype tStart, realtype tEnd, realtype maxDt, realtype tSamp)
 
void Solve (realtype tStart, realtype tEnd, realtype maxDt)
 
void SetMaxSteps (long int numSteps)
 
long int GetMaxSteps ()
 
void SetTolerances (double relTol=1e-5, double absTol=1e-7)
 
double GetRelativeTolerance ()
 
double GetAbsoluteTolerance ()
 
double GetLastStepSize ()
 
bool GetUseAnalyticJacobian () const
 
bool HasAnalyticJacobian () const
 
void ForceUseOfNumericalJacobian (bool useNumericalJacobian=true)
 
- Public Member Functions inherited from AbstractParameterisedSystem< N_Vector >
 AbstractParameterisedSystem (unsigned numberOfStateVariables)
 
N_VectorrGetStateVariables ()
 
N_Vector GetStateVariables ()
 
void SetStateVariables (const N_Vector &rStateVariables)
 
double GetStateVariable (unsigned index) const
 
double GetStateVariable (const std::string &rName) const
 
void SetStateVariable (unsigned index, double newValue)
 
void SetStateVariable (const std::string &rName, double newValue)
 
virtual void VerifyStateVariables ()
 
void SetDefaultInitialConditions (const N_Vector &rInitialConditions)
 
void SetDefaultInitialCondition (unsigned index, double initialCondition)
 
N_Vector GetInitialConditions () const
 
void ResetToInitialConditions ()
 
double GetParameter (unsigned index) const
 
double GetParameter (const std::string &rName) const
 
void SetParameter (const std::string &rName, double value)
 
void SetParameter (unsigned index, double value)
 
double GetAnyVariable (unsigned index, double time=0.0, N_Vector *pDerivedQuantities=NULL)
 
double GetAnyVariable (const std::string &rName, double time=0.0, N_Vector *pDerivedQuantities=NULL)
 
void SetAnyVariable (unsigned index, double value)
 
void SetAnyVariable (const std::string &rName, double value)
 
virtual N_Vector ComputeDerivedQuantities (double time, const N_Vector &rState)
 
N_Vector ComputeDerivedQuantitiesFromCurrentState (double time)
 
- Public Member Functions inherited from AbstractUntemplatedParameterisedSystem
 AbstractUntemplatedParameterisedSystem (unsigned numberOfStateVariables)
 
virtual ~AbstractUntemplatedParameterisedSystem ()
 
boost::shared_ptr< const AbstractOdeSystemInformationGetSystemInformation () const
 
std::string GetSystemName () const
 
unsigned GetNumberOfAttributes () const
 
bool HasAttribute (const std::string &rName) const
 
double GetAttribute (const std::string &rName) const
 
unsigned GetNumberOfStateVariables () const
 
const std::vector< std::string > & rGetStateVariableNames () const
 
const std::vector< std::string > & rGetStateVariableUnits () const
 
unsigned GetStateVariableIndex (const std::string &rName) const
 
bool HasStateVariable (const std::string &rName) const
 
std::string GetStateVariableUnits (unsigned index) const
 
unsigned GetNumberOfParameters () const
 
const std::vector< std::string > & rGetParameterNames () const
 
const std::vector< std::string > & rGetParameterUnits () const
 
unsigned GetParameterIndex (const std::string &rName) const
 
bool HasParameter (const std::string &rName) const
 
std::string GetParameterUnits (unsigned index) const
 
unsigned GetNumberOfDerivedQuantities () const
 
const std::vector< std::string > & rGetDerivedQuantityNames () const
 
const std::vector< std::string > & rGetDerivedQuantityUnits () const
 
unsigned GetDerivedQuantityIndex (const std::string &rName) const
 
bool HasDerivedQuantity (const std::string &rName) const
 
std::string GetDerivedQuantityUnits (unsigned index) const
 
unsigned GetAnyVariableIndex (const std::string &rName) const
 
bool HasAnyVariable (const std::string &rName) const
 
std::string GetAnyVariableUnits (unsigned index) const
 
std::string GetAnyVariableUnits (const std::string &rName) const
 

Protected Attributes

double mMaxDt
 
- Protected Attributes inherited from AbstractCardiacCellInterface
unsigned mVoltageIndex
 
boost::shared_ptr< AbstractIvpOdeSolvermpOdeSolver
 
boost::shared_ptr< AbstractStimulusFunctionmpIntracellularStimulus
 
bool mSetVoltageDerivativeToZero
 
bool mIsUsedInTissue
 
bool mHasDefaultStimulusFromCellML
 
double mFixedVoltage
 
- Protected Attributes inherited from AbstractCvodeSystem
bool mHasAnalyticJacobian
 
bool mUseAnalyticJacobian
 
double mRelTol
 
double mAbsTol
 
void * mpCvodeMem
 
long int mMaxSteps
 
double mLastInternalStepSize
 
- Protected Attributes inherited from AbstractParameterisedSystem< N_Vector >
N_Vector mStateVariables
 
N_Vector mParameters
 
- Protected Attributes inherited from AbstractUntemplatedParameterisedSystem
unsigned mNumberOfStateVariables
 
boost::shared_ptr< AbstractOdeSystemInformationmpSystemInfo
 

Private Member Functions

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

Friends

class boost::serialization::access
 

Additional Inherited Members

- Protected Member Functions inherited from AbstractCvodeSystem
void Init ()
 
- Protected Member Functions inherited from AbstractParameterisedSystem< N_Vector >
std::string DumpState (const std::string &rMessage)
 
std::string DumpState (const std::string &rMessage, N_Vector Y)
 
std::string DumpState (const std::string &rMessage, N_Vector Y, double time)
 
void CheckParametersOnLoad (const std::vector< double > &rParameters, const std::vector< std::string > &rParameterNames)
 

Detailed Description

A cardiac cell that is designed to be simulated using CVODE. It uses CVODE's vector type natively via AbstractCvodeSystem.

Functionality is similar to that provided by AbstractCardiacCell and AbstractOdeSystem, but not identical. It also includes a direct interface to the CVODE solver, via the Solve methods, since the CvodeAdaptor class may be a bit slower.

Assumes that it will be solving stiff systems, so uses BDF/Newton.

Various methods in this class just call methods on AbstractCvodeSystem, to reduce compiler confusion when working with things in an inheritance tree.

Any single cell work should generally use this class rather than AbstractCardiacCell.

This class may also be faster for certain tissue problems (especially when using a relatively large PDE timestep (~0.1ms)). BUT it will come with a memory overhead as every node has to carry a CVODE solver object as it stores the internal state of the solver, not just the state variables as ForwardEuler solver does. So may not be appropriate for very large meshes.

Definition at line 78 of file AbstractCvodeCell.hpp.

Constructor & Destructor Documentation

◆ AbstractCvodeCell()

AbstractCvodeCell::AbstractCvodeCell ( boost::shared_ptr< AbstractIvpOdeSolver pSolver,
unsigned  numberOfStateVariables,
unsigned  voltageIndex,
boost::shared_ptr< AbstractStimulusFunction pIntracellularStimulus 
)

Create a new cardiac cell.

Note
subclasses must call Init() in their constructor after setting mpSystemInfo.
Parameters
pSolvernot used for these cells (they're always solved with CVODE); may be empty
numberOfStateVariablesthe size of the ODE system modelling this cell
voltageIndexthe index of the transmembrane potential within the vector of state variables
pIntracellularStimulusthe intracellular stimulus current

Definition at line 48 of file AbstractCvodeCell.cpp.

◆ ~AbstractCvodeCell()

AbstractCvodeCell::~AbstractCvodeCell ( )
virtual

Free the state variables, if they have been set.

Definition at line 59 of file AbstractCvodeCell.cpp.

Member Function Documentation

◆ Compute()

OdeSolution AbstractCvodeCell::Compute ( double  tStart,
double  tEnd,
double  tSamp = 0.0 
)
virtual

Simulates this cell's behaviour between the time interval [tStart, tEnd], and return state variable values.

The maximum time step to use will be taken as mMaxDt. If this is unset it is the same as tSamp, which defaults to HeartConfig::Instance()->GetPrintingTimeStep().

Parameters
tStartbeginning of the time interval to simulate
tEndend of the time interval to simulate
tSampsampling interval for returned results (defaults to HeartConfig printing time step)
Returns
solution object

Implements AbstractCardiacCellInterface.

Definition at line 102 of file AbstractCvodeCell.cpp.

References DOUBLE_UNSET, HeartConfig::GetPrintingTimeStep(), HeartConfig::Instance(), mMaxDt, SetTimestep(), and AbstractCvodeSystem::Solve().

◆ ComputeExceptVoltage()

void AbstractCvodeCell::ComputeExceptVoltage ( double  tStart,
double  tEnd 
)
virtual

Simulates this cell's behaviour between the time interval [tStart, tEnd], but does not update the voltage.

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

Implements AbstractCardiacCellInterface.

Definition at line 116 of file AbstractCvodeCell.cpp.

References GetVoltage(), SetVoltage(), AbstractCardiacCellInterface::SetVoltageDerivativeToZero(), SolveAndUpdateState(), and AbstractParameterisedSystem< N_Vector >::VerifyStateVariables().

Referenced by AbstractCardiacTissue< SPACE_DIM >::SetConductivityModifier().

◆ GetAnyVariable()

double AbstractCvodeCell::GetAnyVariable ( const std::string &  rName,
double  time = 0.0 
)
virtual

This just calls the method AbstractCvodeSystem::GetAnyVariable

It is here (despite being inherited) because we seem to need to specify explicitly which method in the parent classes we intend to implement to take care of the pure definition in AbstractCardiacCellInterface.

Parameters
rNamevariable name
timethe time at which to evaluate variable (only needed for derived quantities).
Returns
value of the variable at that time

Implements AbstractCardiacCellInterface.

Definition at line 187 of file AbstractCvodeCell.cpp.

References AbstractParameterisedSystem< N_Vector >::GetAnyVariable().

◆ GetNumberOfParameters()

unsigned AbstractCvodeCell::GetNumberOfParameters ( ) const
virtual

This just returns the number of parameters in the cell model.

It is here because we seem to need to specify explicitly which method in the parent classes we intend to implement to take care of the pure definition in AbstractCardiacCellInterface

Returns
the number of parameters

Implements AbstractCardiacCellInterface.

Definition at line 148 of file AbstractCvodeCell.cpp.

References AbstractUntemplatedParameterisedSystem::GetNumberOfParameters().

◆ GetNumberOfStateVariables()

unsigned AbstractCvodeCell::GetNumberOfStateVariables ( ) const
virtual

This just returns the number of state variables in the cell model.

It is here because we seem to need to specify explicitly which method in the parent classes we intend to implement to take care of the pure definition in AbstractCardiacCellInterface

Returns
the number of state variables

Implements AbstractCardiacCellInterface.

Definition at line 143 of file AbstractCvodeCell.cpp.

References AbstractUntemplatedParameterisedSystem::GetNumberOfStateVariables().

Referenced by GetStdVecStateVariables().

◆ GetParameter() [1/2]

double AbstractCvodeCell::GetParameter ( const std::string &  rParameterName)
virtual

This just calls AbstractCvodeSystem::GetParameter

It is here (despite being inherited) because we seem to need to specify explicitly which method in the parent classes we intend to implement to take care of the pure definition in AbstractCardiacCellInterface.

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

Implements AbstractCardiacCellInterface.

Definition at line 192 of file AbstractCvodeCell.cpp.

References AbstractParameterisedSystem< N_Vector >::GetParameter().

◆ GetParameter() [2/2]

double AbstractCvodeCell::GetParameter ( unsigned  parameterIndex)
virtual

This is just here to show there is an alternative to the above!

It just calls the base class method.

Parameters
parameterIndexthe index of the parameter vector entry to return
Returns
the parameter value

Implements AbstractCardiacCellInterface.

Definition at line 197 of file AbstractCvodeCell.cpp.

References AbstractParameterisedSystem< N_Vector >::GetParameter().

◆ GetStdVecStateVariables()

std::vector< double > AbstractCvodeCell::GetStdVecStateVariables ( )
virtual

This just returns the state variables in the cell model.

It is here (despite being inherited) because we seem to need to specify explicitly which method in the parent classes we intend to implement to take care of the pure definition in AbstractCardiacCellInterface.

Returns
the state variables

Implements AbstractCardiacCellInterface.

Definition at line 153 of file AbstractCvodeCell.cpp.

References CopyToStdVector(), GetNumberOfStateVariables(), and AbstractParameterisedSystem< N_Vector >::rGetStateVariables().

◆ GetTimestep()

double AbstractCvodeCell::GetTimestep ( )
Returns
The maximum timestep that is used by CVODE with this cell.

Definition at line 87 of file AbstractCvodeCell.cpp.

References mMaxDt.

◆ GetVoltage()

double AbstractCvodeCell::GetVoltage ( )
virtual
Returns
the current value of the transmembrane potential, as given in our state variable vector.

Implements AbstractCardiacCellInterface.

Definition at line 63 of file AbstractCvodeCell.cpp.

References AbstractParameterisedSystem< N_Vector >::GetAnyVariable(), AbstractParameterisedSystem< N_Vector >::mStateVariables, and AbstractCardiacCellInterface::mVoltageIndex.

Referenced by ComputeExceptVoltage().

◆ rGetStateVariableNames()

const std::vector< std::string > & AbstractCvodeCell::rGetStateVariableNames ( ) const
virtual

Just calls AbstractCvodeSystem::rGetStateVariableNames().

It is here (despite being inherited) because we seem to need to specify explicitly which method in the parent classes we intend to implement to take care of the pure definition in AbstractCardiacCellInterface.

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

Implements AbstractCardiacCellInterface.

Definition at line 160 of file AbstractCvodeCell.cpp.

References AbstractUntemplatedParameterisedSystem::rGetStateVariableNames().

◆ serialize()

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

Archive the member variables.

Parameters
archive
version

Definition at line 90 of file AbstractCvodeCell.hpp.

References mMaxDt.

◆ SetMaxTimestep()

void AbstractCvodeCell::SetMaxTimestep ( double  maxDt)

Set the maximum timestep to use for simulating this cell.

Parameters
maxDtthe maximum timestep

Definition at line 77 of file AbstractCvodeCell.cpp.

References SetTimestep().

◆ SetParameter()

void AbstractCvodeCell::SetParameter ( const std::string &  rParameterName,
double  value 
)
virtual

This just calls AbstractCvodeSystem::SetParameter

It is here (despite being inherited) because we seem to need to specify explicitly which method in the parent classes we intend to implement to take care of the pure definition in AbstractCardiacCellInterface.

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

Implements AbstractCardiacCellInterface.

Definition at line 202 of file AbstractCvodeCell.cpp.

References AbstractParameterisedSystem< N_Vector >::SetParameter().

Referenced by AbstractCvodeCellWithDataClamp::TurnOffDataClamp(), and AbstractCvodeCellWithDataClamp::TurnOnDataClamp().

◆ SetStateVariable() [1/2]

void AbstractCvodeCell::SetStateVariable ( const std::string &  rName,
double  newValue 
)
virtual

This just calls the method AbstractCvodeSystem::SetStateVariable

It is here (despite being inherited) because we seem to need to specify explicitly which method in the parent classes we intend to implement to take care of the pure definition in AbstractCardiacCellInterface.

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

Implements AbstractCardiacCellInterface.

Definition at line 182 of file AbstractCvodeCell.cpp.

References AbstractParameterisedSystem< N_Vector >::SetStateVariable().

◆ SetStateVariable() [2/2]

void AbstractCvodeCell::SetStateVariable ( unsigned  index,
double  newValue 
)
virtual

This just calls the method AbstractCvodeSystem::SetStateVariable

It is here (despite being inherited) because we seem to need to specify explicitly which method in the parent classes we intend to implement to take care of the pure definition in AbstractCardiacCellInterface.

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

Implements AbstractCardiacCellInterface.

Definition at line 177 of file AbstractCvodeCell.cpp.

References AbstractParameterisedSystem< N_Vector >::SetStateVariable().

◆ SetStateVariables() [1/2]

void AbstractCvodeCell::SetStateVariables ( const N_Vector rVariables)

This is also needed now just to show there is an alternative to the above method!

Parameters
rVariablesthe state variables (to take a copy of).

Definition at line 172 of file AbstractCvodeCell.cpp.

References AbstractParameterisedSystem< N_Vector >::SetStateVariables().

◆ SetStateVariables() [2/2]

void AbstractCvodeCell::SetStateVariables ( const std::vector< double > &  rVariables)
virtual

This just sets the state variables in the cell model.

It is here (despite being inherited) because we seem to need to specify explicitly which method in the parent classes we intend to implement to take care of the pure definition in AbstractCardiacCellInterface.

Parameters
rVariablesthe state variables (to take a copy of).

Implements AbstractCardiacCellInterface.

Definition at line 165 of file AbstractCvodeCell.cpp.

References DeleteVector(), MakeNVector(), and AbstractParameterisedSystem< N_Vector >::SetStateVariables().

◆ SetTimestep()

void AbstractCvodeCell::SetTimestep ( double  maxDt)
virtual

Set the maximum timestep to use for simulating this cell.

As CVODE adaptively alters the timestep used, this method just sets the maximum timestep allowed (despite its name). It is required as our base class AbstractCardiacCellInterface declares it as a pure virtual method. Users using this class directly should call SetMaxTimestep instead of this method, for clearer code.

Parameters
maxDtthe maximum timestep

Implements AbstractCardiacCellInterface.

Definition at line 82 of file AbstractCvodeCell.cpp.

References mMaxDt.

Referenced by Compute(), SetMaxTimestep(), and SolveAndUpdateState().

◆ SetVoltage()

void AbstractCvodeCell::SetVoltage ( double  voltage)
virtual

◆ SetVoltageDerivativeToZero()

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

Set whether to clamp the voltage by setting its derivative to zero. We need to ensure CVODE is re-initialised if this setting changes.

Parameters
clampwhether to clamp

Reimplemented from AbstractCardiacCellInterface.

Definition at line 134 of file AbstractCvodeCell.cpp.

References AbstractCardiacCellInterface::mSetVoltageDerivativeToZero, AbstractCvodeSystem::ResetSolver(), and AbstractCardiacCellInterface::SetVoltageDerivativeToZero().

◆ SolveAndUpdateState()

void AbstractCvodeCell::SolveAndUpdateState ( double  tStart,
double  tEnd 
)
virtual

Simulate this cell's behaviour between the time interval [tStart, tEnd], updating the internal state variable values.

The maximum time step to use is given by mMaxDt, which defaults to HeartConfig::Instance()->GetPrintingTimeStep() if unset.

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

Implements AbstractCardiacCellInterface.

Definition at line 93 of file AbstractCvodeCell.cpp.

References DOUBLE_UNSET, HeartConfig::Instance(), mMaxDt, SetTimestep(), and AbstractCvodeSystem::Solve().

Referenced by ComputeExceptVoltage().

Friends And Related Symbol Documentation

◆ boost::serialization::access

friend class boost::serialization::access
friend

Needed for serialization.

Definition at line 82 of file AbstractCvodeCell.hpp.

Member Data Documentation

◆ mMaxDt

double AbstractCvodeCell::mMaxDt
protected

The maximum timestep to use.

Definition at line 100 of file AbstractCvodeCell.hpp.

Referenced by Compute(), GetTimestep(), serialize(), SetTimestep(), and SolveAndUpdateState().


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