Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
AbstractCvodeSystem Class Referenceabstract

#include <AbstractCvodeSystem.hpp>

+ Inheritance diagram for AbstractCvodeSystem:
+ Collaboration diagram for AbstractCvodeSystem:

Public Member Functions

 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 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_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)
 

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< AbstractOdeSystemInformationmpSystemInfo
 

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, const double &rStartTime, const double &rEndTime)
 

Private Attributes

N_Vector mLastSolutionState
 
double mLastSolutionTime
 
bool mForceReset
 
bool mForceMinimalReset
 

Friends

class TestAbstractCvodeSystem
 
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, and uses an inbuilt CVODE solver.

N.B. The CvodeAdaptor solver is not used by AbstractCvodeSystems, it is only used to translate when the ODE system uses the normal std::vector<double> state variables. For Chaste developers - this means there is some duplication of code between this class and CvodeAdaptor as both set up and use CVODE to solve systems. If you change one you should probably change both!

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 chaste_codegen.

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 141 of file AbstractCvodeSystem.hpp.

Constructor & Destructor Documentation

◆ AbstractCvodeSystem()

AbstractCvodeSystem::AbstractCvodeSystem ( unsigned  numberOfStateVariables)

Constructor.

Parameters
numberOfStateVariablesthe number of state variables in the ODE system

Definition at line 158 of file AbstractCvodeSystem.cpp.

References SetTolerances().

◆ ~AbstractCvodeSystem()

AbstractCvodeSystem::~AbstractCvodeSystem ( )
virtual

Member Function Documentation

◆ CvodeError()

void AbstractCvodeSystem::CvodeError ( int  flag,
const char *  msg,
const double rTime,
const double rStartTime,
const double rEndTime 
)
private

Report an error from CVODE.

Parameters
flagCVODE error code
msgOur description of the error
rTimeThe time that the solver got to
rStartTimeThe time that the solver started at.
rEndTimeThe time that the solver was supposed to solve up to.

Definition at line 566 of file AbstractCvodeSystem.cpp.

References EXCEPTION, FreeCvodeMemory(), MakeStdVec(), mpCvodeMem, AbstractParameterisedSystem< N_Vector >::mStateVariables, and AbstractUntemplatedParameterisedSystem::rGetStateVariableNames().

Referenced by Solve(), and Solve().

◆ EvaluateAnalyticJacobian()

virtual void AbstractCvodeSystem::EvaluateAnalyticJacobian ( realtype  time,
N_Vector  y,
N_Vector  ydot,
CHASTE_CVODE_DENSE_MATRIX  jacobian,
N_Vector  tmp1,
N_Vector  tmp2,
N_Vector  tmp3 
)
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.

