Chaste Release::3.1
AbstractCvodeSystem Class Reference

#include <AbstractCvodeSystem.hpp>

Inheritance diagram for AbstractCvodeSystem:
Collaboration diagram for AbstractCvodeSystem:

List of all members.

Public Member Functions

 AbstractCvodeSystem (unsigned numberOfStateVariables)
virtual ~AbstractCvodeSystem ()
virtual void EvaluateYDerivatives (realtype time, N_Vector y, N_Vector ydot)=0
void SetAutoReset (bool autoReset)
void SetMinimalReset (bool minimalReset)
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 ()

Protected Member Functions

void Init ()

Protected Attributes

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

Private Member Functions

template<class Archive >
void save (Archive &archive, const unsigned int version) const
template<class Archive >
void load (Archive &archive, const unsigned int version)
void SetupCvode (N_Vector initialConditions, realtype tStart, realtype maxDt)
void RecordStoppingPoint (double stopTime)
void FreeCvodeMemory ()
void CvodeError (int flag, const char *msg)

Private Attributes

N_Vector mLastSolutionState
double mLastSolutionTime
bool mAutoReset
bool mForceMinimalReset

Friends

class boost::serialization::access

Detailed Description

Abstract OdeSystem class for Cvode systems (N_Vector instead of std::vector)

Sets up variables and functions for a general CVODE system.

ODE systems are specified primarily by the EvaluateYDerivatives() method, which calculates the right-hand side of the system.

Instances can store their state internally in the mStateVariables vector in our base class AbstractParameterisedSystem (see also GetNumberOfStateVariables(), SetStateVariables() and rGetStateVariables()), although this is not essential - the vector may be empty, in which case AbstractIvpOdeSolver::SolveAndUpdateStateVariable may not be used to solve the system.

CVODE systems may also have a vector of parameters, which can be accessed through the GetParameter() and SetParameter() methods of our base class.

Information about what the parameters and state variables represent is provided by a subclass of AbstractOdeSystemInformation. Various wrapper methods (e.g. rGetStateVariableNames()) are provided in our base class to access this information.

Also, subclasses may define a condition at which ODE solvers should stop prematurely. For this class CVODE solvers are being used, so CalculateRootFunction() should be used to detect the stopping time.

Note that the default tolerances for the solver are set by SetTolerances(), these can make quite a difference to the time it takes to solve the ODE system.

Repeated calls to Solve will set up and delete CVODE memory, unless the following methods are called:

SetAutoReset(false) - try not to reset unless time or state variables change

SetMinimalReset(true) - calls SetAutoReset(false) and also ignores changes in state vars.

Definition at line 101 of file AbstractCvodeSystem.hpp.


Constructor & Destructor Documentation

AbstractCvodeSystem::AbstractCvodeSystem ( unsigned  numberOfStateVariables)

Constructor.

Parameters:
numberOfStateVariablesthe number of state variables in the ODE system

Definition at line 80 of file AbstractCvodeSystem.cpp.

References SetTolerances().

AbstractCvodeSystem::~AbstractCvodeSystem ( ) [virtual]

Member Function Documentation

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

Report an error from CVODE.

Parameters:
flagCVODE error code
msgOur description of the error

Definition at line 357 of file AbstractCvodeSystem.cpp.

References EXCEPTION.

Referenced by Solve().

virtual void AbstractCvodeSystem::EvaluateYDerivatives ( realtype  time,
N_Vector  y,
N_Vector  ydot 
) [pure virtual]

Method to evaluate the derivatives of the system.

Parameters:
timethe current time
ythe current values of the state variables
ydotstorage for the derivatives of the system; will be filled in on return
void AbstractCvodeSystem::FreeCvodeMemory ( ) [private]

Free CVODE memory when finished with.

Definition at line 347 of file AbstractCvodeSystem.cpp.

References mpCvodeMem.

Referenced by Solve(), and ~AbstractCvodeSystem().

double AbstractCvodeSystem::GetAbsoluteTolerance ( )

Get the absolute tolerance.

Definition at line 229 of file AbstractCvodeSystem.cpp.

References mAbsTol.

double AbstractCvodeSystem::GetLastStepSize ( )

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

Definition at line 234 of file AbstractCvodeSystem.cpp.

References mLastInternalStepSize.

long int AbstractCvodeSystem::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 212 of file AbstractCvodeSystem.cpp.

References mMaxSteps.

double AbstractCvodeSystem::GetRelativeTolerance ( )

Get the relative tolerance.

Definition at line 224 of file AbstractCvodeSystem.cpp.

References mRelTol.

void AbstractCvodeSystem::Init ( ) [protected]
template<class Archive >
void AbstractCvodeSystem::load ( Archive &  archive,
const unsigned int  version 
) [inline, private]
void AbstractCvodeSystem::RecordStoppingPoint ( double  stopTime) [private]

Record where the last solve got to so we know whether to re-initialise.

Parameters:
stopTimethe finishing time

Definition at line 332 of file AbstractCvodeSystem.cpp.

References CreateVectorIfEmpty(), AbstractUntemplatedParameterisedSystem::GetNumberOfStateVariables(), GetVectorComponent(), mAutoReset, mLastSolutionState, mLastSolutionTime, AbstractParameterisedSystem< N_Vector >::mStateVariables, and SetVectorComponent().

Referenced by Solve().

