Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
AbstractCvodeSystem.cpp
1/*
2
3Copyright (c) 2005-2024, 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 <cassert>
39#include <sstream>
40
41#include "AbstractCvodeSystem.hpp"
42#include "CvodeAdaptor.hpp" // For CvodeErrorHandler
43#include "Exception.hpp"
44#include "MathsCustomFunctions.hpp" // For tolerance comparison
45#include "TimeStepper.hpp"
47
48// CVODE headers
49#include <cvode/cvode.h>
50#include <sundials/sundials_nvector.h>
51
52#if CHASTE_SUNDIALS_VERSION >= 30000
53#include <cvode/cvode_direct.h> /* access to CVDls interface */
54#include <sundials/sundials_types.h> /* defs. of realtype, sunindextype */
55#include <sunlinsol/sunlinsol_dense.h> /* access to dense SUNLinearSolver */
56#include <sunmatrix/sunmatrix_dense.h> /* access to dense SUNMatrix */
57#else
58#include <cvode/cvode_dense.h>
59#endif
60
61#if CHASTE_SUNDIALS_VERSION >= 60000
62#include "CvodeContextManager.hpp" // access to shared SUNContext object required by Sundials 6.0+
63#endif
64
65//#include "Debug.hpp"
66//void DebugSteps(void* pCvodeMem, AbstractCvodeSystem* pSys)
67//{
68// long int num_jac_evals, nniters, num_steps;
69// CVDenseGetNumJacEvals(pCvodeMem, &num_jac_evals);
71// CVodeGetNumNonlinSolvIters(pCvodeMem, &nniters);
72// CVodeGetNumSteps(pCvodeMem, &num_steps);
73// double num_newton_iters = nniters;
74// PRINT_3_VARIABLES(pSys->GetSystemName(), num_newton_iters/num_steps, num_jac_evals/num_newton_iters);
75//}
85int AbstractCvodeSystemRhsAdaptor(realtype t, N_Vector y, N_Vector ydot, void* pData)
86{
87 assert(pData != nullptr);
88 AbstractCvodeSystem* p_ode_system = (AbstractCvodeSystem*)pData;
89 try
90 {
91 p_ode_system->EvaluateYDerivatives(t, y, ydot);
92 }
93 catch (const Exception& e)
94 {
95#if CHASTE_SUNDIALS_VERSION <= 20300
96 // Really old CVODE used to solve past the requested time points and could trigger this exception unnecessarily...
97 if (e.CheckShortMessageContains("is outside the times stored in the data clamp") == "")
98 {
99 return 1; // This may be a recoverable error!
100 }
101#endif
102
103 std::cerr << "CVODE RHS Exception: " << e.GetMessage()
104 << std::endl
105 << std::flush;
106 return -1;
107 }
108
109 // Something like this might help CVODE when things are a bit unstable...
110 // We tried this again in Jan 2023 but now with VerifyStateVariables() returning a recoverable error when values look off. However this didn't improve the situation.
111 // See: https://github.com/Chaste/Chaste/issues/46
112 // try
113 // {
114 // p_ode_system->VerifyStateVariables();
115 // }
116 // catch (const Exception &e)
117 // {
118 // std::cout << "t = " << t << ":\t" << e.GetMessage() << std::endl << std::flush;
119 // return 1; // A positive return flag to CVODE tells it there's been an error but it might be recoverable.
120 // }
121
122 return 0;
123}
124
125/*
126 * Absolute chaos here with four different possible interfaces to the jacobian.
127 */
128#if CHASTE_SUNDIALS_VERSION >= 30000
129// Sundials 3.0 - has taken away the argument N at the top...
130int AbstractCvodeSystemJacAdaptor(realtype t, N_Vector y, N_Vector ydot, CHASTE_CVODE_DENSE_MATRIX jacobian,
131#elif CHASTE_SUNDIALS_VERSION >= 20500
132// Sundials 2.5
133int AbstractCvodeSystemJacAdaptor(long int N, realtype t, N_Vector y, N_Vector ydot, CHASTE_CVODE_DENSE_MATRIX jacobian,
134#elif CHASTE_SUNDIALS_VERSION >= 20400
135// Sundials 2.4
136int AbstractCvodeSystemJacAdaptor(int N, realtype t, N_Vector y, N_Vector ydot, DlsMat jacobian,
137#else
138// Sundials 2.3 and below (not sure how far below, but this is 2006 so old enough).
139int AbstractCvodeSystemJacAdaptor(long int N, DenseMat jacobian, realtype t, N_Vector y, N_Vector ydot,
140#endif
141 void* pData, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
142{
143 assert(pData != nullptr);
144 AbstractCvodeSystem* p_ode_system = (AbstractCvodeSystem*)pData;
145 try
146 {
147 p_ode_system->EvaluateAnalyticJacobian(t, y, ydot, jacobian, tmp1, tmp2, tmp3);
148 }
149 catch (const Exception& e)
150 {
151 std::cerr << "CVODE Jacobian Exception: " << e.GetMessage() << std::endl
152 << std::flush;
153 return -1;
154 }
155 return 0;
156}
157
158AbstractCvodeSystem::AbstractCvodeSystem(unsigned numberOfStateVariables)
159 : AbstractParameterisedSystem<N_Vector>(numberOfStateVariables),
160 mLastSolutionState(nullptr),
161 mLastSolutionTime(0.0),
162#if CHASTE_SUNDIALS_VERSION >= 20400
163 mForceReset(false),
164#else
165 // Old Sundials don't seem to 'go back' when something has changed
166 // properly, and give more inaccurate answers.
167 mForceReset(true),
168#endif
169 mForceMinimalReset(false),
170#if CHASTE_SUNDIALS_VERSION >= 30000
171 mpSundialsDenseMatrix(nullptr),
172 mpSundialsLinearSolver(nullptr),
173#endif
174 mHasAnalyticJacobian(false),
175 mUseAnalyticJacobian(false),
176 mpCvodeMem(nullptr),
177 mMaxSteps(0),
178 mLastInternalStepSize(0)
179{
180 SetTolerances(); // Set the tolerances to the defaults.
181}
182
184{
188#if CHASTE_SUNDIALS_VERSION >= 60000
189 mParameters = N_VNew_Serial(rGetParameterNames().size(), CvodeContextManager::Instance()->GetSundialsContext());
190#else
191 mParameters = N_VNew_Serial(rGetParameterNames().size());
192#endif
193 for (int i = 0; i < NV_LENGTH_S(mParameters); i++)
194 {
195 NV_Ith_S(mParameters, i) = 0.0;
196 }
197}
198
206
207//
208//double AbstractCvodeSystem::CalculateRootFunction(double time, const std::vector<double>& rY)
209//{
210// bool stop = CalculateStoppingEvent(time, rY);
211// return stop ? 0.0 : 1.0;
212//}
213
215 realtype tEnd,
216 realtype maxDt,
217 realtype tSamp)
218{
219 assert(tEnd >= tStart);
220 assert(tSamp > 0.0);
221
222 SetupCvode(mStateVariables, tStart, maxDt);
223
224 TimeStepper stepper(tStart, tEnd, tSamp);
225
226 // Set up ODE solution
227 OdeSolution solutions;
228 solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
229 solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
230 solutions.rGetTimes().push_back(tStart);
232
233 // Main time sampling loop
234 while (!stepper.IsTimeAtEnd())
235 {
236 // This should stop CVODE going past the end of where we wanted and interpolating back.
237 int ierr = CVodeSetStopTime(mpCvodeMem, stepper.GetNextTime());
238 assert(ierr == CV_SUCCESS);
239 UNUSED_OPT(ierr); // avoid unused var warning
240
241 // // This parameter governs how many times we allow a recoverable right hand side failure
242 // int ierr = CVodeSetMaxConvFails(mpCvodeMem, 1000);
243 // assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr); // avoid unused var warning
244
245 double cvode_stopped_at = stepper.GetTime();
246 ierr = CVode(mpCvodeMem, stepper.GetNextTime(), mStateVariables,
247 &cvode_stopped_at, CV_NORMAL);
248 if (ierr < 0)
249 {
250 // DebugSteps(mpCvodeMem, this);
251 CvodeError(ierr, "CVODE failed to solve system", cvode_stopped_at, stepper.GetTime(), stepper.GetNextTime());
252 }
253 // Not root finding, so should have reached requested time
254 assert(fabs(cvode_stopped_at - stepper.GetNextTime()) < DBL_EPSILON);
255#ifndef NDEBUG
257#endif
258 // Store solution
259 solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
260 solutions.rGetTimes().push_back(cvode_stopped_at);
261 stepper.AdvanceOneTimeStep();
262 }
263
264 // stepper.EstimateTimeSteps may have been an overestimate...
265 solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
266
267 int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
268 assert(ierr == CV_SUCCESS);
269 UNUSED_OPT(ierr); // avoid unused var warning
270
272
273 return solutions;
274}
275
276void AbstractCvodeSystem::Solve(realtype tStart,
277 realtype tEnd,
278 realtype maxDt)
279{
280 assert(tEnd >= tStart);
281
282 SetupCvode(mStateVariables, tStart, maxDt);
283
284 // This should stop CVODE going past the end of where we wanted and interpolating back.
285 int ierr = CVodeSetStopTime(mpCvodeMem, tEnd);
286 assert(ierr == CV_SUCCESS);
287 UNUSED_OPT(ierr); // avoid unused var warning
288
289 double cvode_stopped_at = tStart;
290 ierr = CVode(mpCvodeMem, tEnd, mStateVariables, &cvode_stopped_at, CV_NORMAL);
291 if (ierr < 0)
292 {
293 // DebugSteps(mpCvodeMem, this);
294 CvodeError(ierr, "CVODE failed to solve system", cvode_stopped_at, tStart, tEnd);
295 }
296 // Not root finding, so should have reached requested time
297 assert(fabs(cvode_stopped_at - tEnd) < DBL_EPSILON);
298
299 ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
300 assert(ierr == CV_SUCCESS);
301 UNUSED_OPT(ierr); // avoid unused var warning
302
303 RecordStoppingPoint(cvode_stopped_at);
304
305//
306// long int nst, nfe, nsetups, nje, nfeLS, nni, ncfn, netf, nge;
307//
308//
309// CVodeGetNumSteps(mpCvodeMem, &nst);
310// CVodeGetNumRhsEvals(mpCvodeMem, &nfe);
311// CVodeGetNumLinSolvSetups(mpCvodeMem, &nsetups);
312// CVodeGetNumErrTestFails(mpCvodeMem, &netf);
313// CVodeGetNumNonlinSolvIters(mpCvodeMem, &nni);
314// CVodeGetNumNonlinSolvConvFails(mpCvodeMem, &ncfn);
315// CVDlsGetNumJacEvals(mpCvodeMem, &nje);
316// CVDlsGetNumRhsEvals(mpCvodeMem, &nfeLS);
317// CVodeGetNumGEvals(mpCvodeMem, &nge);
318//
319// printf("\nFinal Statistics:\n");
320// printf("nst = %-6ld nfe = %-6ld nsetups = %-6ld nfeLS = %-6ld nje = %ld\n",
321// nst, nfe, nsetups, nfeLS, nje);
322// printf("nni = %-6ld ncfn = %-6ld netf = %-6ld nge = %ld\n \n",
323// nni, ncfn, netf, nge);
324// std::cout << std::flush;
325#ifndef NDEBUG
327#endif
328}
329
330void AbstractCvodeSystem::SetMaxSteps(long int numSteps)
331{
332 mMaxSteps = numSteps;
333}
334
336{
337 return mMaxSteps;
338}
339
340void AbstractCvodeSystem::SetTolerances(double relTol, double absTol)
341{
342 mRelTol = relTol;
343 mAbsTol = absTol;
344 ResetSolver();
345}
346
351
356
361
363{
364 mForceReset = autoReset;
365 if (mForceReset)
366 {
367 ResetSolver();
368 }
369}
370
375
380
382{
383 mForceMinimalReset = minimalReset;
385 {
386 SetForceReset(false);
387 }
388}
389
394
396 realtype tStart,
397 realtype maxDt)
398{
399 assert((unsigned)NV_LENGTH_S(initialConditions) == GetNumberOfStateVariables());
400 assert(maxDt >= 0.0);
401
402 // Find out if we need to (re-)initialise
403 //std::cout << "!mpCvodeMem = " << !mpCvodeMem << ", mForceReset = " << mForceReset << ", !mLastSolutionState = " << !mLastSolutionState << ", comp doubles = " << !CompareDoubles::WithinAnyTolerance(tStart, mLastSolutionTime) << "\n";
405 if (!reinit && !mForceMinimalReset)
406 {
407 const unsigned size = GetNumberOfStateVariables();
408 for (unsigned i = 0; i < size; i++)
409 {
411 {
412 reinit = true;
413 break;
414 }
415 }
416 }
417
418 if (!mpCvodeMem)
419 {
420 //std::cout << "New CVODE solver\n";
421#if CHASTE_SUNDIALS_VERSION >= 60000
422 mpCvodeMem = CVodeCreate(CV_BDF, CvodeContextManager::Instance()->GetSundialsContext());
423#elif CHASTE_SUNDIALS_VERSION >= 40000
424 // v4.0.0 release notes: instead of specifying the nonlinear iteration type when creating the CVODE(S) memory structure,
425 // CVODE(S) uses the SUNNONLINSOL_NEWTON module implementation of a Newton iteration by default.
426 mpCvodeMem = CVodeCreate(CV_BDF);
427#else
428 mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
429#endif
430 if (mpCvodeMem == nullptr)
431 EXCEPTION("Failed to SetupCvode CVODE"); // LCOV_EXCL_LINE
432
433 // Set error handler
434 CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, nullptr);
435// Set the user data
436#if CHASTE_SUNDIALS_VERSION >= 20400
437 CVodeSetUserData(mpCvodeMem, (void*)(this));
438#else
439 CVodeSetFdata(mpCvodeMem, (void*)(this));
440#endif
441// Setup CVODE
442#if CHASTE_SUNDIALS_VERSION >= 20400
443 CVodeInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions);
444 CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
445#else
446 CVodeMalloc(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
447 CV_SS, mRelTol, &mAbsTol);
448#endif
449
450#if CHASTE_SUNDIALS_VERSION >= 60000
451 /* Create dense matrix SUNDenseMatrix for use in linear solves */
452 mpSundialsDenseMatrix = SUNDenseMatrix(NV_LENGTH_S(initialConditions), NV_LENGTH_S(initialConditions), CvodeContextManager::Instance()->GetSundialsContext());
453#elif CHASTE_SUNDIALS_VERSION >= 30000
454 /* Create dense matrix SUNDenseMatrix for use in linear solves */
455 mpSundialsDenseMatrix = SUNDenseMatrix(NV_LENGTH_S(initialConditions), NV_LENGTH_S(initialConditions));
456#endif
457
458#if CHASTE_SUNDIALS_VERSION >= 60000
459 /* Create dense SUNLinSol_Dense object for use by CVode */
460 mpSundialsLinearSolver = SUNLinSol_Dense(initialConditions, mpSundialsDenseMatrix, CvodeContextManager::Instance()->GetSundialsContext());
461
462 /* Call CVodeSetLinearSolver to attach the matrix and linear solver to CVode */
463 CVodeSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
464#elif CHASTE_SUNDIALS_VERSION >= 40000
465 /* Create dense SUNLinSol_Dense object for use by CVode */
466 mpSundialsLinearSolver = SUNLinSol_Dense(initialConditions, mpSundialsDenseMatrix);
467
468 /* Call CVodeSetLinearSolver to attach the matrix and linear solver to CVode */
469 CVodeSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
470#elif CHASTE_SUNDIALS_VERSION >= 30000
471 /* Create dense SUNDenseLinearSolver object for use by CVode */
472 mpSundialsLinearSolver = SUNDenseLinearSolver(initialConditions, mpSundialsDenseMatrix);
473
474 /* Call CVDlsSetLinearSolver to attach the matrix and linear solver to CVode */
475 CVDlsSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
476#else
477 // CVODE < v3.0.0
478 // Attach a linear solver for Newton iteration
479 CVDense(mpCvodeMem, NV_LENGTH_S(initialConditions));
480#endif
481
483 {
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);
490#else
491 CVDenseSetJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor, (void*)(this));
492#endif
493 }
494 }
495 else if (reinit)
496 {
497//std::cout << "Resetting CVODE solver\n";
498#if CHASTE_SUNDIALS_VERSION >= 20400
499 CVodeReInit(mpCvodeMem, tStart, initialConditions);
500 //CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol); - "all solver inputs remain in effect" so we don't need this.
501#else
502 CVodeReInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
503 CV_SS, mRelTol, &mAbsTol);
504#endif
505 }
506
507 // Set max dt and change max steps if wanted
508 if (maxDt > 0)
509 {
510 CVodeSetMaxStep(mpCvodeMem, maxDt);
511 }
512
513 if (mMaxSteps > 0)
514 {
515 CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
516 CVodeSetMaxErrTestFails(mpCvodeMem, 15);
517 }
518}
519
521{
522 // DebugSteps(mpCvodeMem, this);
523
524 // If we're forcing a reset then we don't record the stopping time
525 // as a result it won't match and we will force a reset in SetupCvode() on
526 // the next solve call.
527 if (mForceReset)
528 return;
529
530 // Otherwise we will store the state variables and time for comparison on the
531 // next solve call, to work out whether we need to reset.
532 const unsigned size = GetNumberOfStateVariables();
534 for (unsigned i = 0; i < size; i++)
535 {
537 }
538 mLastSolutionTime = stopTime;
539}
540
542{
543 if (mpCvodeMem)
544 {
545 CVodeFree(&mpCvodeMem);
546 }
547 mpCvodeMem = nullptr;
548
549#if CHASTE_SUNDIALS_VERSION >= 30000
550 if (mpSundialsLinearSolver)
551 {
552 /* Free the linear solver memory */
553 SUNLinSolFree(mpSundialsLinearSolver);
554 }
555 mpSundialsLinearSolver = nullptr;
556
557 if (mpSundialsDenseMatrix)
558 {
559 /* Free the matrix memory */
560 SUNMatDestroy(mpSundialsDenseMatrix);
561 }
562 mpSundialsDenseMatrix = nullptr;
563#endif
564}
565
566void AbstractCvodeSystem::CvodeError(int flag, const char* msg,
567 const double& rTime, const double& rStartTime, const double& rEndTime)
568{
569 std::stringstream err;
570 char* p_flag_name = CVodeGetReturnFlagName(flag);
571 err << msg << ": " << p_flag_name;
572 free(p_flag_name);
573 if (flag == CV_LSETUP_FAIL)
574 {
575#if CHASTE_SUNDIALS_VERSION >= 20500
576 long int ls_flag;
577#else
578 int ls_flag;
579#endif
580 char* p_ls_flag_name;
581
582#if CHASTE_SUNDIALS_VERSION >= 40000
583 CVodeGetLastLinFlag(mpCvodeMem, &ls_flag);
584 p_ls_flag_name = CVodeGetLinReturnFlagName(ls_flag);
585#elif CHASTE_SUNDIALS_VERSION >= 20400
586 CVDlsGetLastFlag(mpCvodeMem, &ls_flag);
587 p_ls_flag_name = CVDlsGetReturnFlagName(ls_flag);
588#else
589 CVDenseGetLastFlag(mpCvodeMem, &ls_flag);
590 p_ls_flag_name = CVDenseGetReturnFlagName(ls_flag);
591#endif
592 err << " (LS flag=" << ls_flag << ":" << p_ls_flag_name << ")";
593 free(p_ls_flag_name);
594 }
595
596 err << "\nGot from time " << rStartTime << " to time " << rTime << ", was supposed to finish at time " << rEndTime << "\n";
597 err << "\nState variables are now:\n";
598 std::vector<double> state_vars = MakeStdVec(mStateVariables);
599 std::vector<std::string> state_var_names = rGetStateVariableNames();
600 for (unsigned i = 0; i < state_vars.size(); i++)
601 {
602 err << "\t" << state_var_names[i] << "\t:\t" << state_vars[i] << std::endl;
603 }
604
606 std::cerr << err.str() << std::endl
607 << std::flush;
608 EXCEPTION(err.str());
609}
610
615
620
622{
623 if (!useNumericalJacobian && !mHasAnalyticJacobian)
624 {
625 EXCEPTION("Analytic Jacobian requested, but this ODE system doesn't have one. You can check this with HasAnalyticJacobian().");
626 }
627
628 if (mUseAnalyticJacobian == useNumericalJacobian)
629 {
630 mUseAnalyticJacobian = !useNumericalJacobian;
631 // We need to re-initialise the solver completely to change this.
632 this->FreeCvodeMemory();
633 }
634}
635
636//#include "MathsCustomFunctions.hpp"
637//#include <algorithm>
638//void AbstractCvodeSystem::CheckAnalyticJacobian(realtype time, N_Vector y, N_Vector ydot,
639// CHASTE_CVODE_DENSE_MATRIX jacobian,
640// N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
641//{
642// N_Vector nudge_ydot = tmp1;
643// N_Vector numeric_jth_col = tmp2;
644// N_Vector ewt = tmp3;
645// const unsigned size = GetNumberOfStateVariables();
646// const double rel_tol = 1e-1;
647// const double abs_tol = 1e-6;
648// realtype* p_y = N_VGetArrayPointer(y);
649// realtype* p_numeric_jth_col = N_VGetArrayPointer(numeric_jth_col);
650//
651// // CVODE internal data for computing the numeric J
652// realtype h;
653// CVodeGetLastStep(mpCvodeMem, &h);
654// CVodeGetErrWeights(mpCvodeMem, ewt);
655// realtype* p_ewt = N_VGetArrayPointer(ewt);
656// // Compute minimum nudge
657// realtype srur = sqrt(DBL_EPSILON);
658// realtype fnorm = N_VWrmsNorm(ydot, ewt);
659// realtype min_nudge = (fnorm != 0.0) ?
660// (1000.0 * fabs(h) * DBL_EPSILON * size * fnorm) : 1.0;
661//
662// for (unsigned j=0; j<size; j++)
663// {
664// // Check the j'th column of the Jacobian
665// realtype yjsaved = p_y[j];
666// realtype nudge = std::max(srur*fabs(yjsaved), min_nudge/p_ewt[j]);
667// p_y[j] += nudge;
668// EvaluateYDerivatives(time, y, nudge_ydot);
669// p_y[j] = yjsaved;
670// realtype nudge_inv = 1.0 / nudge;
671// N_VLinearSum(nudge_inv, nudge_ydot, -nudge_inv, ydot, numeric_jth_col);
672// realtype* p_analytic_jth_col = DENSE_COL(jacobian, j);
673//
674// for (unsigned i=0; i<size; i++)
675// {
676// if (!CompareDoubles::WithinAnyTolerance(p_numeric_jth_col[i], p_analytic_jth_col[i], rel_tol, abs_tol))
677// {
678// EXCEPTION("Analytic Jacobian appears dodgy at time " << time << " entry (" << i << "," << j << ").\n"
679// << "Analytic=" << p_analytic_jth_col[i] << "; numeric=" << p_numeric_jth_col[i] << "."
680// << DumpState("", y, time));
681// }
682// }
683// }
684//}
685
686#endif // CHASTE_CVODE
#define EXCEPTION(message)
#define UNUSED_OPT(var)
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)
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)
void RecordStoppingPoint(double stopTime)
void SetMaxSteps(long int numSteps)
AbstractCvodeSystem(unsigned numberOfStateVariables)
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
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 SetMinimalReset(bool minimalReset)
bool GetForceReset()
Get whether we will force a solver reset on every call to Solve()
const std::vector< std::string > & rGetParameterNames() const
boost::shared_ptr< AbstractOdeSystemInformation > mpSystemInfo
const std::vector< std::string > & rGetStateVariableNames() 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)
bool IsTimeAtEnd() const
unsigned GetTotalTimeStepsTaken() const
double GetTime() const
void AdvanceOneTimeStep()
double GetNextTime() const
unsigned EstimateTimeSteps() const