Parameters
timethe current time (used by ODE systems like y' = f(t,y) only I guess)
ythe current state variables y for y' = f(t,y)
ydotthe current set of derivatives y' = f(t,y)
jacobiana pointer to a jacobian, populated by this method.
tmp1working memory of the correct size provided by CVODE for temporary calculations
tmp2working memory of the correct size provided by CVODE for temporary calculations
tmp3working memory of the correct size provided by CVODE for temporary calculations

Definition at line 356 of file AbstractCvodeSystem.hpp.

References EXCEPTION.

◆ EvaluateYDerivatives()

virtual void AbstractCvodeSystem::EvaluateYDerivatives ( realtype  time,
const 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

◆ ForceUseOfNumericalJacobian()

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.

Parameters
useNumericalJacobianWhether to use a numerical instead of the analytic Jacobian.

Definition at line 621 of file AbstractCvodeSystem.cpp.

References EXCEPTION, FreeCvodeMemory(), mHasAnalyticJacobian, and mUseAnalyticJacobian.

◆ FreeCvodeMemory()

void AbstractCvodeSystem::FreeCvodeMemory ( )
private

Free CVODE memory when finished with.

Definition at line 541 of file AbstractCvodeSystem.cpp.

References mpCvodeMem.

Referenced by ~AbstractCvodeSystem(), CvodeError(), and ForceUseOfNumericalJacobian().

◆ GetAbsoluteTolerance()

double AbstractCvodeSystem::GetAbsoluteTolerance ( )
Returns
the absolute tolerance.

Definition at line 352 of file AbstractCvodeSystem.cpp.

References mAbsTol.

◆ GetForceReset()

bool AbstractCvodeSystem::GetForceReset ( )

Get whether we will force a solver reset on every call to Solve()

Returns
true We will reinitialise the solver at the start of every Solve() call.
false We will not reinitialise the solver at the start of every Solve() call.

Definition at line 376 of file AbstractCvodeSystem.cpp.

References mForceReset.

◆ GetLastStepSize()

double AbstractCvodeSystem::GetLastStepSize ( )
Returns
the last step size used internally by CVODE in the last Solve call.

Definition at line 357 of file AbstractCvodeSystem.cpp.

References mLastInternalStepSize.

◆ GetMaxSteps()

long int AbstractCvodeSystem::GetMaxSteps ( )
Returns
the maximum number of steps to be taken by the solver in its attempt to reach the next output time.

Definition at line 335 of file AbstractCvodeSystem.cpp.

References mMaxSteps.

◆ GetMinimalReset()

bool AbstractCvodeSystem::GetMinimalReset ( )

Get whether we want to run with minimal reset or not (no reinitialisation of the solver if variables change)

Returns
true We are not resetting the solver if time or variables change between solve calls.
false We will reset the solver if time or variables change between solve calls.

Definition at line 371 of file AbstractCvodeSystem.cpp.

References mForceMinimalReset.

◆ GetRelativeTolerance()

double AbstractCvodeSystem::GetRelativeTolerance ( )
Returns
the relative tolerance.

Definition at line 347 of file AbstractCvodeSystem.cpp.

References mRelTol.

◆ GetUseAnalyticJacobian()

bool AbstractCvodeSystem::GetUseAnalyticJacobian ( ) const
Returns
whether an analytic Jacobian is used (mUseAnalyticJacobian).

Definition at line 616 of file AbstractCvodeSystem.cpp.

References mUseAnalyticJacobian.

◆ HasAnalyticJacobian()

bool AbstractCvodeSystem::HasAnalyticJacobian ( ) const
Returns
whether the ODE system has an analytic Jacobian (mHasAnalyticJacobian).

Definition at line 611 of file AbstractCvodeSystem.cpp.

References mHasAnalyticJacobian.

◆ Init()

void AbstractCvodeSystem::Init ( )
protected

◆ load()

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

◆ RecordStoppingPoint()

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 520 of file AbstractCvodeSystem.cpp.

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

Referenced by Solve(), and Solve().

◆ ResetSolver()

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 390 of file AbstractCvodeSystem.cpp.

References DeleteVector(), and mLastSolutionState.

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

◆ save()

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

◆ SetForceReset()

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

Parameters
autoResetwhether to reset on every Solve

Definition at line 362 of file AbstractCvodeSystem.cpp.

References mForceReset, and ResetSolver().

Referenced by SetMinimalReset().

◆ SetMaxSteps()

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 330 of file AbstractCvodeSystem.cpp.

References mMaxSteps.

◆ SetMinimalReset()

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.

Parameters
minimalResetwhether to avoid checking for changes in state variables

Definition at line 381 of file AbstractCvodeSystem.cpp.

References mForceMinimalReset, and SetForceReset().

◆ SetTolerances()

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 340 of file AbstractCvodeSystem.cpp.

References mAbsTol, mRelTol, and ResetSolver().

Referenced by AbstractCvodeSystem().

◆ SetupCvode()

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 395 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(), and Solve().

◆ Solve() [1/2]

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 276 of file AbstractCvodeSystem.cpp.

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

◆ Solve() [2/2]

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
Returns
the solution object with results at each sampling step

Definition at line 214 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().

Friends And Related Symbol Documentation

◆ boost::serialization::access

friend class boost::serialization::access
friend

Definition at line 146 of file AbstractCvodeSystem.hpp.

◆ TestAbstractCvodeSystem

friend class TestAbstractCvodeSystem
friend

Definition at line 144 of file AbstractCvodeSystem.hpp.

Member Data Documentation

◆ mAbsTol

double AbstractCvodeSystem::mAbsTol
protected

Absolute tolerance for solver.

Definition at line 296 of file AbstractCvodeSystem.hpp.

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

◆ mForceMinimalReset

bool AbstractCvodeSystem::mForceMinimalReset
private

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

Definition at line 276 of file AbstractCvodeSystem.hpp.

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

◆ mForceReset

bool AbstractCvodeSystem::mForceReset
private

Whether to automatically reset CVODE on each Solve call.

Definition at line 273 of file AbstractCvodeSystem.hpp.

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

◆ mHasAnalyticJacobian

bool AbstractCvodeSystem::mHasAnalyticJacobian
protected

Whether we have an analytic Jacobian.

Definition at line 287 of file AbstractCvodeSystem.hpp.

Referenced by ForceUseOfNumericalJacobian(), HasAnalyticJacobian(), load(), and save().

◆ mLastInternalStepSize

double AbstractCvodeSystem::mLastInternalStepSize
protected

The size of the previous timestep.

Definition at line 308 of file AbstractCvodeSystem.hpp.

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

◆ mLastSolutionState

N_Vector AbstractCvodeSystem::mLastSolutionState
private

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

Definition at line 267 of file AbstractCvodeSystem.hpp.

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

◆ mLastSolutionTime

double AbstractCvodeSystem::mLastSolutionTime
private

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

Definition at line 270 of file AbstractCvodeSystem.hpp.

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

◆ mMaxSteps

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 305 of file AbstractCvodeSystem.hpp.

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

◆ mpCvodeMem

void* AbstractCvodeSystem::mpCvodeMem
protected

CVODE's internal data.

Definition at line 299 of file AbstractCvodeSystem.hpp.

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

◆ mRelTol

double AbstractCvodeSystem::mRelTol
protected

Relative tolerance for solver.

Definition at line 293 of file AbstractCvodeSystem.hpp.

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

◆ mUseAnalyticJacobian

bool AbstractCvodeSystem::mUseAnalyticJacobian
protected

Whether to use an analytic Jacobian.

Definition at line 290 of file AbstractCvodeSystem.hpp.

Referenced by ForceUseOfNumericalJacobian(), GetUseAnalyticJacobian(), load(), save(), and SetupCvode().


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