Chaste  Release::2017.1
AbstractCvodeSystem.cpp
1 /*
2 
3 Copyright (c) 2005-2017, 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 <sstream>
39 #include <cassert>
40 
41 #include "AbstractCvodeSystem.hpp"
42 #include "Exception.hpp"
44 #include "MathsCustomFunctions.hpp" // For tolerance comparison
45 #include "TimeStepper.hpp"
46 #include "CvodeAdaptor.hpp" // For CvodeErrorHandler
47 
48 // CVODE headers
49 #include <cvode/cvode.h>
50 #include <sundials/sundials_nvector.h>
51 #include <cvode/cvode_dense.h>
52 
53 
54 //#include "Debug.hpp"
55 //void DebugSteps(void* pCvodeMem, AbstractCvodeSystem* pSys)
56 //{
57 // long int num_jac_evals, nniters, num_steps;
58 // CVDenseGetNumJacEvals(pCvodeMem, &num_jac_evals);
60 // CVodeGetNumNonlinSolvIters(pCvodeMem, &nniters);
61 // CVodeGetNumSteps(pCvodeMem, &num_steps);
62 // double num_newton_iters = nniters;
63 // PRINT_3_VARIABLES(pSys->GetSystemName(), num_newton_iters/num_steps, num_jac_evals/num_newton_iters);
64 //}
74 int AbstractCvodeSystemRhsAdaptor(realtype t, N_Vector y, N_Vector ydot, void *pData)
75 {
76  assert(pData != nullptr);
77  AbstractCvodeSystem* p_ode_system = (AbstractCvodeSystem*) pData;
78  try
79  {
80  p_ode_system->EvaluateYDerivatives(t, y, ydot);
81  }
82  catch (const Exception &e)
83  {
84 #if CHASTE_SUNDIALS_VERSION <= 20300
85  // Really old CVODE used to solve past the requested time points and could trigger this exception unnecessarily...
86  if (e.CheckShortMessageContains("is outside the times stored in the data clamp")=="")
87  {
88  return 1; // This may be a recoverable error!
89  }
90 #endif
91 
92  std::cerr << "CVODE RHS Exception: " << e.GetMessage()
93  << std::endl << std::flush;
94  return -1;
95  }
96 
97 // // Something like this might help CVODE when things are a bit unstable...
98 // try
99 // {
100 // p_ode_system->VerifyStateVariables();
101 // }
102 // catch (const Exception &e)
103 // {
104 // std::cout << "t = " << t << ":\t" << e.GetMessage() << std::endl << std::flush;
105 // return 1; // A positive return flag to CVODE tells it there's been an error but it might be recoverable.
106 // }
107 
108  return 0;
109 }
110 
111 
112 /*
113  * Absolute chaos here with three different possible interfaces to the jacobian.
114  */
115 #if CHASTE_SUNDIALS_VERSION >= 20500
116  // Sundials 2.5
117  int AbstractCvodeSystemJacAdaptor(long int N, realtype t, N_Vector y, N_Vector ydot, DlsMat jacobian,
118 #elif CHASTE_SUNDIALS_VERSION >= 20400
119  // Sundials 2.4
120  int AbstractCvodeSystemJacAdaptor(int N, realtype t, N_Vector y, N_Vector ydot, DlsMat jacobian,
121 #else
122  // Sundials 2.3 and below (not sure how far below, but this is 2006 so old enough).
123  int AbstractCvodeSystemJacAdaptor(long int N, DenseMat jacobian, realtype t, N_Vector y, N_Vector ydot,
124 #endif
125  void *pData, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
126 {
127  assert(pData != nullptr);
128  AbstractCvodeSystem* p_ode_system = (AbstractCvodeSystem*) pData;
129  try
130  {
131  p_ode_system->EvaluateAnalyticJacobian((long)(N), t, y, ydot, jacobian, tmp1, tmp2, tmp3);
132  }
133  catch (const Exception &e)
134  {
135  std::cerr << "CVODE Jacobian Exception: " << e.GetMessage() << std::endl << std::flush;
136  return -1;
137  }
138  return 0;
139 }
140 
141 AbstractCvodeSystem::AbstractCvodeSystem(unsigned numberOfStateVariables)
142  : AbstractParameterisedSystem<N_Vector>(numberOfStateVariables),
143  mLastSolutionState(nullptr),
144  mLastSolutionTime(0.0),
145 #if CHASTE_SUNDIALS_VERSION >=20400
146  mForceReset(false),
147 #else
148  // Old Sundials don't seem to 'go back' when something has changed
149  // properly, and give more inaccurate answers.
150  mForceReset(true),
151 #endif
152  mForceMinimalReset(false),
153  mHasAnalyticJacobian(false),
154  mUseAnalyticJacobian(false),
155  mpCvodeMem(nullptr),
156  mMaxSteps(0),
157  mLastInternalStepSize(0)
158 {
159  SetTolerances(); // Set the tolerances to the defaults.
160 }
161 
163 {
167  mParameters = N_VNew_Serial(rGetParameterNames().size());
168  for (int i=0; i<NV_LENGTH_S(mParameters); i++)
169  {
170  NV_Ith_S(mParameters, i) = 0.0;
171  }
172 }
173 
175 {
176  FreeCvodeMemory();
180 }
181 
182 //
183 //double AbstractCvodeSystem::CalculateRootFunction(double time, const std::vector<double>& rY)
184 //{
185 // bool stop = CalculateStoppingEvent(time, rY);
186 // return stop ? 0.0 : 1.0;
187 //}
188 
190  realtype tEnd,
191  realtype maxDt,
192  realtype tSamp)
193 {
194  assert(tEnd >= tStart);
195  assert(tSamp > 0.0);
196 
197  SetupCvode(mStateVariables, tStart, maxDt);
198 
199  TimeStepper stepper(tStart, tEnd, tSamp);
200 
201  // Set up ODE solution
202  OdeSolution solutions;
203  solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
204  solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
205  solutions.rGetTimes().push_back(tStart);
207 
208  // Main time sampling loop
209  while (!stepper.IsTimeAtEnd())
210  {
211  // This should stop CVODE going past the end of where we wanted and interpolating back.
212  int ierr = CVodeSetStopTime(mpCvodeMem, stepper.GetNextTime());
213  assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr); // avoid unused var warning
214 
215 // // This parameter governs how many times we allow a recoverable right hand side failure
216 // int ierr = CVodeSetMaxConvFails(mpCvodeMem, 1000);
217 // assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr); // avoid unused var warning
218 
219  double cvode_stopped_at = stepper.GetTime();
220  ierr = CVode(mpCvodeMem, stepper.GetNextTime(), mStateVariables,
221  &cvode_stopped_at, CV_NORMAL);
222  if (ierr<0)
223  {
224 // DebugSteps(mpCvodeMem, this);
225  CvodeError(ierr, "CVODE failed to solve system", cvode_stopped_at, stepper.GetTime(), stepper.GetNextTime());
226  }
227  // Not root finding, so should have reached requested time
228  assert(fabs(cvode_stopped_at - stepper.GetNextTime()) < DBL_EPSILON);
229 #ifndef NDEBUG
231 #endif
232  // Store solution
233  solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
234  solutions.rGetTimes().push_back(cvode_stopped_at);
235  stepper.AdvanceOneTimeStep();
236  }
237 
238  // stepper.EstimateTimeSteps may have been an overestimate...
239  solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
240 
241  int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
242  assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr); // avoid unused var warning
243 
244  RecordStoppingPoint(tEnd);
245 
246  return solutions;
247 }
248 
249 void AbstractCvodeSystem::Solve(realtype tStart,
250  realtype tEnd,
251  realtype maxDt)
252 {
253  assert(tEnd >= tStart);
254 
255  SetupCvode(mStateVariables, tStart, maxDt);
256 
257  // This should stop CVODE going past the end of where we wanted and interpolating back.
258  int ierr = CVodeSetStopTime(mpCvodeMem, tEnd);
259  assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr); // avoid unused var warning
260 
261  double cvode_stopped_at = tStart;
262  ierr = CVode(mpCvodeMem, tEnd, mStateVariables, &cvode_stopped_at, CV_NORMAL);
263  if (ierr<0)
264  {
265 // DebugSteps(mpCvodeMem, this);
266  CvodeError(ierr, "CVODE failed to solve system", cvode_stopped_at, tStart, tEnd);
267  }
268  // Not root finding, so should have reached requested time
269  assert(fabs(cvode_stopped_at - tEnd) < DBL_EPSILON);
270 
271  ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
272  assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr); // avoid unused var warning
273 
274  RecordStoppingPoint(cvode_stopped_at);
275 
276 //
277 // long int nst, nfe, nsetups, nje, nfeLS, nni, ncfn, netf, nge;
278 //
279 //
280 // CVodeGetNumSteps(mpCvodeMem, &nst);
281 // CVodeGetNumRhsEvals(mpCvodeMem, &nfe);
282 // CVodeGetNumLinSolvSetups(mpCvodeMem, &nsetups);
283 // CVodeGetNumErrTestFails(mpCvodeMem, &netf);
284 // CVodeGetNumNonlinSolvIters(mpCvodeMem, &nni);
285 // CVodeGetNumNonlinSolvConvFails(mpCvodeMem, &ncfn);
286 // CVDlsGetNumJacEvals(mpCvodeMem, &nje);
287 // CVDlsGetNumRhsEvals(mpCvodeMem, &nfeLS);
288 // CVodeGetNumGEvals(mpCvodeMem, &nge);
289 //
290 // printf("\nFinal Statistics:\n");
291 // printf("nst = %-6ld nfe = %-6ld nsetups = %-6ld nfeLS = %-6ld nje = %ld\n",
292 // nst, nfe, nsetups, nfeLS, nje);
293 // printf("nni = %-6ld ncfn = %-6ld netf = %-6ld nge = %ld\n \n",
294 // nni, ncfn, netf, nge);
295 // std::cout << std::flush;
296 #ifndef NDEBUG
298 #endif
299 }
300 
301 
302 void AbstractCvodeSystem::SetMaxSteps(long int numSteps)
303 {
304  mMaxSteps = numSteps;
305 }
306 
308 {
309  return mMaxSteps;
310 }
311 
312 void AbstractCvodeSystem::SetTolerances(double relTol, double absTol)
313 {
314  mRelTol = relTol;
315  mAbsTol = absTol;
316  ResetSolver();
317 }
318 
320 {
321  return mRelTol;
322 }
323 
325 {
326  return mAbsTol;
327 }
328 
330 {
331  return mLastInternalStepSize;
332 }
333 
335 {
336  mForceReset = autoReset;
337  if (mForceReset)
338  {
339  ResetSolver();
340  }
341 }
342 
344 {
345  mForceMinimalReset = minimalReset;
346  if (mForceMinimalReset)
347  {
348  SetForceReset(false);
349  }
350 }
351 
352 
354 {
356 }
357 
358 
360  realtype tStart,
361  realtype maxDt)
362 {
363  assert((unsigned)NV_LENGTH_S(initialConditions) == GetNumberOfStateVariables());
364  assert(maxDt >= 0.0);
365 
366  // Find out if we need to (re-)initialise
367  //std::cout << "!mpCvodeMem = " << !mpCvodeMem << ", mForceReset = " << mForceReset << ", !mLastSolutionState = " << !mLastSolutionState << ", comp doubles = " << !CompareDoubles::WithinAnyTolerance(tStart, mLastSolutionTime) << "\n";
369  if (!reinit && !mForceMinimalReset)
370  {
371  const unsigned size = GetNumberOfStateVariables();
372  for (unsigned i=0; i<size; i++)
373  {
375  {
376  reinit = true;
377  break;
378  }
379  }
380  }
381 
382  if (!mpCvodeMem)
383  {
384  //std::cout << "New CVODE solver\n";
385  mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
386  if (mpCvodeMem == nullptr) EXCEPTION("Failed to SetupCvode CVODE"); // in one line to avoid coverage problem!
387 
388  // Set error handler
389  CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, nullptr);
390  // Set the user data
391 #if CHASTE_SUNDIALS_VERSION >= 20400
392  CVodeSetUserData(mpCvodeMem, (void*)(this));
393 #else
394  CVodeSetFdata(mpCvodeMem, (void*)(this));
395 #endif
396  // Setup CVODE
397 #if CHASTE_SUNDIALS_VERSION >= 20400
398  CVodeInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions);
399  CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
400 #else
401  CVodeMalloc(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
402  CV_SS, mRelTol, &mAbsTol);
403 #endif
404  // Attach a linear solver for Newton iteration
405  CVDense(mpCvodeMem, NV_LENGTH_S(initialConditions));
406 
408  {
409 #if CHASTE_SUNDIALS_VERSION >= 20400
410  CVDlsSetDenseJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor);
411 #else
412  CVDenseSetJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor, (void*)(this));
413 #endif
414  }
415  }
416  else if (reinit)
417  {
418  //std::cout << "Resetting CVODE solver\n";
419 #if CHASTE_SUNDIALS_VERSION >= 20400
420  CVodeReInit(mpCvodeMem, tStart, initialConditions);
421  CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
422 #else
423  CVodeReInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
424  CV_SS, mRelTol, &mAbsTol);
425 #endif
426  }
427 
428  // Set max dt and change max steps if wanted
429  CVodeSetMaxStep(mpCvodeMem, maxDt);
430  if (mMaxSteps > 0)
431  {
432  CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
433  CVodeSetMaxErrTestFails(mpCvodeMem, 15);
434  }
435 }
436 
437 
439 {
440 // DebugSteps(mpCvodeMem, this);
441 
442  // If we're forcing a reset then we don't record the stopping time
443  // as a result it won't match and we will force a reset in SetupCvode() on
444  // the next solve call.
445  if (mForceReset) return;
446 
447  // Otherwise we will store the state variables and time for comparison on the
448  // next solve call, to work out whether we need to reset.
449  const unsigned size = GetNumberOfStateVariables();
451  for (unsigned i=0; i<size; i++)
452  {
454  }
455  mLastSolutionTime = stopTime;
456 }
457 
458 
460 {
461  if (mpCvodeMem)
462  {
463  CVodeFree(&mpCvodeMem);
464  }
465  mpCvodeMem = nullptr;
466 }
467 
468 
469 void AbstractCvodeSystem::CvodeError(int flag, const char * msg,
470  const double& rTime, const double& rStartTime, const double& rEndTime)
471 {
472  std::stringstream err;
473  char* p_flag_name = CVodeGetReturnFlagName(flag);
474  err << msg << ": " << p_flag_name;
475  free(p_flag_name);
476  if (flag == CV_LSETUP_FAIL)
477  {
478 #if CHASTE_SUNDIALS_VERSION >= 20500
479  long int ls_flag;
480 #else
481  int ls_flag;
482 #endif
483  char* p_ls_flag_name;
484 #if CHASTE_SUNDIALS_VERSION >= 20400
485  CVDlsGetLastFlag(mpCvodeMem, &ls_flag);
486  p_ls_flag_name = CVDlsGetReturnFlagName(ls_flag);
487 #else
488  CVDenseGetLastFlag(mpCvodeMem, &ls_flag);
489  p_ls_flag_name = CVDenseGetReturnFlagName(ls_flag);
490 #endif
491  err << " (LS flag=" << ls_flag << ":" << p_ls_flag_name << ")";
492  free(p_ls_flag_name);
493  }
494 
495  err << "\nGot from time " << rStartTime << " to time " << rTime << ", was supposed to finish at time " << rEndTime << "\n";
496  err << "\nState variables are now:\n";
497  std::vector<double> state_vars = MakeStdVec(mStateVariables);
498  std::vector<std::string> state_var_names = rGetStateVariableNames();
499  for (unsigned i=0; i<state_vars.size(); i++)
500  {
501  err << "\t" << state_var_names[i] << "\t:\t" << state_vars[i] << std::endl;
502  }
503 
504  FreeCvodeMemory();
505  std::cerr << err.str() << std::endl << std::flush;
506  EXCEPTION(err.str());
507 }
508 
510 {
511  return mHasAnalyticJacobian;
512 }
513 
515 {
516  return mUseAnalyticJacobian;
517 }
518 
519 
521 {
522  if (!useNumericalJacobian && !mHasAnalyticJacobian)
523  {
524  EXCEPTION("Analytic Jacobian requested, but this ODE system doesn't have one. You can check this with HasAnalyticJacobian().");
525  }
526 
527  if (mUseAnalyticJacobian == useNumericalJacobian)
528  {
529  mUseAnalyticJacobian = !useNumericalJacobian;
530  // We need to re-initialise the solver completely to change this.
531  this->FreeCvodeMemory();
532  }
533 }
534 
535 
536 //#include "MathsCustomFunctions.hpp"
537 //#include <algorithm>
538 //void AbstractCvodeSystem::CheckAnalyticJacobian(realtype time, N_Vector y, N_Vector ydot,
539 // CHASTE_CVODE_DENSE_MATRIX jacobian,
540 // N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
541 //{
542 // N_Vector nudge_ydot = tmp1;
543 // N_Vector numeric_jth_col = tmp2;
544 // N_Vector ewt = tmp3;
545 // const unsigned size = GetNumberOfStateVariables();
546 // const double rel_tol = 1e-1;
547 // const double abs_tol = 1e-6;
548 // realtype* p_y = N_VGetArrayPointer(y);
549 // realtype* p_numeric_jth_col = N_VGetArrayPointer(numeric_jth_col);
550 //
551 // // CVODE internal data for computing the numeric J
552 // realtype h;
553 // CVodeGetLastStep(mpCvodeMem, &h);
554 // CVodeGetErrWeights(mpCvodeMem, ewt);
555 // realtype* p_ewt = N_VGetArrayPointer(ewt);
556 // // Compute minimum nudge
557 // realtype srur = sqrt(DBL_EPSILON);
558 // realtype fnorm = N_VWrmsNorm(ydot, ewt);
559 // realtype min_nudge = (fnorm != 0.0) ?
560 // (1000.0 * fabs(h) * DBL_EPSILON * size * fnorm) : 1.0;
561 //
562 // for (unsigned j=0; j<size; j++)
563 // {
564 // // Check the j'th column of the Jacobian
565 // realtype yjsaved = p_y[j];
566 // realtype nudge = std::max(srur*fabs(yjsaved), min_nudge/p_ewt[j]);
567 // p_y[j] += nudge;
568 // EvaluateYDerivatives(time, y, nudge_ydot);
569 // p_y[j] = yjsaved;
570 // realtype nudge_inv = 1.0 / nudge;
571 // N_VLinearSum(nudge_inv, nudge_ydot, -nudge_inv, ydot, numeric_jth_col);
572 // realtype* p_analytic_jth_col = DENSE_COL(jacobian, j);
573 //
574 // for (unsigned i=0; i<size; i++)
575 // {
576 // if (!CompareDoubles::WithinAnyTolerance(p_numeric_jth_col[i], p_analytic_jth_col[i], rel_tol, abs_tol))
577 // {
578 // EXCEPTION("Analytic Jacobian appears dodgy at time " << time << " entry (" << i << "," << j << ").\n"
579 // << "Analytic=" << p_analytic_jth_col[i] << "; numeric=" << p_numeric_jth_col[i] << "."
580 // << DumpState("", y, time));
581 // }
582 // }
583 // }
584 //}
585 
586 
587 #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
bool HasAnalyticJacobian() const
virtual void EvaluateAnalyticJacobian(long int N, realtype time, N_Vector y, N_Vector ydot, CHASTE_CVODE_DENSE_MATRIX jacobian, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
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)