Chaste  Release::2024.1
AbstractCvodeSystem.cpp
1 /*
2 
3 Copyright (c) 2005-2021, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, 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 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF 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 //#include "Debug.hpp"
62 //void DebugSteps(void* pCvodeMem, AbstractCvodeSystem* pSys)
63 //{
64 // long int num_jac_evals, nniters, num_steps;
65 // CVDenseGetNumJacEvals(pCvodeMem, &num_jac_evals);
67 // CVodeGetNumNonlinSolvIters(pCvodeMem, &nniters);
68 // CVodeGetNumSteps(pCvodeMem, &num_steps);
69 // double num_newton_iters = nniters;
70 // PRINT_3_VARIABLES(pSys->GetSystemName(), num_newton_iters/num_steps, num_jac_evals/num_newton_iters);
71 //}
81 int AbstractCvodeSystemRhsAdaptor(realtype t, N_Vector y, N_Vector ydot, void* pData)
82 {
83  assert(pData != nullptr);
84  AbstractCvodeSystem* p_ode_system = (AbstractCvodeSystem*)pData;
85  try
86  {
87  p_ode_system->EvaluateYDerivatives(t, y, ydot);
88  }
89  catch (const Exception& e)
90  {
91 #if CHASTE_SUNDIALS_VERSION <= 20300
92  // Really old CVODE used to solve past the requested time points and could trigger this exception unnecessarily...
93  if (e.CheckShortMessageContains("is outside the times stored in the data clamp") == "")
94  {
95  return 1; // This may be a recoverable error!
96  }
97 #endif
98 
99  std::cerr << "CVODE RHS Exception: " << e.GetMessage()
100  << std::endl
101  << std::flush;
102  return -1;
103  }
104 
105  // // Something like this might help CVODE when things are a bit unstable...
106  // try
107  // {
108  // p_ode_system->VerifyStateVariables();
109  // }
110  // catch (const Exception &e)
111  // {
112  // std::cout << "t = " << t << ":\t" << e.GetMessage() << std::endl << std::flush;
113  // return 1; // A positive return flag to CVODE tells it there's been an error but it might be recoverable.
114  // }
115 
116  return 0;
117 }
118 
119 /*
120  * Absolute chaos here with four different possible interfaces to the jacobian.
121  */
122 #if CHASTE_SUNDIALS_VERSION >= 30000
123 // Sundials 3.0 - has taken away the argument N at the top...
124 int AbstractCvodeSystemJacAdaptor(realtype t, N_Vector y, N_Vector ydot, CHASTE_CVODE_DENSE_MATRIX jacobian,
125 #elif CHASTE_SUNDIALS_VERSION >= 20500
126 // Sundials 2.5
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
129 // Sundials 2.4
130 int AbstractCvodeSystemJacAdaptor(int N, realtype t, N_Vector y, N_Vector ydot, DlsMat jacobian,
131 #else
132 // Sundials 2.3 and below (not sure how far below, but this is 2006 so old enough).
133 int AbstractCvodeSystemJacAdaptor(long int N, DenseMat jacobian, realtype t, N_Vector y, N_Vector ydot,
134 #endif
135  void* pData, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
136 {
137  assert(pData != nullptr);
138  AbstractCvodeSystem* p_ode_system = (AbstractCvodeSystem*)pData;
139  try
140  {
141  p_ode_system->EvaluateAnalyticJacobian(t, y, ydot, jacobian, tmp1, tmp2, tmp3);
142  }
143  catch (const Exception& e)
144  {
145  std::cerr << "CVODE Jacobian Exception: " << e.GetMessage() << std::endl
146  << std::flush;
147  return -1;
148  }
149  return 0;
150 }
151 
152 AbstractCvodeSystem::AbstractCvodeSystem(unsigned numberOfStateVariables)
153  : AbstractParameterisedSystem<N_Vector>(numberOfStateVariables),
154  mLastSolutionState(nullptr),
155  mLastSolutionTime(0.0),
156 #if CHASTE_SUNDIALS_VERSION >= 20400
157  mForceReset(false),
158 #else
159  // Old Sundials don't seem to 'go back' when something has changed
160  // properly, and give more inaccurate answers.
161  mForceReset(true),
162 #endif
163  mForceMinimalReset(false),
164 #if CHASTE_SUNDIALS_VERSION >= 30000
165  mpSundialsDenseMatrix(nullptr),
166  mpSundialsLinearSolver(nullptr),
167 #endif
168  mHasAnalyticJacobian(false),
169  mUseAnalyticJacobian(false),
170  mpCvodeMem(nullptr),
171  mMaxSteps(0),
172  mLastInternalStepSize(0)
173 {
174  SetTolerances(); // Set the tolerances to the defaults.
175 }
176 
178 {
182  mParameters = N_VNew_Serial(rGetParameterNames().size());
183  for (int i = 0; i < NV_LENGTH_S(mParameters); i++)
184  {
185  NV_Ith_S(mParameters, i) = 0.0;
186  }
187 }
188 
190 {
191  FreeCvodeMemory();
195 }
196 
197 //
198 //double AbstractCvodeSystem::CalculateRootFunction(double time, const std::vector<double>& rY)
199 //{
200 // bool stop = CalculateStoppingEvent(time, rY);
201 // return stop ? 0.0 : 1.0;
202 //}
203 
205  realtype tEnd,
206  realtype maxDt,
207  realtype tSamp)
208 {
209  assert(tEnd >= tStart);
210  assert(tSamp > 0.0);
211 
212  SetupCvode(mStateVariables, tStart, maxDt);
213 
214  TimeStepper stepper(tStart, tEnd, tSamp);
215 
216  // Set up ODE solution
217  OdeSolution solutions;
218  solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
219  solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
220  solutions.rGetTimes().push_back(tStart);
222 
223  // Main time sampling loop
224  while (!stepper.IsTimeAtEnd())
225  {
226  // This should stop CVODE going past the end of where we wanted and interpolating back.
227  int ierr = CVodeSetStopTime(mpCvodeMem, stepper.GetNextTime());
228  assert(ierr == CV_SUCCESS);
229  UNUSED_OPT(ierr); // avoid unused var warning
230 
231  // // This parameter governs how many times we allow a recoverable right hand side failure
232  // int ierr = CVodeSetMaxConvFails(mpCvodeMem, 1000);
233  // assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr); // avoid unused var warning
234 
235  double cvode_stopped_at = stepper.GetTime();
236  ierr = CVode(mpCvodeMem, stepper.GetNextTime(), mStateVariables,
237  &cvode_stopped_at, CV_NORMAL);
238  if (ierr < 0)
239  {
240  // DebugSteps(mpCvodeMem, this);
241  CvodeError(ierr, "CVODE failed to solve system", cvode_stopped_at, stepper.GetTime(), stepper.GetNextTime());
242  }
243  // Not root finding, so should have reached requested time
244  assert(fabs(cvode_stopped_at - stepper.GetNextTime()) < DBL_EPSILON);
245 #ifndef NDEBUG
247 #endif
248  // Store solution
249  solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
250  solutions.rGetTimes().push_back(cvode_stopped_at);
251  stepper.AdvanceOneTimeStep();
252  }
253 
254  // stepper.EstimateTimeSteps may have been an overestimate...
255  solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
256 
257  int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
258  assert(ierr == CV_SUCCESS);
259  UNUSED_OPT(ierr); // avoid unused var warning
260 
261  RecordStoppingPoint(tEnd);
262 
263  return solutions;
264 }
265 
266 void AbstractCvodeSystem::Solve(realtype tStart,
267  realtype tEnd,
268  realtype maxDt)
269 {
270  assert(tEnd >= tStart);
271 
272  SetupCvode(mStateVariables, tStart, maxDt);
273 
274  // This should stop CVODE going past the end of where we wanted and interpolating back.
275  int ierr = CVodeSetStopTime(mpCvodeMem, tEnd);
276  assert(ierr == CV_SUCCESS);
277  UNUSED_OPT(ierr); // avoid unused var warning
278 
279  double cvode_stopped_at = tStart;
280  ierr = CVode(mpCvodeMem, tEnd, mStateVariables, &cvode_stopped_at, CV_NORMAL);
281  if (ierr < 0)
282  {
283  // DebugSteps(mpCvodeMem, this);
284  CvodeError(ierr, "CVODE failed to solve system", cvode_stopped_at, tStart, tEnd);
285  }
286  // Not root finding, so should have reached requested time
287  assert(fabs(cvode_stopped_at - tEnd) < DBL_EPSILON);
288 
289  ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
290  assert(ierr == CV_SUCCESS);
291  UNUSED_OPT(ierr); // avoid unused var warning
292 
293  RecordStoppingPoint(cvode_stopped_at);
294 
295 //
296 // long int nst, nfe, nsetups, nje, nfeLS, nni, ncfn, netf, nge;
297 //
298 //
299 // CVodeGetNumSteps(mpCvodeMem, &nst);
300 // CVodeGetNumRhsEvals(mpCvodeMem, &nfe);
301 // CVodeGetNumLinSolvSetups(mpCvodeMem, &nsetups);
302 // CVodeGetNumErrTestFails(mpCvodeMem, &netf);
303 // CVodeGetNumNonlinSolvIters(mpCvodeMem, &nni);
304 // CVodeGetNumNonlinSolvConvFails(mpCvodeMem, &ncfn);
305 // CVDlsGetNumJacEvals(mpCvodeMem, &nje);
306 // CVDlsGetNumRhsEvals(mpCvodeMem, &nfeLS);
307 // CVodeGetNumGEvals(mpCvodeMem, &nge);
308 //
309 // printf("\nFinal Statistics:\n");
310 // printf("nst = %-6ld nfe = %-6ld nsetups = %-6ld nfeLS = %-6ld nje = %ld\n",
311 // nst, nfe, nsetups, nfeLS, nje);
312 // printf("nni = %-6ld ncfn = %-6ld netf = %-6ld nge = %ld\n \n",
313 // nni, ncfn, netf, nge);
314 // std::cout << std::flush;
315 #ifndef NDEBUG
317 #endif
318 }
319 
320 void AbstractCvodeSystem::SetMaxSteps(long int numSteps)
321 {
322  mMaxSteps = numSteps;
323 }
324 
326 {
327  return mMaxSteps;
328 }
329 
330 void AbstractCvodeSystem::SetTolerances(double relTol, double absTol)
331 {
332  mRelTol = relTol;
333  mAbsTol = absTol;
334  ResetSolver();
335 }
336 
338 {
339  return mRelTol;
340 }
341 
343 {
344  return mAbsTol;
345 }
346 
348 {
349  return mLastInternalStepSize;
350 }
351 
353 {
354  mForceReset = autoReset;
355  if (mForceReset)
356  {
357  ResetSolver();
358  }
359 }
360 
362 {
363  return mForceMinimalReset;
364 }
365 
367 {
368  return mForceReset;
369 }
370 
372 {
373  mForceMinimalReset = minimalReset;
374  if (mForceMinimalReset)
375  {
376  SetForceReset(false);
377  }
378 }
379 
381 {
383 }
384 
386  realtype tStart,
387  realtype maxDt)
388 {
389  assert((unsigned)NV_LENGTH_S(initialConditions) == GetNumberOfStateVariables());
390  assert(maxDt >= 0.0);
391 
392  // Find out if we need to (re-)initialise
393  //std::cout << "!mpCvodeMem = " << !mpCvodeMem << ", mForceReset = " << mForceReset << ", !mLastSolutionState = " << !mLastSolutionState << ", comp doubles = " << !CompareDoubles::WithinAnyTolerance(tStart, mLastSolutionTime) << "\n";
395  if (!reinit && !mForceMinimalReset)
396  {
397  const unsigned size = GetNumberOfStateVariables();
398  for (unsigned i = 0; i < size; i++)
399  {
401  {
402  reinit = true;
403  break;
404  }
405  }
406  }
407 
408  if (!mpCvodeMem)
409  {
410  //std::cout << "New CVODE solver\n";
411 #if CHASTE_SUNDIALS_VERSION >= 40000
412  // v4.0.0 release notes: instead of specifying the nonlinear iteration type when creating the CVODE(S) memory structure,
413  // CVODE(S) uses the SUNNONLINSOL_NEWTON module implementation of a Newton iteration by default.
414  mpCvodeMem = CVodeCreate(CV_BDF);
415 #else
416  mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
417 #endif
418  if (mpCvodeMem == nullptr)
419  EXCEPTION("Failed to SetupCvode CVODE"); // LCOV_EXCL_LINE
420 
421  // Set error handler
422  CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, nullptr);
423 // Set the user data
424 #if CHASTE_SUNDIALS_VERSION >= 20400
425  CVodeSetUserData(mpCvodeMem, (void*)(this));
426 #else
427  CVodeSetFdata(mpCvodeMem, (void*)(this));
428 #endif
429 // Setup CVODE
430 #if CHASTE_SUNDIALS_VERSION >= 20400
431  CVodeInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions);
432  CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
433 #else
434  CVodeMalloc(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
435  CV_SS, mRelTol, &mAbsTol);
436 #endif
437 
438 #if CHASTE_SUNDIALS_VERSION >= 30000
439  /* Create dense matrix SUNDenseMatrix for use in linear solves */
440  mpSundialsDenseMatrix = SUNDenseMatrix(NV_LENGTH_S(initialConditions), NV_LENGTH_S(initialConditions));
441 #endif
442 
443 #if CHASTE_SUNDIALS_VERSION >= 40000
444  /* Create dense SUNLinSol_Dense object for use by CVode */
445  mpSundialsLinearSolver = SUNLinSol_Dense(initialConditions, mpSundialsDenseMatrix);
446 
447  /* Call CVodeSetLinearSolver to attach the matrix and linear solver to CVode */
448  CVodeSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
449 #elif CHASTE_SUNDIALS_VERSION >= 30000
450  /* Create dense SUNDenseLinearSolver object for use by CVode */
451  mpSundialsLinearSolver = SUNDenseLinearSolver(initialConditions, mpSundialsDenseMatrix);
452 
453  /* Call CVDlsSetLinearSolver to attach the matrix and linear solver to CVode */
454  CVDlsSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
455 #else
456  // CVODE < v3.0.0
457  // Attach a linear solver for Newton iteration
458  CVDense(mpCvodeMem, NV_LENGTH_S(initialConditions));
459 #endif
460 
462  {
463 #if CHASTE_SUNDIALS_VERSION >= 40000
464  CVodeSetJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor);
465 #elif CHASTE_SUNDIALS_VERSION >= 30000
466  CVDlsSetJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor);
467 #elif CHASTE_SUNDIALS_VERSION >= 20400
468  CVDlsSetDenseJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor);
469 #else
470  CVDenseSetJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor, (void*)(this));
471 #endif
472  }
473  }
474  else if (reinit)
475  {
476 //std::cout << "Resetting CVODE solver\n";
477 #if CHASTE_SUNDIALS_VERSION >= 20400
478  CVodeReInit(mpCvodeMem, tStart, initialConditions);
479  //CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol); - "all solver inputs remain in effect" so we don't need this.
480 #else
481  CVodeReInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
482  CV_SS, mRelTol, &mAbsTol);
483 #endif
484  }
485 
486  // Set max dt and change max steps if wanted
487  if (maxDt > 0)
488  {
489  CVodeSetMaxStep(mpCvodeMem, maxDt);
490  }
491 
492  if (mMaxSteps > 0)
493  {
494  CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
495  CVodeSetMaxErrTestFails(mpCvodeMem, 15);
496  }
497 }
498 
500 {
501  // DebugSteps(mpCvodeMem, this);
502 
503  // If we're forcing a reset then we don't record the stopping time
504  // as a result it won't match and we will force a reset in SetupCvode() on
505  // the next solve call.
506  if (mForceReset)
507  return;
508 
509  // Otherwise we will store the state variables and time for comparison on the
510  // next solve call, to work out whether we need to reset.
511  const unsigned size = GetNumberOfStateVariables();
513  for (unsigned i = 0; i < size; i++)
514  {
516  }
517  mLastSolutionTime = stopTime;
518 }
519 
521 {
522  if (mpCvodeMem)
523  {
524  CVodeFree(&mpCvodeMem);
525  }
526  mpCvodeMem = nullptr;
527 
528 #if CHASTE_SUNDIALS_VERSION >= 30000
529  if (mpSundialsLinearSolver)
530  {
531  /* Free the linear solver memory */
532  SUNLinSolFree(mpSundialsLinearSolver);
533  }
534  mpSundialsLinearSolver = nullptr;
535 
536  if (mpSundialsDenseMatrix)
537  {
538  /* Free the matrix memory */
539  SUNMatDestroy(mpSundialsDenseMatrix);
540  }
541  mpSundialsDenseMatrix = nullptr;
542 #endif
543 }
544 
545 void AbstractCvodeSystem::CvodeError(int flag, const char* msg,
546  const double& rTime, const double& rStartTime, const double& rEndTime)
547 {
548  std::stringstream err;
549  char* p_flag_name = CVodeGetReturnFlagName(flag);
550  err << msg << ": " << p_flag_name;
551  free(p_flag_name);
552  if (flag == CV_LSETUP_FAIL)
553  {
554 #if CHASTE_SUNDIALS_VERSION >= 20500
555  long int ls_flag;
556 #else
557  int ls_flag;
558 #endif
559  char* p_ls_flag_name;
560 
561 #if CHASTE_SUNDIALS_VERSION >= 40000
562  CVodeGetLastLinFlag(mpCvodeMem, &ls_flag);
563  p_ls_flag_name = CVodeGetLinReturnFlagName(ls_flag);
564 #elif CHASTE_SUNDIALS_VERSION >= 20400
565  CVDlsGetLastFlag(mpCvodeMem, &ls_flag);
566  p_ls_flag_name = CVDlsGetReturnFlagName(ls_flag);
567 #else
568  CVDenseGetLastFlag(mpCvodeMem, &ls_flag);
569  p_ls_flag_name = CVDenseGetReturnFlagName(ls_flag);
570 #endif
571  err << " (LS flag=" << ls_flag << ":" << p_ls_flag_name << ")";
572  free(p_ls_flag_name);
573  }
574 
575  err << "\nGot from time " << rStartTime << " to time " << rTime << ", was supposed to finish at time " << rEndTime << "\n";
576  err << "\nState variables are now:\n";
577  std::vector<double> state_vars = MakeStdVec(mStateVariables);
578  std::vector<std::string> state_var_names = rGetStateVariableNames();
579  for (unsigned i = 0; i < state_vars.size(); i++)
580  {
581  err << "\t" << state_var_names[i] << "\t:\t" << state_vars[i] << std::endl;
582  }
583 
584  FreeCvodeMemory();
585  std::cerr << err.str() << std::endl
586  << std::flush;
587  EXCEPTION(err.str());
588 }
589 
591 {
592  return mHasAnalyticJacobian;
593 }
594 
596 {
597  return mUseAnalyticJacobian;
598 }
599 
601 {
602  if (!useNumericalJacobian && !mHasAnalyticJacobian)
603  {
604  EXCEPTION("Analytic Jacobian requested, but this ODE system doesn't have one. You can check this with HasAnalyticJacobian().");
605  }
606 
607  if (mUseAnalyticJacobian == useNumericalJacobian)
608  {
609  mUseAnalyticJacobian = !useNumericalJacobian;
610  // We need to re-initialise the solver completely to change this.
611  this->FreeCvodeMemory();
612  }
613 }
614 
615 //#include "MathsCustomFunctions.hpp"
616 //#include <algorithm>
617 //void AbstractCvodeSystem::CheckAnalyticJacobian(realtype time, N_Vector y, N_Vector ydot,
618 // CHASTE_CVODE_DENSE_MATRIX jacobian,
619 // N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
620 //{
621 // N_Vector nudge_ydot = tmp1;
622 // N_Vector numeric_jth_col = tmp2;
623 // N_Vector ewt = tmp3;
624 // const unsigned size = GetNumberOfStateVariables();
625 // const double rel_tol = 1e-1;
626 // const double abs_tol = 1e-6;
627 // realtype* p_y = N_VGetArrayPointer(y);
628 // realtype* p_numeric_jth_col = N_VGetArrayPointer(numeric_jth_col);
629 //
630 // // CVODE internal data for computing the numeric J
631 // realtype h;
632 // CVodeGetLastStep(mpCvodeMem, &h);
633 // CVodeGetErrWeights(mpCvodeMem, ewt);
634 // realtype* p_ewt = N_VGetArrayPointer(ewt);
635 // // Compute minimum nudge
636 // realtype srur = sqrt(DBL_EPSILON);
637 // realtype fnorm = N_VWrmsNorm(ydot, ewt);
638 // realtype min_nudge = (fnorm != 0.0) ?
639 // (1000.0 * fabs(h) * DBL_EPSILON * size * fnorm) : 1.0;
640 //
641 // for (unsigned j=0; j<size; j++)
642 // {
643 // // Check the j'th column of the Jacobian
644 // realtype yjsaved = p_y[j];
645 // realtype nudge = std::max(srur*fabs(yjsaved), min_nudge/p_ewt[j]);
646 // p_y[j] += nudge;
647 // EvaluateYDerivatives(time, y, nudge_ydot);
648 // p_y[j] = yjsaved;
649 // realtype nudge_inv = 1.0 / nudge;
650 // N_VLinearSum(nudge_inv, nudge_ydot, -nudge_inv, ydot, numeric_jth_col);
651 // realtype* p_analytic_jth_col = DENSE_COL(jacobian, j);
652 //
653 // for (unsigned i=0; i<size; i++)
654 // {
655 // if (!CompareDoubles::WithinAnyTolerance(p_numeric_jth_col[i], p_analytic_jth_col[i], rel_tol, abs_tol))
656 // {
657 // EXCEPTION("Analytic Jacobian appears dodgy at time " << time << " entry (" << i << "," << j << ").\n"
658 // << "Analytic=" << p_analytic_jth_col[i] << "; numeric=" << p_numeric_jth_col[i] << "."
659 // << DumpState("", y, time));
660 // }
661 // }
662 // }
663 //}
664 
665 #endif // CHASTE_CVODE
const std::vector< std::string > & rGetParameterNames() const
bool IsTimeAtEnd() const
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()
void RecordStoppingPoint(double stopTime)
#define EXCEPTION(message)
Definition: Exception.hpp:143
std::string GetMessage() const
Definition: Exception.cpp:83
double GetNextTime() const
void SetupCvode(N_Vector initialConditions, realtype tStart, realtype maxDt)
unsigned GetTotalTimeStepsTaken() const
const std::vector< std::string > & rGetStateVariableNames() const
void AdvanceOneTimeStep()
double GetTime() const
void CreateVectorIfEmpty(VECTOR &rVec, unsigned size)
void SetOdeSystemInformation(boost::shared_ptr< const AbstractOdeSystemInformation > pOdeSystemInfo)
Definition: OdeSolution.cpp:66
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 GetForceReset()
Get whether we will force a solver reset on every call to Solve()
bool GetUseAnalyticJacobian() const
void SetForceReset(bool autoReset)
virtual void EvaluateYDerivatives(realtype time, const N_Vector y, N_Vector ydot)=0
unsigned EstimateTimeSteps() const
void SetVectorComponent(VECTOR &rVec, unsigned index, double value)
void SetMaxSteps(long int numSteps)
std::vector< double > & rGetTimes()
bool GetMinimalReset()
Get whether we want to run with minimal reset or not (no reinitialisation of the solver if variables ...
#define UNUSED_OPT(var)
Definition: Exception.hpp:216
boost::shared_ptr< AbstractOdeSystemInformation > mpSystemInfo
OdeSolution Solve(realtype tStart, realtype tEnd, realtype maxDt, realtype tSamp)
double GetVectorComponent(const VECTOR &rVec, unsigned index)
std::string CheckShortMessageContains(std::string expected) const
Definition: Exception.cpp:103
void ForceUseOfNumericalJacobian(bool useNumericalJacobian=true)
AbstractCvodeSystem(unsigned numberOfStateVariables)
void SetNumberOfTimeSteps(unsigned numTimeSteps)
Definition: OdeSolution.cpp:58
void DeleteVector(VECTOR &rVec)