AbstractCvodeCell Class Reference

#include <AbstractCvodeCell.hpp>

Inheritance diagram for AbstractCvodeCell:

Inheritance graph
[legend]
Collaboration diagram for AbstractCvodeCell:

Collaboration graph
[legend]

List of all members.

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 ResetToInitialConditions ()
N_Vector GetInitialConditions ()
void SetStateVariables (N_Vector stateVars)
void SetStateVariablesUsingACopyOfThisVector (N_Vector stateVars)
N_Vector GetStateVariables ()
virtual void EvaluateRhs (realtype t, N_Vector y, N_Vector ydot)=0
OdeSolution Compute (double tStart, double tEnd, double tSamp=0.0)
void ComputeExceptVoltage (double tStart, double tEnd)
void SetVoltageDerivativeToZero (bool clamp=true)
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-4, double absTol=1e-6)
double GetRelativeTolerance ()
double GetAbsoluteTolerance ()
double GetLastStepSize ()

Protected Member Functions

std::string DumpState (const std::string &rMessage, N_Vector Y=NULL)
void Init ()
N_Vector CopyVector (N_Vector originalVec)

Protected Attributes

double mRelTol
double mAbsTol
void * mpCvodeMem
long int mMaxSteps
double mLastInternalStepSize

Private Member Functions

void SetupCvode (N_Vector initialConditions, realtype tStart, realtype maxDt)
void FreeCvodeMemory ()
void CvodeError (int flag, const char *msg)
std::vector< double > MakeStdVec (N_Vector v)


Detailed Description

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

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 doesn't work for us.

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

Note that a call to Solve will initialise the CVODE solver, and free its working memory when done. There is thus a non-trivial overhead involved.

Todo:
#890 Add an option to just initialise once, and assume subsequent Solve calls are continuing from where we left off.
Todo:
#889 Integrate this better into the main hierarchy?

Definition at line 69 of file AbstractCvodeCell.hpp.


