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>
81 int AbstractCvodeSystemRhsAdaptor(realtype t,
N_Vector y,
N_Vector ydot,
void* pData)
83 assert(pData !=
nullptr);
91 #if CHASTE_SUNDIALS_VERSION <= 20300
99 std::cerr <<
"CVODE RHS Exception: " << e.
GetMessage()
122 #if CHASTE_SUNDIALS_VERSION >= 30000
124 int AbstractCvodeSystemJacAdaptor(realtype t,
N_Vector y,
N_Vector ydot, CHASTE_CVODE_DENSE_MATRIX jacobian,
125 #elif CHASTE_SUNDIALS_VERSION >= 20500
127 int AbstractCvodeSystemJacAdaptor(
long int N, realtype t,
N_Vector y,
N_Vector ydot, CHASTE_CVODE_DENSE_MATRIX jacobian,
128 #elif CHASTE_SUNDIALS_VERSION >= 20400
130 int AbstractCvodeSystemJacAdaptor(
int N, realtype t,
N_Vector y,
N_Vector ydot, DlsMat jacobian,
133 int AbstractCvodeSystemJacAdaptor(
long int N, DenseMat jacobian, realtype t,
N_Vector y,
N_Vector ydot,
137 assert(pData !=
nullptr);
145 std::cerr <<
"CVODE Jacobian Exception: " << e.
GetMessage() << std::endl
154 mLastSolutionState(nullptr),
155 mLastSolutionTime(0.0),
156 #if CHASTE_SUNDIALS_VERSION >= 20400
163 mForceMinimalReset(false),
164 #if CHASTE_SUNDIALS_VERSION >= 30000
165 mpSundialsDenseMatrix(nullptr),
166 mpSundialsLinearSolver(nullptr),
168 mHasAnalyticJacobian(false),
169 mUseAnalyticJacobian(false),
172 mLastInternalStepSize(0)
209 assert(tEnd >= tStart);
228 assert(ierr == CV_SUCCESS);
235 double cvode_stopped_at = stepper.
GetTime();
237 &cvode_stopped_at, CV_NORMAL);
244 assert(fabs(cvode_stopped_at - stepper.
GetNextTime()) < DBL_EPSILON);
250 solutions.
rGetTimes().push_back(cvode_stopped_at);
258 assert(ierr == CV_SUCCESS);
270 assert(tEnd >= tStart);
275 int ierr = CVodeSetStopTime(
mpCvodeMem, tEnd);
276 assert(ierr == CV_SUCCESS);
279 double cvode_stopped_at = tStart;
284 CvodeError(ierr,
"CVODE failed to solve system", cvode_stopped_at, tStart, tEnd);
287 assert(fabs(cvode_stopped_at - tEnd) < DBL_EPSILON);
290 assert(ierr == CV_SUCCESS);
380 assert(maxDt >= 0.0);
388 for (
unsigned i = 0; i < size; i++)
405 CVodeSetErrHandlerFn(
mpCvodeMem, CvodeErrorHandler,
nullptr);
407 #if CHASTE_SUNDIALS_VERSION >= 20400
413 #if CHASTE_SUNDIALS_VERSION >= 20400
414 CVodeInit(
mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions);
417 CVodeMalloc(
mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
421 #if CHASTE_SUNDIALS_VERSION >= 30000
423 mpSundialsDenseMatrix = SUNDenseMatrix(NV_LENGTH_S(initialConditions), NV_LENGTH_S(initialConditions));
426 mpSundialsLinearSolver = SUNDenseLinearSolver(initialConditions, mpSundialsDenseMatrix);
429 CVDlsSetLinearSolver(
mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
432 CVDense(
mpCvodeMem, NV_LENGTH_S(initialConditions));
437 #if CHASTE_SUNDIALS_VERSION >= 30000
438 CVDlsSetJacFn(
mpCvodeMem, AbstractCvodeSystemJacAdaptor);
439 #elif CHASTE_SUNDIALS_VERSION >= 20400
440 CVDlsSetDenseJacFn(
mpCvodeMem, AbstractCvodeSystemJacAdaptor);
442 CVDenseSetJacFn(
mpCvodeMem, AbstractCvodeSystemJacAdaptor, (
void*)(
this));
449 #if CHASTE_SUNDIALS_VERSION >= 20400
450 CVodeReInit(
mpCvodeMem, tStart, initialConditions);
453 CVodeReInit(
mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
481 for (
unsigned i = 0; i < size; i++)
496 #if CHASTE_SUNDIALS_VERSION >= 30000
497 if (mpSundialsLinearSolver)
500 SUNLinSolFree(mpSundialsLinearSolver);
502 mpSundialsLinearSolver =
nullptr;
504 if (mpSundialsDenseMatrix)
507 SUNMatDestroy(mpSundialsDenseMatrix);
509 mpSundialsDenseMatrix =
nullptr;
514 const double& rTime,
const double& rStartTime,
const double& rEndTime)
516 std::stringstream err;
517 char* p_flag_name = CVodeGetReturnFlagName(flag);
518 err << msg <<
": " << p_flag_name;
520 if (flag == CV_LSETUP_FAIL)
522 #if CHASTE_SUNDIALS_VERSION >= 20500
527 char* p_ls_flag_name;
528 #if CHASTE_SUNDIALS_VERSION >= 20400
530 p_ls_flag_name = CVDlsGetReturnFlagName(ls_flag);
533 p_ls_flag_name = CVDenseGetReturnFlagName(ls_flag);
535 err <<
" (LS flag=" << ls_flag <<
":" << p_ls_flag_name <<
")";
536 free(p_ls_flag_name);
539 err <<
"\nGot from time " << rStartTime <<
" to time " << rTime <<
", was supposed to finish at time " << rEndTime <<
"\n";
540 err <<
"\nState variables are now:\n";
543 for (
unsigned i = 0; i < state_vars.size(); i++)
545 err <<
"\t" << state_var_names[i] <<
"\t:\t" << state_vars[i] << std::endl;
549 std::cerr << err.str() << std::endl
568 EXCEPTION(
"Analytic Jacobian requested, but this ODE system doesn't have one. You can check this with HasAnalyticJacobian().");
629 #endif // CHASTE_CVODE
double GetAbsoluteTolerance()
void SetMinimalReset(bool minimalReset)
void SetTolerances(double relTol=1e-5, double absTol=1e-7)
std::vector< double > MakeStdVec(N_Vector v)
void CvodeError(int flag, const char *msg, const double &rTime, const double &rStartTime, const double &rEndTime)
static bool WithinAnyTolerance(double number1, double number2, double relTol=DBL_EPSILON, double absTol=DBL_EPSILON, bool printError=false)
std::vector< std::vector< double > > & rGetSolutions()
bool GetUseAnalyticJacobian() const
void RecordStoppingPoint(double stopTime)
#define EXCEPTION(message)
double mLastInternalStepSize
void SetupCvode(N_Vector initialConditions, realtype tStart, realtype maxDt)
const std::vector< std::string > & rGetStateVariableNames() const
void AdvanceOneTimeStep()
unsigned GetNumberOfStateVariables() const
N_Vector GetInitialConditions() const
void CreateVectorIfEmpty(VECTOR &rVec, unsigned size)
void SetOdeSystemInformation(boost::shared_ptr< const AbstractOdeSystemInformation > pOdeSystemInfo)
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 mUseAnalyticJacobian
bool HasAnalyticJacobian() const
void SetForceReset(bool autoReset)
virtual void EvaluateYDerivatives(realtype time, const N_Vector y, N_Vector ydot)=0
std::string GetMessage() const
void SetVectorComponent(VECTOR &rVec, unsigned index, double value)
unsigned GetTotalTimeStepsTaken() const
void SetMaxSteps(long int numSteps)
virtual void VerifyStateVariables()
virtual ~AbstractCvodeSystem()
std::vector< double > & rGetTimes()
unsigned EstimateTimeSteps() const
double GetRelativeTolerance()
boost::shared_ptr< AbstractOdeSystemInformation > mpSystemInfo
OdeSolution Solve(realtype tStart, realtype tEnd, realtype maxDt, realtype tSamp)
N_Vector mLastSolutionState
std::string CheckShortMessageContains(std::string expected) const
double GetNextTime() const
double GetVectorComponent(const VECTOR &rVec, unsigned index)
void ForceUseOfNumericalJacobian(bool useNumericalJacobian=true)
AbstractCvodeSystem(unsigned numberOfStateVariables)
bool mHasAnalyticJacobian
void SetNumberOfTimeSteps(unsigned numTimeSteps)
const std::vector< std::string > & rGetParameterNames() const
void DeleteVector(VECTOR &rVec)