Chaste  Release::2018.1
CvodeAdaptor.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 "CvodeAdaptor.hpp"
39 
40 #include "Exception.hpp"
41 #include "MathsCustomFunctions.hpp" // For tolerance check.
42 #include "TimeStepper.hpp"
44 
45 #include <iostream>
46 #include <sstream>
47 
48 // CVODE headers
49 #include <sundials/sundials_nvector.h>
50 
51 #if CHASTE_SUNDIALS_VERSION >= 30000
52 #include <cvode/cvode_direct.h> /* access to CVDls interface */
53 #include <sundials/sundials_types.h> /* defs. of realtype, sunindextype */
54 #include <sunlinsol/sunlinsol_dense.h> /* access to dense SUNLinearSolver */
55 #include <sunmatrix/sunmatrix_dense.h> /* access to dense SUNMatrix */
56 #else
57 #include <cvode/cvode_dense.h>
58 #endif
59 
76 int CvodeRhsAdaptor(realtype t, N_Vector y, N_Vector ydot, void* pData)
77 {
78  assert(pData != nullptr);
79  CvodeData* p_data = (CvodeData*)pData;
80  // Get y, ydot into std::vector<>s
81  static std::vector<realtype> ydot_vec;
82  CopyToStdVector(y, *p_data->pY);
83  CopyToStdVector(ydot, ydot_vec);
84  // Call our function
85  try
86  {
87  p_data->pSystem->EvaluateYDerivatives(t, *(p_data->pY), ydot_vec);
88  }
89  catch (const Exception& e)
90  {
91  std::cerr << "CVODE RHS Exception: " << e.GetMessage() << std::endl
92  << std::flush;
93  return -1;
94  }
95  // Copy derivative back
96  CopyFromStdVector(ydot_vec, ydot);
97  return 0;
98 }
99 
120 int CvodeRootAdaptor(realtype t, N_Vector y, realtype* pGOut, void* pData)
121 {
122  assert(pData != nullptr);
123  CvodeData* p_data = (CvodeData*)pData;
124  // Get y into a std::vector
125  CopyToStdVector(y, *p_data->pY);
126  // Call our function
127  try
128  {
129  *pGOut = p_data->pSystem->CalculateRootFunction(t, *p_data->pY);
130  }
131  catch (const Exception& e)
132  {
133  std::cerr << "CVODE Root Exception: " << e.GetMessage() << std::endl
134  << std::flush;
135  return -1;
136  }
137  return 0;
138 }
139 
140 // /**
141 // * Jacobian computation adaptor function.
142 // *
143 // * If solving an AbstractOdeSystemWithAnalyticJacobian, this function
144 // * can be used to allow CVODE to compute the Jacobian analytically.
145 // *
146 // * Note to self: can test using pSystem->GetUseAnalyticJacobian().
147 // */
148 // int CvodeDenseJacobianAdaptor(long int numberOfStateVariables, DenseMat J,
149 // realtype t, N_Vector y, N_Vector fy,
150 // void* pData,
151 // N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
152 // {
153 // AbstractOdeSystemWithAnalyticJacobian* pSystem
154 // = (AbstractOdeSystemWithAnalyticJacobian*) pData;
155 // // Get the current time step
156 // double dt;
157 // CVodeGetCurrentStep(CvodeMem, &dt);
158 // // Get std::vector<> for y and double** for J
159 // std::vector<realtype>& y_vec = *NV_DATA_STL(y);
160 // double** ppJ = J->data; // organised column-wise: J_{i,j} = ppJ[j][i]
161 // // Call our function
162 // try
163 // {
164 // pSystem->AnalyticJacobian(y_vec, ppJ, t, dt);
165 // }
166 // catch (const Exception &e)
167 // {
168 // std::cerr << "CVODE J Exception: " << e.GetMessage() << std::endl << std::flush;
169 // return -1;
170 // }
171 // // Update J (if needed)
172 // return 0;
173 // }
174 
175 void CvodeErrorHandler(int errorCode, const char* module, const char* function,
176  char* message, void* pData)
177 {
178  std::stringstream err;
179  err << "CVODE Error " << errorCode << " in module " << module
180  << " function " << function << ": " << message;
181  std::cerr << "*" << err.str() << std::endl
182  << std::flush;
183  // Throwing the exception here causes termination on Maverick (g++ 4.4)
184  //EXCEPTION(err.str());
185 }
186 
188  std::vector<double>& rInitialY,
189  double startTime,
190  double maxStep)
191 {
192  assert(rInitialY.size() == pOdeSystem->GetNumberOfStateVariables());
193  assert(maxStep > 0.0);
194 
196  N_Vector initial_values = N_VMake_Serial(rInitialY.size(), &(rInitialY[0]));
197  assert(NV_DATA_S(initial_values) == &(rInitialY[0]));
198  assert(!NV_OWN_DATA_S(initial_values));
199  // std::cout << " Initial values: "; N_VPrint_Serial(initial_values);
200  // std::cout << " Rtol: " << mRelTol << ", Atol: " << mAbsTol << std::endl;
201  // std::cout << " Start: " << startTime << " max dt=" << maxStep << std::endl << std::flush;
202 
203  // Find out if we need to (re-)initialise
205  if (!reinit && !mForceMinimalReset)
206  {
207  const unsigned size = GetVectorSize(rInitialY);
208  for (unsigned i = 0; i < size; i++)
209  {
211  {
212  reinit = true;
213  break;
214  }
215  }
216  }
217 
218  if (!mpCvodeMem) // First run of this solver, set up CVODE memory
219  {
220  // Set up CVODE's memory.
221  mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
222  if (mpCvodeMem == nullptr) EXCEPTION("Failed to SetupCvode CVODE"); // In one line to avoid coverage problem!
223  // Set error handler
224  CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, nullptr);
225  // Set the user data
226  mData.pSystem = pOdeSystem;
227  mData.pY = &rInitialY;
228 #if CHASTE_SUNDIALS_VERSION >= 20400
229  CVodeSetUserData(mpCvodeMem, (void*)(&mData));
230 #else
231  CVodeSetFdata(mpCvodeMem, (void*)(&mData));
232 #endif
233 
234 // Setup CVODE
235 #if CHASTE_SUNDIALS_VERSION >= 20400
236  CVodeInit(mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values);
237  CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
238 #else
239  CVodeMalloc(mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values,
240  CV_SS, mRelTol, &mAbsTol);
241 #endif
242 
243  // Set the rootfinder function if wanted
244  if (mCheckForRoots)
245  {
246 #if CHASTE_SUNDIALS_VERSION >= 20400
247  CVodeRootInit(mpCvodeMem, 1, CvodeRootAdaptor);
248 #else
249  CVodeRootInit(mpCvodeMem, 1, CvodeRootAdaptor, (void*)(&mData));
250 #endif
251  }
252 
253 #if CHASTE_SUNDIALS_VERSION >= 30000
254  /* Create dense SUNMatrix for use in linear solves */
255  mpSundialsDenseMatrix = SUNDenseMatrix(rInitialY.size(), rInitialY.size());
256 
257  /* Create dense SUNLinearSolver object for use by CVode */
258  mpSundialsLinearSolver = SUNDenseLinearSolver(initial_values, mpSundialsDenseMatrix);
259 
260  /* Call CVDlsSetLinearSolver to attach the matrix and linear solver to CVode */
261  CVDlsSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
262 #else
263  // Attach a linear solver for Newton iteration
264  CVDense(mpCvodeMem, rInitialY.size());
265 #endif
266  }
267  else if (reinit) // Could be new ODE system, or new Y values
268  {
269  // Set the user data
270  mData.pSystem = pOdeSystem; // stays the same on a re-initialize
271  mData.pY = &rInitialY; // changes on a re-initialize
272 #if CHASTE_SUNDIALS_VERSION >= 20400
273  CVodeSetUserData(mpCvodeMem, (void*)(&mData));
274 #else
275  CVodeSetFdata(mpCvodeMem, (void*)(&mData));
276 #endif
277 
278 #if CHASTE_SUNDIALS_VERSION >= 20400
279  CVodeReInit(mpCvodeMem, startTime, initial_values);
280  CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
281 #else
282  CVodeReInit(mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values,
283  CV_SS, mRelTol, &mAbsTol);
284 #endif
285 
286 #if CHASTE_SUNDIALS_VERSION >= 30000
287  if (mpSundialsLinearSolver)
288  {
289  /* Free the linear solver memory */
290  SUNLinSolFree(mpSundialsLinearSolver);
291  }
292  if (mpSundialsDenseMatrix)
293  {
294  /* Free the matrix memory */
295  SUNMatDestroy(mpSundialsDenseMatrix);
296  }
297 
298  /* Create dense SUNMatrix for use in linear solves */
299  mpSundialsDenseMatrix = SUNDenseMatrix(rInitialY.size(), rInitialY.size());
300 
301  /* Create dense SUNLinearSolver object for use by CVode */
302  mpSundialsLinearSolver = SUNDenseLinearSolver(initial_values, mpSundialsDenseMatrix);
303 
304  /* Call CVDlsSetLinearSolver to attach the matrix and linear solver to CVode */
305  CVDlsSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
306 #else
307  // Attach a linear solver for Newton iteration
308  CVDense(mpCvodeMem, rInitialY.size());
309 #endif
310  }
311 
312  CVodeSetMaxStep(mpCvodeMem, maxStep);
313  // Change max steps if wanted
314  if (mMaxSteps > 0)
315  {
316  CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
317  CVodeSetMaxErrTestFails(mpCvodeMem, 15);
318  }
319  DeleteVector(initial_values);
320 }
321 
323 {
324  if (mpCvodeMem)
325  {
326  CVodeFree(&mpCvodeMem);
327  }
328  mpCvodeMem = nullptr;
329 
330 #if CHASTE_SUNDIALS_VERSION >= 30000
331  if (mpSundialsLinearSolver)
332  {
333  /* Free the linear solver memory */
334  SUNLinSolFree(mpSundialsLinearSolver);
335  }
336  mpSundialsLinearSolver = nullptr;
337 
338  if (mpSundialsDenseMatrix)
339  {
340  /* Free the matrix memory */
341  SUNMatDestroy(mpSundialsDenseMatrix);
342  }
343  mpSundialsDenseMatrix = nullptr;
344 #endif
345 }
346 
347 void CvodeAdaptor::SetForceReset(bool autoReset)
348 {
349  mForceReset = autoReset;
350  if (mForceReset)
351  {
352  SetMinimalReset(false);
353  ResetSolver();
354  }
355 }
356 
357 void CvodeAdaptor::SetMinimalReset(bool minimalReset)
358 {
359  mForceMinimalReset = minimalReset;
360  if (mForceMinimalReset)
361  {
362  SetForceReset(false);
363  }
364 }
365 
367 {
368  // Doing this makes the SetupCvode think it needs to reset things.
370 }
371 
372 void CvodeAdaptor::CvodeError(int flag, const char* msg)
373 {
374  std::stringstream err;
375  char* p_flag_name = CVodeGetReturnFlagName(flag);
376  err << msg << ": " << p_flag_name;
377  free(p_flag_name);
378  std::cerr << err.str() << std::endl
379  << std::flush;
380  EXCEPTION(err.str());
381 }
382 
384  std::vector<double>& rYValues,
385  double startTime,
386  double endTime,
387  double maxStep,
388  double timeSampling)
389 {
390  assert(endTime > startTime);
391  assert(timeSampling > 0.0);
392 
393  mStoppingEventOccurred = false;
394  if (mCheckForRoots && pOdeSystem->CalculateStoppingEvent(startTime, rYValues) == true)
395  {
396  EXCEPTION("(Solve with sampling) Stopping event is true for initial condition");
397  }
398 
399  SetupCvode(pOdeSystem, rYValues, startTime, maxStep);
400 
401  TimeStepper stepper(startTime, endTime, timeSampling);
402  N_Vector yout = N_VMake_Serial(rYValues.size(), &(rYValues[0]));
403 
404  // Set up ODE solution
405  OdeSolution solutions;
406  solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
407  solutions.rGetSolutions().push_back(rYValues);
408  solutions.rGetTimes().push_back(startTime);
409  solutions.SetOdeSystemInformation(pOdeSystem->GetSystemInformation());
410 
411  // Main time sampling loop
412  while (!stepper.IsTimeAtEnd() && !mStoppingEventOccurred)
413  {
414  // This should stop CVODE going past the end of where we wanted and interpolating back.
415  int ierr = CVodeSetStopTime(mpCvodeMem, stepper.GetNextTime());
416  assert(ierr == CV_SUCCESS);
417  UNUSED_OPT(ierr); // avoid unused var warning
418 
419  double tend;
420  ierr = CVode(mpCvodeMem, stepper.GetNextTime(), yout, &tend, CV_NORMAL);
421  if (ierr < 0)
422  {
423  FreeCvodeMemory();
424  DeleteVector(yout);
425  CvodeError(ierr, "CVODE failed to solve system");
426  }
427  // Store solution
428  solutions.rGetSolutions().push_back(rYValues);
429  solutions.rGetTimes().push_back(tend);
430  if (ierr == CV_ROOT_RETURN)
431  {
432  // Stopping event occurred
433  mStoppingEventOccurred = true;
434  mStoppingTime = tend;
435  }
436  mLastSolutionTime = tend;
437  stepper.AdvanceOneTimeStep();
438  }
439 
440  // stepper.EstimateTimeSteps may have been an overestimate...
441  solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
442 
443  int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
444  assert(ierr == CV_SUCCESS);
445  UNUSED_OPT(ierr); // avoid unused var warning
447  DeleteVector(yout);
448 
449  return solutions;
450 }
451 
453  std::vector<double>& rYValues,
454  double startTime,
455  double endTime,
456  double maxStep)
457 {
458  assert(endTime > startTime);
459 
460  mStoppingEventOccurred = false;
461  if (mCheckForRoots && pOdeSystem->CalculateStoppingEvent(startTime, rYValues) == true)
462  {
463  EXCEPTION("(Solve) Stopping event is true for initial condition");
464  }
465 
466  SetupCvode(pOdeSystem, rYValues, startTime, maxStep);
467 
468  N_Vector yout = N_VMake_Serial(rYValues.size(), &(rYValues[0]));
469 
470  // This should stop CVODE going past the end of where we wanted and interpolating back.
471  int ierr = CVodeSetStopTime(mpCvodeMem, endTime);
472  assert(ierr == CV_SUCCESS);
473  UNUSED_OPT(ierr); // avoid unused var warning
474 
475  double tend;
476  ierr = CVode(mpCvodeMem, endTime, yout, &tend, CV_NORMAL);
477  if (ierr < 0)
478  {
479  FreeCvodeMemory();
480  DeleteVector(yout);
481  CvodeError(ierr, "CVODE failed to solve system");
482  }
483  if (ierr == CV_ROOT_RETURN)
484  {
485  // Stopping event occurred
486  mStoppingEventOccurred = true;
487  mStoppingTime = tend;
488  }
489  assert(NV_DATA_S(yout) == &(rYValues[0]));
490  assert(!NV_OWN_DATA_S(yout));
491 
492  // long int steps;
493  // CVodeGetNumSteps(mpCvodeMem, &steps);
494  // std::cout << " Solved to " << endTime << " in " << steps << " steps.\n";
495 
496  ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
497  assert(ierr == CV_SUCCESS);
498  UNUSED_OPT(ierr); // avoid unused var warning
499  RecordStoppingPoint(tend, yout);
500  DeleteVector(yout);
501 }
502 
503 CvodeAdaptor::CvodeAdaptor(double relTol, double absTol)
505  mpCvodeMem(nullptr),
506  mRelTol(relTol),
507  mAbsTol(absTol),
508  mLastInternalStepSize(-0.0),
509  mMaxSteps(0),
510  mCheckForRoots(false),
511  mLastSolutionState(nullptr),
512  mLastSolutionTime(0.0),
513 #if CHASTE_SUNDIALS_VERSION >= 20400
514  mForceReset(false),
515 #else
516  mForceReset(true),
517 #endif
518  mForceMinimalReset(false)
519 #if CHASTE_SUNDIALS_VERSION >= 30000
520  ,
521  mpSundialsDenseMatrix(nullptr),
522  mpSundialsLinearSolver(nullptr)
523 #endif
524 {
525 }
526 
528 {
529  FreeCvodeMemory();
531 }
532 
533 void CvodeAdaptor::RecordStoppingPoint(double stopTime, N_Vector yEnd)
534 {
535  if (!mForceReset)
536  {
537  const unsigned size = GetVectorSize(yEnd);
539  for (unsigned i = 0; i < size; i++)
540  {
542  }
543  mLastSolutionTime = stopTime;
544  }
545 }
546 
547 void CvodeAdaptor::SetTolerances(double relTol, double absTol)
548 {
549  mRelTol = relTol;
550  mAbsTol = absTol;
551 }
552 
554 {
555  return mRelTol;
556 }
557 
559 {
560  return mAbsTol;
561 }
562 
564 {
565  return mLastInternalStepSize;
566 }
567 
569 {
570  mCheckForRoots = true;
571 }
572 
573 void CvodeAdaptor::SetMaxSteps(long int numSteps)
574 {
575  mMaxSteps = numSteps;
576 }
577 
579 {
580  return mMaxSteps;
581 }
582 
583 // Serialization for Boost >= 1.36
586 
587 #endif // CHASTE_CVODE
AbstractOdeSystem * pSystem
void SetupCvode(AbstractOdeSystem *pOdeSystem, std::vector< double > &rInitialY, double startTime, double maxStep)
static bool WithinAnyTolerance(double number1, double number2, double relTol=DBL_EPSILON, double absTol=DBL_EPSILON, bool printError=false)
CvodeAdaptor(double relTol=1e-4, double absTol=1e-6)
std::vector< std::vector< double > > & rGetSolutions()
void ResetSolver()
void CheckForStoppingEvents()
#define EXCEPTION(message)
Definition: Exception.hpp:143
N_Vector mLastSolutionState
virtual void EvaluateYDerivatives(double time, const std::vector< double > &rY, std::vector< double > &rDY)=0
void FreeCvodeMemory()
long int GetMaxSteps()
double mLastInternalStepSize
void AdvanceOneTimeStep()
bool IsTimeAtEnd() const
double GetAbsoluteTolerance()
void CreateVectorIfEmpty(VECTOR &rVec, unsigned size)
void SetOdeSystemInformation(boost::shared_ptr< const AbstractOdeSystemInformation > pOdeSystemInfo)
Definition: OdeSolution.cpp:66
void CvodeError(int flag, const char *msg)
void SetMinimalReset(bool minimalReset)
std::string GetMessage() const
Definition: Exception.cpp:83
bool mForceMinimalReset
virtual double CalculateRootFunction(double time, const std::vector< double > &rY)
void SetVectorComponent(VECTOR &rVec, unsigned index, double value)
std::vector< realtype > * pY
unsigned GetTotalTimeStepsTaken() const
boost::shared_ptr< const AbstractOdeSystemInformation > GetSystemInformation() const
std::vector< double > & rGetTimes()
unsigned EstimateTimeSteps() const
void SetTolerances(double relTol=1e-4, double absTol=1e-6)
#define UNUSED_OPT(var)
Definition: Exception.hpp:216
#define CHASTE_CLASS_EXPORT(T)
void RecordStoppingPoint(double stopTime, N_Vector yEnd)
void SetForceReset(bool autoReset)
void CopyFromStdVector(const std::vector< double > &rSrc, VECTOR &rDest)
double mLastSolutionTime
CvodeData mData
void CopyToStdVector(const VECTOR &rSrc, std::vector< double > &rDest)
long int mMaxSteps
double GetNextTime() const
double GetVectorComponent(const VECTOR &rVec, unsigned index)
void SetMaxSteps(long int numSteps)
double GetLastStepSize()
double GetRelativeTolerance()
OdeSolution Solve(AbstractOdeSystem *pOdeSystem, std::vector< double > &rYValues, double startTime, double endTime, double maxStep, double timeSampling)
void SetNumberOfTimeSteps(unsigned numTimeSteps)
Definition: OdeSolution.cpp:58
void DeleteVector(VECTOR &rVec)
unsigned GetVectorSize(const VECTOR &rVec)
virtual bool CalculateStoppingEvent(double time, const std::vector< double > &rY)