Constructor & Destructor Documentation

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:
pSolver not used for these cells (they're always solved with CVODE); may be empty
numberOfStateVariables the size of the ODE system modelling this cell
voltageIndex the index of the transmembrane potential within the vector of state variables
pIntracellularStimulus the intracellular stimulus current

Definition at line 73 of file AbstractCvodeCell.cpp.

References SetTolerances().

AbstractCvodeCell::~AbstractCvodeCell (  )  [virtual]

Free the state variables, if they have been set.

Definition at line 86 of file AbstractCvodeCell.cpp.

References AbstractParameterisedSystem< N_Vector >::mParameters, and AbstractParameterisedSystem< N_Vector >::mStateVariables.


Member Function Documentation

std::string AbstractCvodeCell::DumpState ( const std::string &  rMessage,
N_Vector  Y = NULL 
) [protected]

Used to include extra debugging information in exception messages, e.g. EXCEPTION(DumpState("Gating variable out of range"));

Parameters:
rMessage explanatory message to include in exception
Y current state variable values (defaults to using mStateVariables)

Definition at line 375 of file AbstractCvodeCell.cpp.

References AbstractParameterisedSystem< N_Vector >::mStateVariables, and AbstractParameterisedSystem< N_Vector >::rGetStateVariableNames().

void AbstractCvodeCell::Init (  )  [protected]

Must be called by concrete subclass constructors to initialise the state variables, after setting mpSystemInfo.

Definition at line 101 of file AbstractCvodeCell.cpp.

References GetInitialConditions(), AbstractParameterisedSystem< N_Vector >::mParameters, AbstractParameterisedSystem< N_Vector >::rGetParameterNames(), and SetStateVariables().

N_Vector AbstractCvodeCell::CopyVector ( N_Vector  originalVec  )  [protected]

Create a new N_Vector by copying the given vector. Used by GetStateVariables and SetStateVariablesUsingACopyOfThisVector.

Parameters:
originalVec the vector to copy

Definition at line 391 of file AbstractCvodeCell.cpp.

Referenced by GetStateVariables(), and SetStateVariablesUsingACopyOfThisVector().

double AbstractCvodeCell::GetVoltage (  )  [virtual]

Get the current value of the transmembrane potential, as given in our state variable vector.

Implements AbstractCardiacCellInterface.

Definition at line 117 of file AbstractCvodeCell.cpp.

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

void AbstractCvodeCell::SetVoltage ( double  voltage  )  [virtual]

Set the transmembrane potential

Parameters:
voltage new value

Implements AbstractCardiacCellInterface.

Definition at line 123 of file AbstractCvodeCell.cpp.

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

void AbstractCvodeCell::ResetToInitialConditions (  )  [virtual]

Reset the model's state variables to the default initial conditions.

Implements AbstractCardiacCellInterface.

Definition at line 130 of file AbstractCvodeCell.cpp.

References GetInitialConditions(), and SetStateVariables().

N_Vector AbstractCvodeCell::GetInitialConditions (  ) 

Get the initial conditions for the cell.

Creates and returns a fresh N_Vector, which must be destroyed by the caller when finished with.

Definition at line 135 of file AbstractCvodeCell.cpp.

References AbstractParameterisedSystem< N_Vector >::mpSystemInfo.

Referenced by Init(), and ResetToInitialConditions().

void AbstractCvodeCell::SetStateVariables ( N_Vector  stateVars  ) 

Assign a vector to be used for this cell's state.

The cell takes responsibility for freeing this vector when it is destroyed. If the cell already has state, it will be freed.

Parameters:
stateVars new state variables vector

Todo:
#890 re-init CVODE here?

Definition at line 154 of file AbstractCvodeCell.cpp.

References AbstractParameterisedSystem< N_Vector >::mStateVariables.

Referenced by Init(), and ResetToInitialConditions().

void AbstractCvodeCell::SetStateVariablesUsingACopyOfThisVector ( N_Vector  stateVars  ) 

Assign a vector to be copied for this cell's state.

Caller retains responsibility for freeing the vector.

Parameters:
stateVars new state variables vector

Definition at line 164 of file AbstractCvodeCell.cpp.

References CopyVector(), and AbstractParameterisedSystem< N_Vector >::mStateVariables.

N_Vector AbstractCvodeCell::GetStateVariables (  ) 

Takes a copy of the state variable vector. Doesn't really return a vector (N_Vector is a pointer type) but named like this to match AbstractCardiacCell. Caller takes responsibility for freeing the vector.

Definition at line 147 of file AbstractCvodeCell.cpp.

References CopyVector(), AbstractParameterisedSystem< N_Vector >::mpSystemInfo, and AbstractParameterisedSystem< N_Vector >::mStateVariables.

virtual void AbstractCvodeCell::EvaluateRhs ( realtype  t,
N_Vector  y,
N_Vector  ydot 
) [pure virtual]

RHS evaluation function, to be provided by subclasses.

Parameters:
t the time at which to evaluate the RHS
y the values of the state variables at time t
ydot to be filled in with the derivatives of the state variables

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 the same as the sampling time step, which defaults to HeartConfig::Instance()->GetPrintingTimeStep(). This is a bit of a hack.

Parameters:
tStart beginning of the time interval to simulate
tEnd end of the time interval to simulate
tSamp sampling interval for returned results (defaults to dt)

Implements AbstractCardiacCellInterface.

Definition at line 182 of file AbstractCvodeCell.cpp.

References HeartConfig::GetPrintingTimeStep(), HeartConfig::Instance(), and Solve().

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:
tStart beginning of the time interval to simulate
tEnd end of the time interval to simulate

Implements AbstractCardiacCellInterface.

Definition at line 192 of file AbstractCvodeCell.cpp.

References EXCEPTION.

void AbstractCvodeCell::SetVoltageDerivativeToZero ( bool  clamp = true  ) 

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

Parameters:
clamp whether to clamp

Definition at line 176 of file AbstractCvodeCell.cpp.

References AbstractCardiacCellInterface::mSetVoltageDerivativeToZero.

OdeSolution AbstractCvodeCell::Solve ( realtype  tStart,
realtype  tEnd,
realtype  maxDt,
realtype  tSamp 
)

Simulate the cell, returning a sampling of the state variables.

Uses the current values of the state variables at initial conditions. If the state variables have not been set (either by a prior solve, or a call to SetStateVariables) the initial conditions (given by GetInitialConditions) will be used.

The final values of the state variables will also be stored in this object.

Parameters:
tStart start time of simulation
tEnd end time of simulation
maxDt maximum time step to be taken by the adaptive solver (set this appropriately to avoid missing a stimulus)
tSamp sampling interval at which to store results

Definition at line 198 of file AbstractCvodeCell.cpp.

References TimeStepper::AdvanceOneTimeStep(), CvodeError(), TimeStepper::EstimateTimeSteps(), FreeCvodeMemory(), TimeStepper::GetNextTime(), TimeStepper::GetTotalTimeStepsTaken(), TimeStepper::IsTimeAtEnd(), MakeStdVec(), mLastInternalStepSize, mpCvodeMem, AbstractParameterisedSystem< N_Vector >::mpSystemInfo, AbstractParameterisedSystem< N_Vector >::mStateVariables, OdeSolution::rGetSolutions(), OdeSolution::rGetTimes(), OdeSolution::SetNumberOfTimeSteps(), OdeSolution::SetOdeSystemInformation(), SetupCvode(), and AbstractCardiacCellInterface::VerifyStateVariables().

Referenced by Compute().

void AbstractCvodeCell::Solve ( realtype  tStart,
realtype  tEnd,
realtype  maxDt 
)

Simulate the cell, updating its internal state variables.

Uses the current values of the state variables at initial conditions. If the state variables have not been set (either by a prior solve, or a call to SetStateVariables) the initial conditions (given by GetInitialConditions) will be used.

Parameters:
tStart start time of simulation
tEnd end time of simulation
maxDt maximum time step to be taken by the adaptive solver (set this appropriately to avoid missing a stimulus)

Definition at line 255 of file AbstractCvodeCell.cpp.

References CvodeError(), FreeCvodeMemory(), mLastInternalStepSize, mpCvodeMem, AbstractParameterisedSystem< N_Vector >::mStateVariables, SetupCvode(), and AbstractCardiacCellInterface::VerifyStateVariables().

void AbstractCvodeCell::SetMaxSteps ( long int  numSteps  ) 

Change the maximum number of steps to be taken by the solver in its attempt to reach the next output time. Default is 500 (set by CVODE).

Parameters:
numSteps new maximum

Definition at line 286 of file AbstractCvodeCell.cpp.

References mMaxSteps.

long int AbstractCvodeCell::GetMaxSteps (  ) 

Get the maximum number of steps to be taken by the solver in its attempt to reach the next output time.

Definition at line 291 of file AbstractCvodeCell.cpp.

References mMaxSteps.

void AbstractCvodeCell::SetTolerances ( double  relTol = 1e-4,
double  absTol = 1e-6 
)

Set relative and absolute tolerances; both scalars. If no parameters are given, tolerances will be reset to default values.

Parameters:
relTol the relative tolerance for the solver (defaults to 1e-4)
absTol the absolute tolerance for the solver (defaults to 1e-6)

Definition at line 296 of file AbstractCvodeCell.cpp.

References mAbsTol, and mRelTol.

Referenced by AbstractCvodeCell().

double AbstractCvodeCell::GetRelativeTolerance (  ) 

Get the relative tolerance.

Definition at line 302 of file AbstractCvodeCell.cpp.

References mRelTol.

double AbstractCvodeCell::GetAbsoluteTolerance (  ) 

Get the absolute tolerance.

Definition at line 307 of file AbstractCvodeCell.cpp.

References mAbsTol.

double AbstractCvodeCell::GetLastStepSize (  ) 

Get the last step size used internally by CVODE in the last Solve call.

Definition at line 312 of file AbstractCvodeCell.cpp.

References mLastInternalStepSize.

void AbstractCvodeCell::SetupCvode ( N_Vector  initialConditions,
realtype  tStart,
realtype  maxDt 
) [private]

Set up the CVODE data structures needed to solve the given system.

Parameters:
initialConditions initial conditions
tStart start time of simulation
maxDt maximum time step to take

Definition at line 318 of file AbstractCvodeCell.cpp.

References EXCEPTION, AbstractParameterisedSystem< N_Vector >::GetNumberOfStateVariables(), mAbsTol, mMaxSteps, mpCvodeMem, and mRelTol.

Referenced by Solve().

void AbstractCvodeCell::FreeCvodeMemory (  )  [private]

Free CVODE memory after a solve.

Definition at line 345 of file AbstractCvodeCell.cpp.

References mpCvodeMem.

Referenced by Solve().

void AbstractCvodeCell::CvodeError ( int  flag,
const char *  msg 
) [private]

Report an error from CVODE.

Parameters:
flag CVODE error code
msg Our description of the error

Definition at line 353 of file AbstractCvodeCell.cpp.

References EXCEPTION.

Referenced by Solve().

std::vector< double > AbstractCvodeCell::MakeStdVec ( N_Vector  v  )  [private]

Copy an N_Vector into a std::vector<double>

Parameters:
v vector to copy

Definition at line 363 of file AbstractCvodeCell.cpp.

Referenced by Solve().


Member Data Documentation

double AbstractCvodeCell::mRelTol [protected]

Relative tolerance for solver.

Definition at line 73 of file AbstractCvodeCell.hpp.

Referenced by GetRelativeTolerance(), SetTolerances(), and SetupCvode().

double AbstractCvodeCell::mAbsTol [protected]

Absolute tolerance for solver.

Definition at line 75 of file AbstractCvodeCell.hpp.

Referenced by GetAbsoluteTolerance(), SetTolerances(), and SetupCvode().

void* AbstractCvodeCell::mpCvodeMem [protected]

CVODE's internal data.

Definition at line 78 of file AbstractCvodeCell.hpp.

Referenced by FreeCvodeMemory(), SetupCvode(), and Solve().

long int AbstractCvodeCell::mMaxSteps [protected]

The maximum number of steps to be taken by the solver in its attempt to reach the next output time.

Definition at line 83 of file AbstractCvodeCell.hpp.

Referenced by GetMaxSteps(), SetMaxSteps(), and SetupCvode().

The size of the previous timestep.

Definition at line 86 of file AbstractCvodeCell.hpp.

Referenced by GetLastStepSize(), and Solve().


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

Generated on Mon Nov 1 12:35:41 2010 for Chaste by  doxygen 1.5.5