41#include "AbstractCvodeSystem.hpp"
42#include "CvodeAdaptor.hpp"
45#include "TimeStepper.hpp"
49#include <cvode/cvode.h>
50#include <sundials/sundials_nvector.h>
52#if CHASTE_SUNDIALS_VERSION >= 30000
53#include <cvode/cvode_direct.h>
54#include <sundials/sundials_types.h>
55#include <sunlinsol/sunlinsol_dense.h>
56#include <sunmatrix/sunmatrix_dense.h>
58#include <cvode/cvode_dense.h>
61#if CHASTE_SUNDIALS_VERSION >= 60000
62#include "CvodeContextManager.hpp"
85int AbstractCvodeSystemRhsAdaptor(realtype t,
N_Vector y,
N_Vector ydot,
void* pData)
87 assert(pData !=
nullptr);
95#if CHASTE_SUNDIALS_VERSION <= 20300
97 if (e.CheckShortMessageContains(
"is outside the times stored in the data clamp") ==
"")
103 std::cerr <<
"CVODE RHS Exception: " << e.GetMessage()
128#if CHASTE_SUNDIALS_VERSION >= 30000
130int AbstractCvodeSystemJacAdaptor(realtype t,
N_Vector y,
N_Vector ydot, CHASTE_CVODE_DENSE_MATRIX jacobian,
131#elif CHASTE_SUNDIALS_VERSION >= 20500
133int AbstractCvodeSystemJacAdaptor(
long int N, realtype t,
N_Vector y,
N_Vector ydot, CHASTE_CVODE_DENSE_MATRIX jacobian,
134#elif CHASTE_SUNDIALS_VERSION >= 20400
136int AbstractCvodeSystemJacAdaptor(
int N, realtype t,
N_Vector y,
N_Vector ydot, DlsMat jacobian,
139int AbstractCvodeSystemJacAdaptor(
long int N, DenseMat jacobian, realtype t,
N_Vector y,
N_Vector ydot,
143 assert(pData !=
nullptr);
151 std::cerr <<
"CVODE Jacobian Exception: " << e.GetMessage() << std::endl
160 mLastSolutionState(nullptr),
161 mLastSolutionTime(0.0),
162#if CHASTE_SUNDIALS_VERSION >= 20400
169 mForceMinimalReset(false),
170#if CHASTE_SUNDIALS_VERSION >= 30000
171 mpSundialsDenseMatrix(nullptr),
172 mpSundialsLinearSolver(nullptr),
174 mHasAnalyticJacobian(false),
175 mUseAnalyticJacobian(false),
178 mLastInternalStepSize(0)
188#if CHASTE_SUNDIALS_VERSION >= 60000
219 assert(tEnd >= tStart);
238 assert(ierr == CV_SUCCESS);
245 double cvode_stopped_at = stepper.
GetTime();
247 &cvode_stopped_at, CV_NORMAL);
254 assert(fabs(cvode_stopped_at - stepper.
GetNextTime()) < DBL_EPSILON);
260 solutions.
rGetTimes().push_back(cvode_stopped_at);
268 assert(ierr == CV_SUCCESS);
280 assert(tEnd >= tStart);
285 int ierr = CVodeSetStopTime(
mpCvodeMem, tEnd);
286 assert(ierr == CV_SUCCESS);
289 double cvode_stopped_at = tStart;
294 CvodeError(ierr,
"CVODE failed to solve system", cvode_stopped_at, tStart, tEnd);
297 assert(fabs(cvode_stopped_at - tEnd) < DBL_EPSILON);
300 assert(ierr == CV_SUCCESS);
400 assert(maxDt >= 0.0);
408 for (
unsigned i = 0; i < size; i++)
421#if CHASTE_SUNDIALS_VERSION >= 60000
422 mpCvodeMem = CVodeCreate(CV_BDF, CvodeContextManager::Instance()->GetSundialsContext());
423#elif CHASTE_SUNDIALS_VERSION >= 40000
434 CVodeSetErrHandlerFn(
mpCvodeMem, CvodeErrorHandler,
nullptr);
436#if CHASTE_SUNDIALS_VERSION >= 20400
442#if CHASTE_SUNDIALS_VERSION >= 20400
443 CVodeInit(
mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions);
446 CVodeMalloc(
mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
450#if CHASTE_SUNDIALS_VERSION >= 60000
452 mpSundialsDenseMatrix = SUNDenseMatrix(NV_LENGTH_S(initialConditions), NV_LENGTH_S(initialConditions), CvodeContextManager::Instance()->GetSundialsContext());
453#elif CHASTE_SUNDIALS_VERSION >= 30000
455 mpSundialsDenseMatrix = SUNDenseMatrix(NV_LENGTH_S(initialConditions), NV_LENGTH_S(initialConditions));
458#if CHASTE_SUNDIALS_VERSION >= 60000
460 mpSundialsLinearSolver = SUNLinSol_Dense(initialConditions, mpSundialsDenseMatrix, CvodeContextManager::Instance()->GetSundialsContext());
463 CVodeSetLinearSolver(
mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
464#elif CHASTE_SUNDIALS_VERSION >= 40000
466 mpSundialsLinearSolver = SUNLinSol_Dense(initialConditions, mpSundialsDenseMatrix);
469 CVodeSetLinearSolver(
mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
470#elif CHASTE_SUNDIALS_VERSION >= 30000
472 mpSundialsLinearSolver = SUNDenseLinearSolver(initialConditions, mpSundialsDenseMatrix);
475 CVDlsSetLinearSolver(
mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
479 CVDense(
mpCvodeMem, NV_LENGTH_S(initialConditions));
484#if CHASTE_SUNDIALS_VERSION >= 40000
485 CVodeSetJacFn(
mpCvodeMem, AbstractCvodeSystemJacAdaptor);
486#elif CHASTE_SUNDIALS_VERSION >= 30000
487 CVDlsSetJacFn(
mpCvodeMem, AbstractCvodeSystemJacAdaptor);
488#elif CHASTE_SUNDIALS_VERSION >= 20400
489 CVDlsSetDenseJacFn(
mpCvodeMem, AbstractCvodeSystemJacAdaptor);
491 CVDenseSetJacFn(
mpCvodeMem, AbstractCvodeSystemJacAdaptor, (
void*)(
this));
498#if CHASTE_SUNDIALS_VERSION >= 20400
499 CVodeReInit(
mpCvodeMem, tStart, initialConditions);
502 CVodeReInit(
mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
534 for (
unsigned i = 0; i < size; i++)
549#if CHASTE_SUNDIALS_VERSION >= 30000
550 if (mpSundialsLinearSolver)
553 SUNLinSolFree(mpSundialsLinearSolver);
555 mpSundialsLinearSolver =
nullptr;
557 if (mpSundialsDenseMatrix)
560 SUNMatDestroy(mpSundialsDenseMatrix);
562 mpSundialsDenseMatrix =
nullptr;
567 const double& rTime,
const double& rStartTime,
const double& rEndTime)
569 std::stringstream err;
570 char* p_flag_name = CVodeGetReturnFlagName(flag);
571 err << msg <<
": " << p_flag_name;
573 if (flag == CV_LSETUP_FAIL)
575#if CHASTE_SUNDIALS_VERSION >= 20500
580 char* p_ls_flag_name;
582#if CHASTE_SUNDIALS_VERSION >= 40000
584 p_ls_flag_name = CVodeGetLinReturnFlagName(ls_flag);
585#elif CHASTE_SUNDIALS_VERSION >= 20400
587 p_ls_flag_name = CVDlsGetReturnFlagName(ls_flag);
590 p_ls_flag_name = CVDenseGetReturnFlagName(ls_flag);
592 err <<
" (LS flag=" << ls_flag <<
":" << p_ls_flag_name <<
")";
593 free(p_ls_flag_name);
596 err <<
"\nGot from time " << rStartTime <<
" to time " << rTime <<
", was supposed to finish at time " << rEndTime <<
"\n";
597 err <<
"\nState variables are now:\n";
600 for (
unsigned i = 0; i < state_vars.size(); i++)
602 err <<
"\t" << state_var_names[i] <<
"\t:\t" << state_vars[i] << std::endl;
606 std::cerr << err.str() << std::endl
625 EXCEPTION(
"Analytic Jacobian requested, but this ODE system doesn't have one. You can check this with HasAnalyticJacobian().");
#define EXCEPTION(message)
void DeleteVector(VECTOR &rVec)
void CreateVectorIfEmpty(VECTOR &rVec, unsigned size)
void SetVectorComponent(VECTOR &rVec, unsigned index, double value)
double GetVectorComponent(const VECTOR &rVec, unsigned index)
std::vector< double > MakeStdVec(N_Vector v)
void SetTolerances(double relTol=1e-5, double absTol=1e-7)
bool GetMinimalReset()
Get whether we want to run with minimal reset or not (no reinitialisation of the solver if variables ...
void SetForceReset(bool autoReset)
bool mUseAnalyticJacobian
void SetupCvode(N_Vector initialConditions, realtype tStart, realtype maxDt)
void CvodeError(int flag, const char *msg, const double &rTime, const double &rStartTime, const double &rEndTime)
double GetRelativeTolerance()
void RecordStoppingPoint(double stopTime)
void SetMaxSteps(long int numSteps)
bool HasAnalyticJacobian() const
AbstractCvodeSystem(unsigned numberOfStateVariables)
double GetAbsoluteTolerance()
OdeSolution Solve(realtype tStart, realtype tEnd, realtype maxDt, realtype tSamp)
void ForceUseOfNumericalJacobian(bool useNumericalJacobian=true)
virtual void EvaluateYDerivatives(realtype time, const N_Vector y, N_Vector ydot)=0
bool mHasAnalyticJacobian
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)
bool GetUseAnalyticJacobian() const
double mLastInternalStepSize
void SetMinimalReset(bool minimalReset)
bool GetForceReset()
Get whether we will force a solver reset on every call to Solve()
N_Vector mLastSolutionState
virtual ~AbstractCvodeSystem()
N_Vector GetInitialConditions() const
virtual void VerifyStateVariables()
const std::vector< std::string > & rGetParameterNames() const
boost::shared_ptr< AbstractOdeSystemInformation > mpSystemInfo
const std::vector< std::string > & rGetStateVariableNames() const
unsigned GetNumberOfStateVariables() const
static bool WithinAnyTolerance(double number1, double number2, double relTol=DBL_EPSILON, double absTol=DBL_EPSILON, bool printError=false)
void SetNumberOfTimeSteps(unsigned numTimeSteps)
std::vector< std::vector< double > > & rGetSolutions()
std::vector< double > & rGetTimes()
void SetOdeSystemInformation(boost::shared_ptr< const AbstractOdeSystemInformation > pOdeSystemInfo)
unsigned GetTotalTimeStepsTaken() const
void AdvanceOneTimeStep()
double GetNextTime() const
unsigned EstimateTimeSteps() const