void AbstractCvodeSystem::ResetSolver ( )

Successive calls to Solve will attempt to intelligently determine whether to re-initialise the internal CVODE solver, or whether we are simply extending the previous solution forward in time. This mechanism compares the state vector to its previous value, and the start time to the end of the last solve, which captures most cases where re-initialisation is required. However, changes to the RHS function can also require this, and cannot be automatically detected. In such cases users must call this function to force re-initialisation.

Definition at line 258 of file AbstractCvodeSystem.cpp.

References DeleteVector(), and mLastSolutionState.

Referenced by SetAutoReset(), SetTolerances(), and AbstractCvodeCell::SetVoltageDerivativeToZero().

template<class Archive >
void AbstractCvodeSystem::save ( Archive &  archive,
const unsigned int  version 
) const [inline, private]
void AbstractCvodeSystem::SetAutoReset ( bool  autoReset)

Set whether to automatically re-initialise CVODE on every call to Solve, or whether to attempt to guess when re-initialisation is needed. For example it will re-initialise if the time changes, or any state variables change.

See also ResetSolver and SetMinimalReset

Parameters:
autoResetwhether to reset on every Solve

Definition at line 239 of file AbstractCvodeSystem.cpp.

References mAutoReset, and ResetSolver().

Referenced by SetMinimalReset().

void AbstractCvodeSystem::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:
numStepsnew maximum

Definition at line 207 of file AbstractCvodeSystem.cpp.

References mMaxSteps.

void AbstractCvodeSystem::SetMinimalReset ( bool  minimalReset)

Set whether to reduce the checking done when guessing when re-initialisation is needed, so it ignores changes in the state variables. If call with true argument, will call SetAutoReset(false).

Parameters:
minimalResetwhether to avoid checking for changes in state variables

Definition at line 248 of file AbstractCvodeSystem.cpp.

References mForceMinimalReset, and SetAutoReset().

void AbstractCvodeSystem::SetTolerances ( double  relTol = 1e-5,
double  absTol = 1e-7 
)

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

Parameters:
relTolthe relative tolerance for the solver (defaults to 1e-5)
absTolthe absolute tolerance for the solver (defaults to 1e-7)

Definition at line 217 of file AbstractCvodeSystem.cpp.

References mAbsTol, mRelTol, and ResetSolver().

Referenced by AbstractCvodeSystem().

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

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

Parameters:
initialConditionsinitial conditions
tStartstart time of simulation
maxDtmaximum time step to take

Definition at line 264 of file AbstractCvodeSystem.cpp.

References EXCEPTION, AbstractUntemplatedParameterisedSystem::GetNumberOfStateVariables(), GetVectorComponent(), mAbsTol, mAutoReset, mForceMinimalReset, mLastSolutionState, mLastSolutionTime, mMaxSteps, mpCvodeMem, mRelTol, AbstractParameterisedSystem< N_Vector >::mStateVariables, and CompareDoubles::WithinAnyTolerance().

Referenced by Solve().

OdeSolution AbstractCvodeSystem::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.

Note:
See also the ResetSolver method.
Parameters:
tStartstart time of simulation
tEndend time of simulation
maxDtmaximum time step to be taken by the adaptive solver (set this appropriately to avoid missing a stimulus)
tSampsampling interval at which to store results

Definition at line 128 of file AbstractCvodeSystem.cpp.

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

Referenced by AbstractCvodeCell::Compute(), and AbstractCvodeCell::SolveAndUpdateState().

void AbstractCvodeSystem::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.

Note:
See also the ResetSolver method.
Parameters:
tStartstart time of simulation
tEndend time of simulation
maxDtmaximum time step to be taken by the adaptive solver (set this appropriately to avoid missing a stimulus)

Definition at line 180 of file AbstractCvodeSystem.cpp.

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


Member Data Documentation

Absolute tolerance for solver.

Definition at line 230 of file AbstractCvodeSystem.hpp.

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

Whether to automatically reset CVODE on each Solve call.

Definition at line 216 of file AbstractCvodeSystem.hpp.

Referenced by load(), RecordStoppingPoint(), save(), SetAutoReset(), and SetupCvode().

Whether to ignore changes in the state variables when deciding whether to reset.

Definition at line 219 of file AbstractCvodeSystem.hpp.

Referenced by load(), save(), SetMinimalReset(), and SetupCvode().

The size of the previous timestep.

Definition at line 242 of file AbstractCvodeSystem.hpp.

Referenced by GetLastStepSize(), load(), save(), and Solve().

Remember where the last solve got to so we know whether to re-initialise.

Definition at line 210 of file AbstractCvodeSystem.hpp.

Referenced by RecordStoppingPoint(), ResetSolver(), SetupCvode(), and ~AbstractCvodeSystem().

Remember where the last solve got to so we know whether to re-initialise.

Definition at line 213 of file AbstractCvodeSystem.hpp.

Referenced by load(), RecordStoppingPoint(), save(), and SetupCvode().

long int AbstractCvodeSystem::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 239 of file AbstractCvodeSystem.hpp.

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

CVODE's internal data.

Definition at line 233 of file AbstractCvodeSystem.hpp.

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

Relative tolerance for solver.

Definition at line 227 of file AbstractCvodeSystem.hpp.

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

Whether to use an analytic Jacobian.

Definition at line 224 of file AbstractCvodeSystem.hpp.

Referenced by load(), and save().


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