Chaste
Release::3.4
|
#include <AbstractCvodeSystem.hpp>
Public Member Functions | |
AbstractCvodeSystem (unsigned numberOfStateVariables) | |
virtual | ~AbstractCvodeSystem () |
virtual void | EvaluateYDerivatives (realtype time, const N_Vector y, N_Vector ydot)=0 |
virtual void | EvaluateAnalyticJacobian (long int N, 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) |
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_Vector & | rGetStateVariables () |
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 AbstractOdeSystemInformation > | GetSystemInformation () 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 Member Functions | |
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_VectorY) |
std::string | DumpState (const std::string &rMessage, N_VectorY, double time) |
void | CheckParametersOnLoad (const std::vector< double > &rParameters, const std::vector< std::string > &rParameterNames) |
Protected Attributes | |
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 < AbstractOdeSystemInformation > | mpSystemInfo |
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, const double &rTime) |
Private Attributes | |
N_Vector | mLastSolutionState |
double | mLastSolutionTime |
bool | mForceReset |
bool | mForceMinimalReset |
Friends | |
class | TestAbstractCvodeSystem |
class | boost::serialization::access |
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.
You can also specify an EvaluateAnalyticJacobian() method which will allow CVODE to evaluate the Jacobian matrix precisely instead of using multiple calls to EvaluateYDerivatives() to make an approximation to it. This generally provides extra accuracy and better convergence, and so gives a speed up for a given tolerance. Cardiac action potential model Jacobians are calculated automatically by PyCML (see python/pycml) via Maple for symbolic differentiation.
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 no longer set up and delete CVODE memory, unless the following method is called:
SetForceReset(true) - reset each time Solve() is called
default - reset if state variables change, or we ask to solve from a different time than the last solve call finished.
SetMinimalReset(true) - ignore changes in state vars and just reset if the time is inconsistent.
Definition at line 119 of file AbstractCvodeSystem.hpp.
AbstractCvodeSystem::AbstractCvodeSystem | ( | unsigned | numberOfStateVariables | ) |
Constructor.
numberOfStateVariables | the number of state variables in the ODE system |
Definition at line 141 of file AbstractCvodeSystem.cpp.
References SetTolerances().
|
virtual |
Virtual destructor since we have virtual methods.
Definition at line 174 of file AbstractCvodeSystem.cpp.
References DeleteVector(), FreeCvodeMemory(), mLastSolutionState, AbstractParameterisedSystem< N_Vector >::mParameters, and AbstractParameterisedSystem< N_Vector >::mStateVariables.
|
private |
Report an error from CVODE.
flag | CVODE error code |
msg | Our description of the error |
rTime | The time that the solver got to |
Definition at line 469 of file AbstractCvodeSystem.cpp.
References EXCEPTION, FreeCvodeMemory(), MakeStdVec(), mpCvodeMem, AbstractParameterisedSystem< N_Vector >::mStateVariables, and AbstractUntemplatedParameterisedSystem::rGetStateVariableNames().
Referenced by Solve().
|
inlinevirtual |
This method is called by AbstractCvodeSystemJacAdaptor method in the .cpp file.
It provides an interface between the methods different versions of CVODE are expecting and the Jacobians provided by Chaste CVODE systems.
N | the size of the ODE system |
time | the current time (used by ODE systems like y' = f(t,y) only I guess) |
y | the current state variables y for y' = f(t,y) |
ydot | the current set of derivatives y' = f(t,y) |
jacobian | a pointer to a jacobian, populated by this method. |
tmp1 | working memory of the correct size provided by CVODE for temporary calculations |
tmp2 | working memory of the correct size provided by CVODE for temporary calculations |
tmp3 | working memory of the correct size provided by CVODE for temporary calculations |
Definition at line 327 of file AbstractCvodeSystem.hpp.
References EXCEPTION.
|
pure virtual |
Method to evaluate the derivatives of the system.
time | the current time |
y | the current values of the state variables |
ydot | storage for the derivatives of the system; will be filled in on return |
void AbstractCvodeSystem::ForceUseOfNumericalJacobian | ( | bool | useNumericalJacobian = true | ) |
Force the use of a numerical Jacobian, even if an analytic form is provided. This is needed for a handful of troublesome models.
useNumericalJacobian | Whether to use a numerical instead of the analytic Jacobian. |
Definition at line 519 of file AbstractCvodeSystem.cpp.
References EXCEPTION, FreeCvodeMemory(), mHasAnalyticJacobian, and mUseAnalyticJacobian.
|
private |
Free CVODE memory when finished with.
Definition at line 459 of file AbstractCvodeSystem.cpp.
References mpCvodeMem.
Referenced by CvodeError(), ForceUseOfNumericalJacobian(), and ~AbstractCvodeSystem().
double AbstractCvodeSystem::GetAbsoluteTolerance | ( | ) |
Definition at line 324 of file AbstractCvodeSystem.cpp.
References mAbsTol.
double AbstractCvodeSystem::GetLastStepSize | ( | ) |
Definition at line 329 of file AbstractCvodeSystem.cpp.
References mLastInternalStepSize.
long int AbstractCvodeSystem::GetMaxSteps | ( | ) |
Definition at line 307 of file AbstractCvodeSystem.cpp.
References mMaxSteps.
double AbstractCvodeSystem::GetRelativeTolerance | ( | ) |
Definition at line 319 of file AbstractCvodeSystem.cpp.
References mRelTol.
bool AbstractCvodeSystem::GetUseAnalyticJacobian | ( | ) | const |
Definition at line 513 of file AbstractCvodeSystem.cpp.
References mUseAnalyticJacobian.
bool AbstractCvodeSystem::HasAnalyticJacobian | ( | ) | const |
Definition at line 508 of file AbstractCvodeSystem.cpp.
References mHasAnalyticJacobian.
|
protected |
Must be called by concrete subclass constructors to initialise the state variables, after setting mpSystemInfo.
Definition at line 162 of file AbstractCvodeSystem.cpp.
References DeleteVector(), AbstractParameterisedSystem< N_Vector >::GetInitialConditions(), AbstractParameterisedSystem< N_Vector >::mParameters, AbstractParameterisedSystem< N_Vector >::mStateVariables, and AbstractUntemplatedParameterisedSystem::rGetParameterNames().
|
inlineprivate |
Archive the member variables.
archive | the archive |
version | the current version of this class |
Definition at line 175 of file AbstractCvodeSystem.hpp.
References AbstractParameterisedSystem< N_Vector >::CheckParametersOnLoad(), CopyFromStdVector(), mAbsTol, mForceMinimalReset, mForceReset, mHasAnalyticJacobian, mLastInternalStepSize, mLastSolutionTime, mMaxSteps, AbstractUntemplatedParameterisedSystem::mNumberOfStateVariables, mRelTol, AbstractParameterisedSystem< N_Vector >::mStateVariables, and mUseAnalyticJacobian.
|
private |
Record where the last solve got to so we know whether to re-initialise.
stopTime | the finishing time |
Definition at line 438 of file AbstractCvodeSystem.cpp.
References CreateVectorIfEmpty(), AbstractUntemplatedParameterisedSystem::GetNumberOfStateVariables(), GetVectorComponent(), mForceReset, 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 353 of file AbstractCvodeSystem.cpp.
References DeleteVector(), and mLastSolutionState.
Referenced by SetForceReset(), SetTolerances(), and AbstractCvodeCell::SetVoltageDerivativeToZero().
|
inlineprivate |
Archive the member variables.
archive | the archive |
version | the current version of this class |
Definition at line 132 of file AbstractCvodeSystem.hpp.
References mAbsTol, MakeStdVec(), mForceMinimalReset, mForceReset, mHasAnalyticJacobian, mLastInternalStepSize, mLastSolutionTime, mMaxSteps, AbstractUntemplatedParameterisedSystem::mNumberOfStateVariables, AbstractParameterisedSystem< N_Vector >::mParameters, mRelTol, AbstractParameterisedSystem< N_Vector >::mStateVariables, mUseAnalyticJacobian, and AbstractUntemplatedParameterisedSystem::rGetParameterNames().
void AbstractCvodeSystem::SetForceReset | ( | 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. You can safely change parameters between solve calls with or without resets.
See also ResetSolver and SetMinimalReset
autoReset | whether to reset on every Solve |
Definition at line 334 of file AbstractCvodeSystem.cpp.
References mForceReset, 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).
numSteps | new maximum |
Definition at line 302 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 called with true argument, will call SetForceReset(false). You can safely change parameters between solve calls with or without resets.
minimalReset | whether to avoid checking for changes in state variables |
Definition at line 343 of file AbstractCvodeSystem.cpp.
References mForceMinimalReset, and SetForceReset().
Set relative and absolute tolerances; both scalars. If no parameters are given, tolerances will be reset to default values.
relTol | the relative tolerance for the solver (defaults to 1e-5) |
absTol | the absolute tolerance for the solver (defaults to 1e-7) |
Definition at line 312 of file AbstractCvodeSystem.cpp.
References mAbsTol, mRelTol, and ResetSolver().
Referenced by AbstractCvodeSystem().
|
private |
Set up the CVODE data structures needed to solve the given system from a given point.
initialConditions | initial conditions |
tStart | start time of simulation |
maxDt | maximum time step to take |
Definition at line 359 of file AbstractCvodeSystem.cpp.
References EXCEPTION, AbstractUntemplatedParameterisedSystem::GetNumberOfStateVariables(), GetVectorComponent(), mAbsTol, mForceMinimalReset, mForceReset, mLastSolutionState, mLastSolutionTime, mMaxSteps, mpCvodeMem, mRelTol, AbstractParameterisedSystem< N_Vector >::mStateVariables, mUseAnalyticJacobian, 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.
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 189 of file AbstractCvodeSystem.cpp.
References TimeStepper::AdvanceOneTimeStep(), CvodeError(), TimeStepper::EstimateTimeSteps(), TimeStepper::GetNextTime(), TimeStepper::GetTime(), TimeStepper::GetTotalTimeStepsTaken(), TimeStepper::IsTimeAtEnd(), MakeStdVec(), mLastInternalStepSize, mpCvodeMem, AbstractUntemplatedParameterisedSystem::mpSystemInfo, AbstractParameterisedSystem< N_Vector >::mStateVariables, RecordStoppingPoint(), OdeSolution::rGetSolutions(), OdeSolution::rGetTimes(), OdeSolution::SetNumberOfTimeSteps(), OdeSolution::SetOdeSystemInformation(), SetupCvode(), UNUSED_OPT, 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.
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 249 of file AbstractCvodeSystem.cpp.
References CvodeError(), mLastInternalStepSize, mpCvodeMem, AbstractParameterisedSystem< N_Vector >::mStateVariables, RecordStoppingPoint(), SetupCvode(), UNUSED_OPT, and AbstractParameterisedSystem< N_Vector >::VerifyStateVariables().
|
protected |
Absolute tolerance for solver.
Definition at line 265 of file AbstractCvodeSystem.hpp.
Referenced by GetAbsoluteTolerance(), load(), save(), SetTolerances(), and SetupCvode().
|
private |
Whether to ignore changes in the state variables when deciding whether to reset.
Definition at line 251 of file AbstractCvodeSystem.hpp.
Referenced by load(), save(), SetMinimalReset(), and SetupCvode().
|
private |
Whether to automatically reset CVODE on each Solve call.
Definition at line 248 of file AbstractCvodeSystem.hpp.
Referenced by load(), RecordStoppingPoint(), save(), SetForceReset(), and SetupCvode().
|
protected |
Whether we have an analytic Jacobian.
Definition at line 256 of file AbstractCvodeSystem.hpp.
Referenced by ForceUseOfNumericalJacobian(), HasAnalyticJacobian(), load(), and save().
|
protected |
The size of the previous timestep.
Definition at line 277 of file AbstractCvodeSystem.hpp.
Referenced by GetLastStepSize(), load(), save(), and Solve().
|
private |
Remember where the last solve got to so we know whether to re-initialise.
Definition at line 242 of file AbstractCvodeSystem.hpp.
Referenced by RecordStoppingPoint(), ResetSolver(), SetupCvode(), and ~AbstractCvodeSystem().
|
private |
Remember where the last solve got to so we know whether to re-initialise.
Definition at line 245 of file AbstractCvodeSystem.hpp.
Referenced by load(), RecordStoppingPoint(), save(), and SetupCvode().
|
protected |
The maximum number of steps to be taken by the solver in its attempt to reach the next output time.
Definition at line 274 of file AbstractCvodeSystem.hpp.
Referenced by GetMaxSteps(), load(), save(), SetMaxSteps(), and SetupCvode().
|
protected |
CVODE's internal data.
Definition at line 268 of file AbstractCvodeSystem.hpp.
Referenced by CvodeError(), FreeCvodeMemory(), SetupCvode(), and Solve().
|
protected |
Relative tolerance for solver.
Definition at line 262 of file AbstractCvodeSystem.hpp.
Referenced by GetRelativeTolerance(), load(), save(), SetTolerances(), and SetupCvode().
|
protected |
Whether to use an analytic Jacobian.
Definition at line 259 of file AbstractCvodeSystem.hpp.
Referenced by ForceUseOfNumericalJacobian(), GetUseAnalyticJacobian(), load(), save(), and SetupCvode().