38 #include "CvodeAdaptor.hpp" 41 #include "TimeStepper.hpp" 49 #include <sundials/sundials_nvector.h> 50 #include <cvode/cvode_dense.h> 71 assert(pData !=
nullptr);
74 static std::vector<realtype> ydot_vec;
84 std::cerr <<
"CVODE RHS Exception: " << e.
GetMessage() << std::endl << std::flush;
112 int CvodeRootAdaptor(realtype t,
N_Vector y, realtype* pGOut,
void* pData)
114 assert(pData !=
nullptr);
125 std::cerr <<
"CVODE Root Exception: " << e.
GetMessage() << std::endl << std::flush;
167 void CvodeErrorHandler(
int errorCode,
const char *module,
const char *
function,
168 char *message,
void* pData)
170 std::stringstream err;
171 err <<
"CVODE Error " << errorCode <<
" in module " << module
172 <<
" function " <<
function <<
": " << message;
173 std::cerr <<
"*" << err.str() << std::endl << std::flush;
180 std::vector<double>& rInitialY,
185 assert(maxStep > 0.0);
188 N_Vector initial_values = N_VMake_Serial(rInitialY.size(), &(rInitialY[0]));
189 assert(NV_DATA_S(initial_values) == &(rInitialY[0]));
190 assert(!NV_OWN_DATA_S(initial_values));
200 for (
unsigned i=0; i<size; i++)
216 CVodeSetErrHandlerFn(
mpCvodeMem, CvodeErrorHandler,
nullptr);
220 #if CHASTE_SUNDIALS_VERSION >= 20400 227 #if CHASTE_SUNDIALS_VERSION >= 20400 228 CVodeInit(
mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values);
231 CVodeMalloc(
mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values,
238 #if CHASTE_SUNDIALS_VERSION >= 20400 239 CVodeRootInit(
mpCvodeMem, 1, CvodeRootAdaptor);
252 #if CHASTE_SUNDIALS_VERSION >= 20400 258 #if CHASTE_SUNDIALS_VERSION >= 20400 259 CVodeReInit(
mpCvodeMem, startTime, initial_values);
262 CVodeReInit(
mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values,
317 std::stringstream err;
318 char* p_flag_name = CVodeGetReturnFlagName(flag);
319 err << msg <<
": " << p_flag_name;
321 std::cerr << err.str() << std::endl << std::flush;
327 std::vector<double>& rYValues,
333 assert(endTime > startTime);
334 assert(timeSampling > 0.0);
339 EXCEPTION(
"(Solve with sampling) Stopping event is true for initial condition");
342 SetupCvode(pOdeSystem, rYValues, startTime, maxStep);
344 TimeStepper stepper(startTime, endTime, timeSampling);
345 N_Vector yout = N_VMake_Serial(rYValues.size(), &(rYValues[0]));
351 solutions.
rGetTimes().push_back(startTime);
367 CvodeError(ierr,
"CVODE failed to solve system");
372 if (ierr == CV_ROOT_RETURN)
395 std::vector<double>& rYValues,
400 assert(endTime > startTime);
405 EXCEPTION(
"(Solve) Stopping event is true for initial condition");
408 SetupCvode(pOdeSystem, rYValues, startTime, maxStep);
410 N_Vector yout = N_VMake_Serial(rYValues.size(), &(rYValues[0]));
413 int ierr = CVodeSetStopTime(
mpCvodeMem, endTime);
417 ierr = CVode(
mpCvodeMem, endTime, yout, &tend, CV_NORMAL);
422 CvodeError(ierr,
"CVODE failed to solve system");
424 if (ierr == CV_ROOT_RETURN)
430 assert(NV_DATA_S(yout) == &(rYValues[0]));
431 assert(!NV_OWN_DATA_S(yout));
453 #if CHASTE_SUNDIALS_VERSION >= 20400
474 for (
unsigned i=0; i<size; i++)
522 #endif // CHASTE_CVODE
AbstractOdeSystem * pSystem
bool mStoppingEventOccurred
void SetupCvode(AbstractOdeSystem *pOdeSystem, std::vector< double > &rInitialY, double startTime, double maxStep)
static bool WithinAnyTolerance(double number1, double number2, double relTol=DBL_EPSILON, double absTol=DBL_EPSILON, bool printError=false)
CvodeAdaptor(double relTol=1e-4, double absTol=1e-6)
std::vector< std::vector< double > > & rGetSolutions()
void CheckForStoppingEvents()
#define EXCEPTION(message)
N_Vector mLastSolutionState
virtual void EvaluateYDerivatives(double time, const std::vector< double > &rY, std::vector< double > &rDY)=0
double mLastInternalStepSize
void AdvanceOneTimeStep()
double GetAbsoluteTolerance()
unsigned GetNumberOfStateVariables() const
void CreateVectorIfEmpty(VECTOR &rVec, unsigned size)
void SetOdeSystemInformation(boost::shared_ptr< const AbstractOdeSystemInformation > pOdeSystemInfo)
void CvodeError(int flag, const char *msg)
void SetMinimalReset(bool minimalReset)
std::string GetMessage() const
virtual double CalculateRootFunction(double time, const std::vector< double > &rY)
void SetVectorComponent(VECTOR &rVec, unsigned index, double value)
std::vector< realtype > * pY
unsigned GetTotalTimeStepsTaken() const
boost::shared_ptr< const AbstractOdeSystemInformation > GetSystemInformation() const
std::vector< double > & rGetTimes()
unsigned EstimateTimeSteps() const
void SetTolerances(double relTol=1e-4, double absTol=1e-6)
#define CHASTE_CLASS_EXPORT(T)
void RecordStoppingPoint(double stopTime, N_Vector yEnd)
void SetForceReset(bool autoReset)
void CopyFromStdVector(const std::vector< double > &rSrc, VECTOR &rDest)
void CopyToStdVector(const VECTOR &rSrc, std::vector< double > &rDest)
double GetNextTime() const
double GetVectorComponent(const VECTOR &rVec, unsigned index)
void SetMaxSteps(long int numSteps)
double GetRelativeTolerance()
OdeSolution Solve(AbstractOdeSystem *pOdeSystem, std::vector< double > &rYValues, double startTime, double endTime, double maxStep, double timeSampling)
void SetNumberOfTimeSteps(unsigned numTimeSteps)
void DeleteVector(VECTOR &rVec)
unsigned GetVectorSize(const VECTOR &rVec)
virtual bool CalculateStoppingEvent(double time, const std::vector< double > &rY)