Chaste  Release::2018.1
AbstractCvodeSystem.cpp
1 /*
2 
3 Copyright (c) 2005-2018, 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  mForceMinimalReset = minimalReset;
364  if (mForceMinimalReset)
365  {
366  SetForceReset(false);
367  }
368 }
369 
371 {
373 }
374 
376  realtype tStart,
377  realtype maxDt)
378 {
379  assert((unsigned)NV_LENGTH_S(initialConditions) == GetNumberOfStateVariables());
380  assert(maxDt >= 0.0);
381 
382  // Find out if we need to (re-)initialise
383  //std::cout << "!mpCvodeMem = " << !mpCvodeMem << ", mForceReset = " << mForceReset << ", !mLastSolutionState = " << !mLastSolutionState << ", comp doubles = " << !CompareDoubles::WithinAnyTolerance(tStart, mLastSolutionTime) << "\n";
385  if (!reinit && !mForceMinimalReset)
386  {
387  const unsigned size = GetNumberOfStateVariables();
388  for (unsigned i = 0; i < size; i++)
389  {
391  {
392  reinit = true;
393  break;
394  }
395  }
396  }
397 
398  if (!mpCvodeMem)
399  {
400  //std::cout << "New CVODE solver\n";
401  mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
402  if (mpCvodeMem == nullptr) EXCEPTION("Failed to SetupCvode CVODE"); // in one line to avoid coverage problem!
403 
404  // Set error handler
405  CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, nullptr);
406 // Set the user data
407 #if CHASTE_SUNDIALS_VERSION >= 20400
408  CVodeSetUserData(mpCvodeMem, (void*)(this));
409 #else
410  CVodeSetFdata(mpCvodeMem, (void*)(this));
411 #endif
412 // Setup CVODE
413 #if CHASTE_SUNDIALS_VERSION >= 20400
414  CVodeInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions);
415  CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
416 #else
417  CVodeMalloc(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
418  CV_SS, mRelTol, &mAbsTol);
419 #endif
420 
421 #if CHASTE_SUNDIALS_VERSION >= 30000
422  /* Create dense SUNMatrix for use in linear solves */
423  mpSundialsDenseMatrix = SUNDenseMatrix(NV_LENGTH_S(initialConditions), NV_LENGTH_S(initialConditions));
424 
425  /* Create dense SUNLinearSolver object for use by CVode */
426  mpSundialsLinearSolver = SUNDenseLinearSolver(initialConditions, mpSundialsDenseMatrix);
427 
428  /* Call CVDlsSetLinearSolver to attach the matrix and linear solver to CVode */
429  CVDlsSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
430 #else
431  // Attach a linear solver for Newton iteration
432  CVDense(mpCvodeMem, NV_LENGTH_S(initialConditions));
433 #endif
434 
436  {
437 #if CHASTE_SUNDIALS_VERSION >= 30000
438  CVDlsSetJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor);
439 #elif CHASTE_SUNDIALS_VERSION >= 20400
440  CVDlsSetDenseJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor);
441 #else
442  CVDenseSetJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor, (void*)(this));
443 #endif
444  }
445  }
446  else if (reinit)
447  {
448 //std::cout << "Resetting CVODE solver\n";
449 #if CHASTE_SUNDIALS_VERSION >= 20400
450  CVodeReInit(mpCvodeMem, tStart, initialConditions);
451  CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
452 #else
453  CVodeReInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
454  CV_SS, mRelTol, &mAbsTol);
455 #endif
456  }
457 
458  // Set max dt and change max steps if wanted
459  CVodeSetMaxStep(mpCvodeMem, maxDt);
460  if (mMaxSteps > 0)
461  {
462  CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
463  CVodeSetMaxErrTestFails(mpCvodeMem, 15);
464  }
465 }
466 
468 {
469  // DebugSteps(mpCvodeMem, this);
470 
471  // If we're forcing a reset then we don't record the stopping time
472  // as a result it won't match and we will force a reset in SetupCvode() on
473  // the next solve call.
474  if (mForceReset)
475  return;
476 
477  // Otherwise we will store the state variables and time for comparison on the
478  // next solve call, to work out whether we need to reset.
479  const unsigned size = GetNumberOfStateVariables();
481  for (unsigned i = 0; i < size; i++)
482  {
484  }
485  mLastSolutionTime = stopTime;
486 }
487 
489 {
490  if (mpCvodeMem)
491  {
492  CVodeFree(&mpCvodeMem);
493  }
494  mpCvodeMem = nullptr;
495 
496 #if CHASTE_SUNDIALS_VERSION >= 30000
497  if (mpSundialsLinearSolver)
498  {
499  /* Free the linear solver memory */
500  SUNLinSolFree(mpSundialsLinearSolver);
501  }
502  mpSundialsLinearSolver = nullptr;
503 
504  if (mpSundialsDenseMatrix)
505  {
506  /* Free the matrix memory */
507  SUNMatDestroy(mpSundialsDenseMatrix);
508  }
509  mpSundialsDenseMatrix = nullptr;
510 #endif
511 }
512 
513 void AbstractCvodeSystem::CvodeError(int flag, const char* msg,
514  const double& rTime, const double& rStartTime, const double& rEndTime)
515 {
516  std::stringstream err;
517  char* p_flag_name = CVodeGetReturnFlagName(flag);
518  err << msg << ": " << p_flag_name;
519  free(p_flag_name);
520  if (flag == CV_LSETUP_FAIL)
521  {
522 #if CHASTE_SUNDIALS_VERSION >= 20500
523  long int ls_flag;
524 #else
525  int ls_flag;
526 #endif
527  char* p_ls_flag_name;
528 #if CHASTE_SUNDIALS_VERSION >= 20400
529  CVDlsGetLastFlag(mpCvodeMem, &ls_flag);
530  p_ls_flag_name = CVDlsGetReturnFlagName(ls_flag);
531 #else
532  CVDenseGetLastFlag(mpCvodeMem, &ls_flag);
533  p_ls_flag_name = CVDenseGetReturnFlagName(ls_flag);
534 #endif
535  err << " (LS flag=" << ls_flag << ":" << p_ls_flag_name << ")";
536  free(p_ls_flag_name);
537  }
538 
539  err << "\nGot from time " << rStartTime << " to time " << rTime << ", was supposed to finish at time " << rEndTime << "\n";
540  err << "\nState variables are now:\n";
541  std::vector<double> state_vars = MakeStdVec(mStateVariables);
542  std::vector<std::string> state_var_names = rGetStateVariableNames();
543  for (unsigned i = 0; i < state_vars.size(); i++)
544  {
545  err << "\t" << state_var_names[i] << "\t:\t" << state_vars[i] << std::endl;
546  }
547 
548  FreeCvodeMemory();
549  std::cerr << err.str() << std::endl
550  << std::flush;
551  EXCEPTION(err.str());
552 }
553 
555 {
556  return mHasAnalyticJacobian;
557 }
558 
560 {
561  return mUseAnalyticJacobian;
562 }
563 
565 {
566  if (!useNumericalJacobian && !mHasAnalyticJacobian)
567  {
568  EXCEPTION("Analytic Jacobian requested, but this ODE system doesn't have one. You can check this with HasAnalyticJacobian().");
569  }
570 
571  if (mUseAnalyticJacobian == useNumericalJacobian)
572  {
573  mUseAnalyticJacobian = !useNumericalJacobian;
574  // We need to re-initialise the solver completely to change this.
575  this->FreeCvodeMemory();
576  }
577 }
578 
579 //#include "MathsCustomFunctions.hpp"
580 //#include <algorithm>
581 //void AbstractCvodeSystem::CheckAnalyticJacobian(realtype time, N_Vector y, N_Vector ydot,
582 // CHASTE_CVODE_DENSE_MATRIX jacobian,
583 // N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
584 //{
585 // N_Vector nudge_ydot = tmp1;
586 // N_Vector numeric_jth_col = tmp2;
587 // N_Vector ewt = tmp3;
588 // const unsigned size = GetNumberOfStateVariables();
589 // const double rel_tol = 1e-1;
590 // const double abs_tol = 1e-6;
591 // realtype* p_y = N_VGetArrayPointer(y);
592 // realtype* p_numeric_jth_col = N_VGetArrayPointer(numeric_jth_col);
593 //
594 // // CVODE internal data for computing the numeric J
595 // realtype h;
596 // CVodeGetLastStep(mpCvodeMem, &h);
597 // CVodeGetErrWeights(mpCvodeMem, ewt);
598 // realtype* p_ewt = N_VGetArrayPointer(ewt);
599 // // Compute minimum nudge
600 // realtype srur = sqrt(DBL_EPSILON);
601 // realtype fnorm = N_VWrmsNorm(ydot, ewt);
602 // realtype min_nudge = (fnorm != 0.0) ?
603 // (1000.0 * fabs(h) * DBL_EPSILON * size * fnorm) : 1.0;
604 //
605 // for (unsigned j=0; j<size; j++)
606 // {
607 // // Check the j'th column of the Jacobian
608 // realtype yjsaved = p_y[j];
609 // realtype nudge = std::max(srur*fabs(yjsaved), min_nudge/p_ewt[j]);
610 // p_y[j] += nudge;
611 // EvaluateYDerivatives(time, y, nudge_ydot);
612 // p_y[j] = yjsaved;
613 // realtype nudge_inv = 1.0 / nudge;
614 // N_VLinearSum(nudge_inv, nudge_ydot, -nudge_inv, ydot, numeric_jth_col);
615 // realtype* p_analytic_jth_col = DENSE_COL(jacobian, j);
616 //
617 // for (unsigned i=0; i<size; i++)
618 // {
619 // if (!CompareDoubles::WithinAnyTolerance(p_numeric_jth_col[i], p_analytic_jth_col[i], rel_tol, abs_tol))
620 // {
621 // EXCEPTION("Analytic Jacobian appears dodgy at time " << time << " entry (" << i << "," << j << ").\n"
622 // << "Analytic=" << p_analytic_jth_col[i] << "; numeric=" << p_numeric_jth_col[i] << "."
623 // << DumpState("", y, time));
624 // }
625 // }
626 // }
627 //}
628 
629 #endif // CHASTE_CVODE
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()
bool GetUseAnalyticJacobian() const
double GetTime() const
void RecordStoppingPoint(double stopTime)
#define EXCEPTION(message)
Definition: Exception.hpp:143
void SetupCvode(N_Vector initialConditions, realtype tStart, realtype maxDt)
const std::vector< std::string > & rGetStateVariableNames() const
void AdvanceOneTimeStep()
bool IsTimeAtEnd() 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 HasAnalyticJacobian() const
void SetForceReset(bool autoReset)
virtual void EvaluateYDerivatives(realtype time, const N_Vector y, N_Vector ydot)=0
std::string GetMessage() const
Definition: Exception.cpp:83
void SetVectorComponent(VECTOR &rVec, unsigned index, double value)
unsigned GetTotalTimeStepsTaken() const
void SetMaxSteps(long int numSteps)
std::vector< double > & rGetTimes()
unsigned EstimateTimeSteps() const
#define UNUSED_OPT(var)
Definition: Exception.hpp:216
boost::shared_ptr< AbstractOdeSystemInformation > mpSystemInfo
OdeSolution Solve(realtype tStart, realtype tEnd, realtype maxDt, realtype tSamp)
std::string CheckShortMessageContains(std::string expected) const
Definition: Exception.cpp:103
double GetNextTime() const
double GetVectorComponent(const VECTOR &rVec, unsigned index)
void ForceUseOfNumericalJacobian(bool useNumericalJacobian=true)
AbstractCvodeSystem(unsigned numberOfStateVariables)
void SetNumberOfTimeSteps(unsigned numTimeSteps)
Definition: OdeSolution.cpp:58
const std::vector< std::string > & rGetParameterNames() const
void DeleteVector(VECTOR &rVec)