Chaste Commit::675f9facbe008c5eacb9006feaeb6423206579ea
CvodeAdaptor.cpp
1/*
2
3Copyright (c) 2005-2025, University of Oxford.
4All rights reserved.
5
6University of Oxford means the Chancellor, Masters and Scholars of the
7University of Oxford, having an administrative office at Wellington
8Square, Oxford OX1 2JD, UK.
9
10This file is part of Chaste.
11
12Redistribution and use in source and binary forms, with or without
13modification, are permitted provided that the following conditions are met:
14 * Redistributions of source code must retain the above copyright notice,
15 this list of conditions and the following disclaimer.
16 * Redistributions in binary form must reproduce the above copyright notice,
17 this list of conditions and the following disclaimer in the documentation
18 and/or other materials provided with the distribution.
19 * Neither the name of the University of Oxford nor the names of its
20 contributors may be used to endorse or promote products derived from this
21 software without specific prior written permission.
22
23THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33
34*/
35
36#ifdef CHASTE_CVODE
37
38#include "CvodeAdaptor.hpp"
39
40#include "Exception.hpp"
41#include "MathsCustomFunctions.hpp" // For tolerance check.
42#include "TimeStepper.hpp"
44
45#include <iostream>
46#include <sstream>
47
48// CVODE headers
49#include <sundials/sundials_nvector.h>
50
51#if CHASTE_SUNDIALS_VERSION >= 30000
52#if CHASTE_SUNDIALS_VERSION < 70000
53#include <cvode/cvode_direct.h> /* access to CVDls interface */
54#endif
55#include <sundials/sundials_types.h> /* defs. of realtype, sunindextype */
56#include <sunlinsol/sunlinsol_dense.h> /* access to dense SUNLinearSolver */
57#include <sunmatrix/sunmatrix_dense.h> /* access to dense SUNMatrix */
58#else
59#include <cvode/cvode_dense.h>
60#endif
61
62#if CHASTE_SUNDIALS_VERSION >= 60000
63#include "CvodeContextManager.hpp" // access to shared SUNContext object required by Sundials 6.0+
64#endif
65
82int CvodeRhsAdaptor(realtype t, N_Vector y, N_Vector ydot, void* pData)
83{
84 assert(pData != nullptr);
85 CvodeData* p_data = (CvodeData*)pData;
86 // Get y, ydot into std::vector<>s
87 static std::vector<realtype> ydot_vec;
88 CopyToStdVector(y, *p_data->pY);
89 CopyToStdVector(ydot, ydot_vec);
90 // Call our function
91 try
92 {
93 p_data->pSystem->EvaluateYDerivatives(t, *(p_data->pY), ydot_vec);
94 }
95 catch (const Exception& e)
96 {
97 std::cerr << "CVODE RHS Exception: " << e.GetMessage() << std::endl
98 << std::flush;
99 return -1;
100 }
101 // Copy derivative back
102 CopyFromStdVector(ydot_vec, ydot);
103 return 0;
104}
105
126int CvodeRootAdaptor(realtype t, N_Vector y, realtype* pGOut, void* pData)
127{
128 assert(pData != nullptr);
129 CvodeData* p_data = (CvodeData*)pData;
130 // Get y into a std::vector
131 CopyToStdVector(y, *p_data->pY);
132 // Call our function
133 try
134 {
135 *pGOut = p_data->pSystem->CalculateRootFunction(t, *p_data->pY);
136 }
137 catch (const Exception& e)
138 {
139 std::cerr << "CVODE Root Exception: " << e.GetMessage() << std::endl
140 << std::flush;
141 return -1;
142 }
143 return 0;
144}
145
146// /**
147// * Jacobian computation adaptor function.
148// *
149// * If solving an AbstractOdeSystemWithAnalyticJacobian, this function
150// * can be used to allow CVODE to compute the Jacobian analytically.
151// *
152// * Note to self: can test using pSystem->GetUseAnalyticJacobian().
153// */
154// int CvodeDenseJacobianAdaptor(long int numberOfStateVariables, DenseMat J,
155// realtype t, N_Vector y, N_Vector fy,
156// void* pData,
157// N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
158// {
159// AbstractOdeSystemWithAnalyticJacobian* pSystem
160// = (AbstractOdeSystemWithAnalyticJacobian*) pData;
161// // Get the current time step
162// double dt;
163// CVodeGetCurrentStep(CvodeMem, &dt);
164// // Get std::vector<> for y and double** for J
165// std::vector<realtype>& y_vec = *NV_DATA_STL(y);
166// double** ppJ = J->data; // organised column-wise: J_{i,j} = ppJ[j][i]
167// // Call our function
168// try
169// {
170// pSystem->AnalyticJacobian(y_vec, ppJ, t, dt);
171// }
172// catch (const Exception &e)
173// {
174// std::cerr << "CVODE J Exception: " << e.GetMessage() << std::endl << std::flush;
175// return -1;
176// }
177// // Update J (if needed)
178// return 0;
179// }
180
181#if CHASTE_SUNDIALS_VERSION >= 70000
182void CvodeErrorHandler(int errorCode, const char* module, const char* function,
183 const char* message, SUNErrCode errCode, void* pData, SUNContext sunContext)
184#else
185void CvodeErrorHandler(int errorCode, const char* module, const char* function,
186 char* message, void* pData)
187#endif
188{
189 std::stringstream err;
190 err << "CVODE Error " << errorCode << " in module " << module
191 << " function " << function << ": " << message;
192 std::cerr << "*" << err.str() << std::endl
193 << std::flush;
194 // Throwing the exception here causes termination on Maverick (g++ 4.4)
195 //EXCEPTION(err.str());
196}
197
199 std::vector<double>& rInitialY,
200 double startTime,
201 double maxStep)
202{
203 assert(rInitialY.size() == pOdeSystem->GetNumberOfStateVariables());
204 assert(maxStep > 0.0);
205
207#if CHASTE_SUNDIALS_VERSION >= 60000
208 N_Vector initial_values = N_VMake_Serial(rInitialY.size(), &(rInitialY[0]), CvodeContextManager::Instance()->GetSundialsContext());
209#else
210 N_Vector initial_values = N_VMake_Serial(rInitialY.size(), &(rInitialY[0]));
211#endif
212 assert(NV_DATA_S(initial_values) == &(rInitialY[0]));
213 assert(!NV_OWN_DATA_S(initial_values));
214 // std::cout << " Initial values: "; N_VPrint_Serial(initial_values);
215 // std::cout << " Rtol: " << mRelTol << ", Atol: " << mAbsTol << std::endl;
216 // std::cout << " Start: " << startTime << " max dt=" << maxStep << std::endl << std::flush;
217
218 // Find out if we need to (re-)initialise
220 if (!reinit && !mForceMinimalReset)
221 {
222 const unsigned size = GetVectorSize(rInitialY);
223 for (unsigned i = 0; i < size; i++)
224 {
226 {
227 reinit = true;
228 break;
229 }
230 }
231 }
232
233 if (!mpCvodeMem) // First run of this solver, set up CVODE memory
234 {
235 // Set up CVODE's memory.
236#if CHASTE_SUNDIALS_VERSION >= 60000
237 mpCvodeMem = CVodeCreate(CV_BDF, CvodeContextManager::Instance()->GetSundialsContext());
238#elif CHASTE_SUNDIALS_VERSION >= 40000
239 // v4.0.0 release notes: instead of specifying the nonlinear iteration type when creating the CVODE(S) memory structure,
240 // CVODE(S) uses the SUNNONLINSOL_NEWTON module implementation of a Newton iteration by default.
241 mpCvodeMem = CVodeCreate(CV_BDF);
242#else
243 mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
244#endif
245
246 if (mpCvodeMem == nullptr)
247 EXCEPTION("Failed to SetupCvode CVODE"); // LCOV_EXCL_LINE
248 // Set error handler
249#if CHASTE_SUNDIALS_VERSION >= 70000
250 SUNContext_PushErrHandler(CvodeContextManager::Instance()->GetSundialsContext(), CvodeErrorHandler, nullptr);
251#else
252 CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, nullptr);
253#endif
254 // Set the user data
255 mData.pSystem = pOdeSystem;
256 mData.pY = &rInitialY;
257#if CHASTE_SUNDIALS_VERSION >= 20400
258 CVodeSetUserData(mpCvodeMem, (void*)(&mData));
259#else
260 CVodeSetFdata(mpCvodeMem, (void*)(&mData));
261#endif
262
263// Setup CVODE
264#if CHASTE_SUNDIALS_VERSION >= 20400
265 CVodeInit(mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values);
266 CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
267#else
268 CVodeMalloc(mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values,
269 CV_SS, mRelTol, &mAbsTol);
270#endif
271
272 // Set the rootfinder function if wanted
273 if (mCheckForRoots)
274 {
275#if CHASTE_SUNDIALS_VERSION >= 20400
276 CVodeRootInit(mpCvodeMem, 1, CvodeRootAdaptor);
277#else
278 CVodeRootInit(mpCvodeMem, 1, CvodeRootAdaptor, (void*)(&mData));
279#endif
280 }
281
282 /* Create dense SUNMatrix for use in linear solves */
283#if CHASTE_SUNDIALS_VERSION >= 60000
284 mpSundialsDenseMatrix = SUNDenseMatrix(rInitialY.size(), rInitialY.size(), CvodeContextManager::Instance()->GetSundialsContext());
285#elif CHASTE_SUNDIALS_VERSION >= 30000
286 mpSundialsDenseMatrix = SUNDenseMatrix(rInitialY.size(), rInitialY.size());
287#endif
288
289#if CHASTE_SUNDIALS_VERSION >= 60000
290 /* Create dense SUNLinearSolver object for use by CVode */
291 mpSundialsLinearSolver = SUNLinSol_Dense(initial_values, mpSundialsDenseMatrix, CvodeContextManager::Instance()->GetSundialsContext());
292
293 /* Call CVodeSetLinearSolver to attach the matrix and linear solver to CVode */
294 CVodeSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
295#elif CHASTE_SUNDIALS_VERSION >= 40000
296 /* Create dense SUNLinearSolver object for use by CVode */
297 mpSundialsLinearSolver = SUNLinSol_Dense(initial_values, mpSundialsDenseMatrix);
298
299 /* Call CVodeSetLinearSolver to attach the matrix and linear solver to CVode */
300 CVodeSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
301#elif CHASTE_SUNDIALS_VERSION >= 30000
302 /* Create dense SUNLinearSolver object for use by CVode */
303 mpSundialsLinearSolver = SUNDenseLinearSolver(initial_values, mpSundialsDenseMatrix);
304
305 /* Call CVDlsSetLinearSolver to attach the matrix and linear solver to CVode */
306 CVDlsSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
307#else
308 // CVODE < v3.0.0
309 // Attach a linear solver for Newton iteration
310 CVDense(mpCvodeMem, rInitialY.size());
311#endif
312 }
313 else if (reinit) // Could be new ODE system, or new Y values
314 {
315 // Set the user data
316 mData.pSystem = pOdeSystem; // stays the same on a re-initialize
317 mData.pY = &rInitialY; // changes on a re-initialize
318#if CHASTE_SUNDIALS_VERSION >= 20400
319 CVodeSetUserData(mpCvodeMem, (void*)(&mData));
320#else
321 CVodeSetFdata(mpCvodeMem, (void*)(&mData));
322#endif
323
324#if CHASTE_SUNDIALS_VERSION >= 20400
325 CVodeReInit(mpCvodeMem, startTime, initial_values);
326 CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
327#else
328 CVodeReInit(mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values,
329 CV_SS, mRelTol, &mAbsTol);
330#endif
331
332#if CHASTE_SUNDIALS_VERSION >= 30000
333 if (mpSundialsLinearSolver)
334 {
335 /* Free the linear solver memory */
336 SUNLinSolFree(mpSundialsLinearSolver);
337 }
338 if (mpSundialsDenseMatrix)
339 {
340 /* Free the matrix memory */
341 SUNMatDestroy(mpSundialsDenseMatrix);
342 }
343
344 /* Create dense matrix of type SUNDenseMatrix for use in linear solves */
345#if CHASTE_SUNDIALS_VERSION >= 60000
346 mpSundialsDenseMatrix = SUNDenseMatrix(rInitialY.size(), rInitialY.size(), CvodeContextManager::Instance()->GetSundialsContext());
347#else
348 mpSundialsDenseMatrix = SUNDenseMatrix(rInitialY.size(), rInitialY.size());
349#endif
350
351#if CHASTE_SUNDIALS_VERSION >= 60000
352 /* Create dense SUNLinSol_Dense object for use by CVode */
353 mpSundialsLinearSolver = SUNLinSol_Dense(initial_values, mpSundialsDenseMatrix, CvodeContextManager::Instance()->GetSundialsContext());
354
355 /* Call CVodeSetLinearSolver to attach the matrix and linear solver to CVode */
356 CVodeSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
357#elif CHASTE_SUNDIALS_VERSION >= 40000
358 /* Create dense SUNLinSol_Dense object for use by CVode */
359 mpSundialsLinearSolver = SUNLinSol_Dense(initial_values, mpSundialsDenseMatrix);
360
361 /* Call CVodeSetLinearSolver to attach the matrix and linear solver to CVode */
362 CVodeSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
363#else
364 /* Create dense SUNLinearSolver object for use by CVode */
365 mpSundialsLinearSolver = SUNDenseLinearSolver(initial_values, mpSundialsDenseMatrix);
366
367 /* Call CVDlsSetLinearSolver to attach the matrix and linear solver to CVode */
368 CVDlsSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
369#endif
370#else
371 // Attach a linear solver for Newton iteration
372 CVDense(mpCvodeMem, rInitialY.size());
373#endif
374 }
375
376 CVodeSetMaxStep(mpCvodeMem, maxStep);
377 // Change max steps if wanted
378 if (mMaxSteps > 0)
379 {
380 CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
381 CVodeSetMaxErrTestFails(mpCvodeMem, 15);
382 }
383 DeleteVector(initial_values);
384}
385
387{
388 if (mpCvodeMem)
389 {
390 CVodeFree(&mpCvodeMem);
391 }
392 mpCvodeMem = nullptr;
393
394#if CHASTE_SUNDIALS_VERSION >= 30000
395 if (mpSundialsLinearSolver)
396 {
397 /* Free the linear solver memory */
398 SUNLinSolFree(mpSundialsLinearSolver);
399 }
400 mpSundialsLinearSolver = nullptr;
401
402 if (mpSundialsDenseMatrix)
403 {
404 /* Free the matrix memory */
405 SUNMatDestroy(mpSundialsDenseMatrix);
406 }
407 mpSundialsDenseMatrix = nullptr;
408#endif
409}
410
411void CvodeAdaptor::SetForceReset(bool autoReset)
412{
413 mForceReset = autoReset;
414 if (mForceReset)
415 {
416 SetMinimalReset(false);
417 ResetSolver();
418 }
419}
420
421void CvodeAdaptor::SetMinimalReset(bool minimalReset)
422{
423 mForceMinimalReset = minimalReset;
425 {
426 SetForceReset(false);
427 }
428}
429
431{
432 // Doing this makes the SetupCvode think it needs to reset things.
434}
435
436void CvodeAdaptor::CvodeError(int flag, const char* msg)
437{
438 std::stringstream err;
439 char* p_flag_name = CVodeGetReturnFlagName(flag);
440 err << msg << ": " << p_flag_name;
441 free(p_flag_name);
442 std::cerr << err.str() << std::endl
443 << std::flush;
444 EXCEPTION(err.str());
445}
446
448 std::vector<double>& rYValues,
449 double startTime,
450 double endTime,
451 double maxStep,
452 double timeSampling)
453{
454 assert(endTime > startTime);
455 assert(timeSampling > 0.0);
456
458 if (mCheckForRoots && pOdeSystem->CalculateStoppingEvent(startTime, rYValues) == true)
459 {
460 EXCEPTION("(Solve with sampling) Stopping event is true for initial condition");
461 }
462
463 SetupCvode(pOdeSystem, rYValues, startTime, maxStep);
464
465 TimeStepper stepper(startTime, endTime, timeSampling);
466#if CHASTE_SUNDIALS_VERSION >= 60000
467 N_Vector yout = N_VMake_Serial(rYValues.size(), &(rYValues[0]), CvodeContextManager::Instance()->GetSundialsContext());
468#else
469 N_Vector yout = N_VMake_Serial(rYValues.size(), &(rYValues[0]));
470#endif
471
472 // Set up ODE solution
473 OdeSolution solutions;
474 solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
475 solutions.rGetSolutions().push_back(rYValues);
476 solutions.rGetTimes().push_back(startTime);
477 solutions.SetOdeSystemInformation(pOdeSystem->GetSystemInformation());
478
479 // Main time sampling loop
480 while (!stepper.IsTimeAtEnd() && !mStoppingEventOccurred)
481 {
482 // This should stop CVODE going past the end of where we wanted and interpolating back.
483 int ierr = CVodeSetStopTime(mpCvodeMem, stepper.GetNextTime());
484 assert(ierr == CV_SUCCESS);
485 UNUSED_OPT(ierr); // avoid unused var warning
486
487 double tend;
488 ierr = CVode(mpCvodeMem, stepper.GetNextTime(), yout, &tend, CV_NORMAL);
489 if (ierr < 0)
490 {
492 DeleteVector(yout);
493 CvodeError(ierr, "CVODE failed to solve system");
494 }
495 // Store solution
496 solutions.rGetSolutions().push_back(rYValues);
497 solutions.rGetTimes().push_back(tend);
498 if (ierr == CV_ROOT_RETURN)
499 {
500 // Stopping event occurred
502 mStoppingTime = tend;
503 }
504 mLastSolutionTime = tend;
505 stepper.AdvanceOneTimeStep();
506 }
507
508 // stepper.EstimateTimeSteps may have been an overestimate...
509 solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
510
511 int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
512 assert(ierr == CV_SUCCESS);
513 UNUSED_OPT(ierr); // avoid unused var warning
515 DeleteVector(yout);
516
517 return solutions;
518}
519
521 std::vector<double>& rYValues,
522 double startTime,
523 double endTime,
524 double maxStep)
525{
526 assert(endTime > startTime);
527
529 if (mCheckForRoots && pOdeSystem->CalculateStoppingEvent(startTime, rYValues) == true)
530 {
531 EXCEPTION("(Solve) Stopping event is true for initial condition");
532 }
533
534 SetupCvode(pOdeSystem, rYValues, startTime, maxStep);
535
536#if CHASTE_SUNDIALS_VERSION >= 60000
537 N_Vector yout = N_VMake_Serial(rYValues.size(), &(rYValues[0]), CvodeContextManager::Instance()->GetSundialsContext());
538#else
539 N_Vector yout = N_VMake_Serial(rYValues.size(), &(rYValues[0]));
540#endif
541
542 // This should stop CVODE going past the end of where we wanted and interpolating back.
543 int ierr = CVodeSetStopTime(mpCvodeMem, endTime);
544 assert(ierr == CV_SUCCESS);
545 UNUSED_OPT(ierr); // avoid unused var warning
546
547 double tend;
548 ierr = CVode(mpCvodeMem, endTime, yout, &tend, CV_NORMAL);
549 if (ierr < 0)
550 {
552 DeleteVector(yout);
553 CvodeError(ierr, "CVODE failed to solve system");
554 }
555 if (ierr == CV_ROOT_RETURN)
556 {
557 // Stopping event occurred
559 mStoppingTime = tend;
560 }
561 assert(NV_DATA_S(yout) == &(rYValues[0]));
562 assert(!NV_OWN_DATA_S(yout));
563
564 // long int steps;
565 // CVodeGetNumSteps(mpCvodeMem, &steps);
566 // std::cout << " Solved to " << endTime << " in " << steps << " steps.\n";
567
568 ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
569 assert(ierr == CV_SUCCESS);
570 UNUSED_OPT(ierr); // avoid unused var warning
571 RecordStoppingPoint(tend, yout);
572 DeleteVector(yout);
573}
574
575CvodeAdaptor::CvodeAdaptor(double relTol, double absTol)
577 mpCvodeMem(nullptr),
578 mRelTol(relTol),
579 mAbsTol(absTol),
580 mLastInternalStepSize(-0.0),
581 mMaxSteps(0),
582 mCheckForRoots(false),
583 mLastSolutionState(nullptr),
584 mLastSolutionTime(0.0),
585#if CHASTE_SUNDIALS_VERSION >= 20400
586 mForceReset(false),
587#else
588 mForceReset(true),
589#endif
590 mForceMinimalReset(false)
591#if CHASTE_SUNDIALS_VERSION >= 30000
592 ,
593 mpSundialsDenseMatrix(nullptr),
594 mpSundialsLinearSolver(nullptr)
595#endif
596{
597}
598
604
606{
607 if (!mForceReset)
608 {
609 const unsigned size = GetVectorSize(yEnd);
611 for (unsigned i = 0; i < size; i++)
612 {
614 }
615 mLastSolutionTime = stopTime;
616 }
617}
618
619void CvodeAdaptor::SetTolerances(double relTol, double absTol)
620{
621 mRelTol = relTol;
622 mAbsTol = absTol;
623}
624
626{
627 return mRelTol;
628}
629
631{
632 return mAbsTol;
633}
634
639
644
645void CvodeAdaptor::SetMaxSteps(long int numSteps)
646{
647 mMaxSteps = numSteps;
648}
649
651{
652 return mMaxSteps;
653}
654
655// Serialization for Boost >= 1.36
658
659#endif // CHASTE_CVODE
#define EXCEPTION(message)
#define UNUSED_OPT(var)
#define CHASTE_CLASS_EXPORT(T)
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)
unsigned GetVectorSize(const VECTOR &rVec)
void CopyFromStdVector(const std::vector< double > &rSrc, VECTOR &rDest)
void CopyToStdVector(const VECTOR &rSrc, std::vector< double > &rDest)
virtual double CalculateRootFunction(double time, const std::vector< double > &rY)
virtual void EvaluateYDerivatives(double time, const std::vector< double > &rY, std::vector< double > &rDY)=0
virtual bool CalculateStoppingEvent(double time, const std::vector< double > &rY)
boost::shared_ptr< const AbstractOdeSystemInformation > GetSystemInformation() const
static bool WithinAnyTolerance(double number1, double number2, double relTol=DBL_EPSILON, double absTol=DBL_EPSILON, bool printError=false)
double GetAbsoluteTolerance()
long int GetMaxSteps()
void CvodeError(int flag, const char *msg)
long int mMaxSteps
void SetMinimalReset(bool minimalReset)
N_Vector mLastSolutionState
void SetMaxSteps(long int numSteps)
void FreeCvodeMemory()
void SetForceReset(bool autoReset)
void CheckForStoppingEvents()
CvodeData mData
OdeSolution Solve(AbstractOdeSystem *pOdeSystem, std::vector< double > &rYValues, double startTime, double endTime, double maxStep, double timeSampling)
double GetLastStepSize()
double mLastInternalStepSize
void RecordStoppingPoint(double stopTime, N_Vector yEnd)
void SetupCvode(AbstractOdeSystem *pOdeSystem, std::vector< double > &rInitialY, double startTime, double maxStep)
double mLastSolutionTime
double GetRelativeTolerance()
void SetTolerances(double relTol=1e-4, double absTol=1e-6)
CvodeAdaptor(double relTol=1e-4, double absTol=1e-6)
void SetNumberOfTimeSteps(unsigned numTimeSteps)
std::vector< std::vector< double > > & rGetSolutions()
std::vector< double > & rGetTimes()
void SetOdeSystemInformation(boost::shared_ptr< const AbstractOdeSystemInformation > pOdeSystemInfo)
bool IsTimeAtEnd() const
unsigned GetTotalTimeStepsTaken() const
void AdvanceOneTimeStep()
double GetNextTime() const
unsigned EstimateTimeSteps() const
AbstractOdeSystem * pSystem
std::vector< realtype > * pY