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