Chaste  Release::2024.1
CvodeAdaptor.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 "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 #if CHASTE_SUNDIALS_VERSION >= 40000
222  // v4.0.0 release notes: instead of specifying the nonlinear iteration type when creating the CVODE(S) memory structure,
223  // CVODE(S) uses the SUNNONLINSOL_NEWTON module implementation of a Newton iteration by default.
224  mpCvodeMem = CVodeCreate(CV_BDF);
225 #else
226  mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
227 #endif
228 
229  if (mpCvodeMem == nullptr)
230  EXCEPTION("Failed to SetupCvode CVODE"); // LCOV_EXCL_LINE
231  // Set error handler
232  CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, nullptr);
233  // Set the user data
234  mData.pSystem = pOdeSystem;
235  mData.pY = &rInitialY;
236 #if CHASTE_SUNDIALS_VERSION >= 20400
237  CVodeSetUserData(mpCvodeMem, (void*)(&mData));
238 #else
239  CVodeSetFdata(mpCvodeMem, (void*)(&mData));
240 #endif
241 
242 // Setup CVODE
243 #if CHASTE_SUNDIALS_VERSION >= 20400
244  CVodeInit(mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values);
245  CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
246 #else
247  CVodeMalloc(mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values,
248  CV_SS, mRelTol, &mAbsTol);
249 #endif
250 
251  // Set the rootfinder function if wanted
252  if (mCheckForRoots)
253  {
254 #if CHASTE_SUNDIALS_VERSION >= 20400
255  CVodeRootInit(mpCvodeMem, 1, CvodeRootAdaptor);
256 #else
257  CVodeRootInit(mpCvodeMem, 1, CvodeRootAdaptor, (void*)(&mData));
258 #endif
259  }
260 
261 #if CHASTE_SUNDIALS_VERSION >= 30000
262  /* Create dense SUNMatrix for use in linear solves */
263  mpSundialsDenseMatrix = SUNDenseMatrix(rInitialY.size(), rInitialY.size());
264 #endif
265 
266 #if CHASTE_SUNDIALS_VERSION >= 40000
267  /* Create dense SUNLinearSolver object for use by CVode */
268  mpSundialsLinearSolver = SUNLinSol_Dense(initial_values, mpSundialsDenseMatrix);
269 
270  /* Call CVodeSetLinearSolver to attach the matrix and linear solver to CVode */
271  CVodeSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
272 #elif CHASTE_SUNDIALS_VERSION >= 30000
273  /* Create dense SUNLinearSolver object for use by CVode */
274  mpSundialsLinearSolver = SUNDenseLinearSolver(initial_values, mpSundialsDenseMatrix);
275 
276  /* Call CVDlsSetLinearSolver to attach the matrix and linear solver to CVode */
277  CVDlsSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
278 #else
279  // CVODE < v3.0.0
280  // Attach a linear solver for Newton iteration
281  CVDense(mpCvodeMem, rInitialY.size());
282 #endif
283  }
284  else if (reinit) // Could be new ODE system, or new Y values
285  {
286  // Set the user data
287  mData.pSystem = pOdeSystem; // stays the same on a re-initialize
288  mData.pY = &rInitialY; // changes on a re-initialize
289 #if CHASTE_SUNDIALS_VERSION >= 20400
290  CVodeSetUserData(mpCvodeMem, (void*)(&mData));
291 #else
292  CVodeSetFdata(mpCvodeMem, (void*)(&mData));
293 #endif
294 
295 #if CHASTE_SUNDIALS_VERSION >= 20400
296  CVodeReInit(mpCvodeMem, startTime, initial_values);
297  CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
298 #else
299  CVodeReInit(mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values,
300  CV_SS, mRelTol, &mAbsTol);
301 #endif
302 
303 #if CHASTE_SUNDIALS_VERSION >= 30000
304  if (mpSundialsLinearSolver)
305  {
306  /* Free the linear solver memory */
307  SUNLinSolFree(mpSundialsLinearSolver);
308  }
309  if (mpSundialsDenseMatrix)
310  {
311  /* Free the matrix memory */
312  SUNMatDestroy(mpSundialsDenseMatrix);
313  }
314 
315  /* Create dense matrix of type SUNDenseMatrix for use in linear solves */
316  mpSundialsDenseMatrix = SUNDenseMatrix(rInitialY.size(), rInitialY.size());
317 
318 #if CHASTE_SUNDIALS_VERSION >= 40000
319  /* Create dense SUNLinSol_Dense object for use by CVode */
320  mpSundialsLinearSolver = SUNLinSol_Dense(initial_values, mpSundialsDenseMatrix);
321 
322  /* Call CVodeSetLinearSolver to attach the matrix and linear solver to CVode */
323  CVodeSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
324 #else
325  /* Create dense SUNLinearSolver object for use by CVode */
326  mpSundialsLinearSolver = SUNDenseLinearSolver(initial_values, mpSundialsDenseMatrix);
327 
328  /* Call CVDlsSetLinearSolver to attach the matrix and linear solver to CVode */
329  CVDlsSetLinearSolver(mpCvodeMem, mpSundialsLinearSolver, mpSundialsDenseMatrix);
330 #endif
331 #else
332  // Attach a linear solver for Newton iteration
333  CVDense(mpCvodeMem, rInitialY.size());
334 #endif
335  }
336 
337  CVodeSetMaxStep(mpCvodeMem, maxStep);
338  // Change max steps if wanted
339  if (mMaxSteps > 0)
340  {
341  CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
342  CVodeSetMaxErrTestFails(mpCvodeMem, 15);
343  }
344  DeleteVector(initial_values);
345 }
346 
348 {
349  if (mpCvodeMem)
350  {
351  CVodeFree(&mpCvodeMem);
352  }
353  mpCvodeMem = nullptr;
354 
355 #if CHASTE_SUNDIALS_VERSION >= 30000
356  if (mpSundialsLinearSolver)
357  {
358  /* Free the linear solver memory */
359  SUNLinSolFree(mpSundialsLinearSolver);
360  }
361  mpSundialsLinearSolver = nullptr;
362 
363  if (mpSundialsDenseMatrix)
364  {
365  /* Free the matrix memory */
366  SUNMatDestroy(mpSundialsDenseMatrix);
367  }
368  mpSundialsDenseMatrix = nullptr;
369 #endif
370 }
371 
372 void CvodeAdaptor::SetForceReset(bool autoReset)
373 {
374  mForceReset = autoReset;
375  if (mForceReset)
376  {
377  SetMinimalReset(false);
378  ResetSolver();
379  }
380 }
381 
382 void CvodeAdaptor::SetMinimalReset(bool minimalReset)
383 {
384  mForceMinimalReset = minimalReset;
385  if (mForceMinimalReset)
386  {
387  SetForceReset(false);
388  }
389 }
390 
392 {
393  // Doing this makes the SetupCvode think it needs to reset things.
395 }
396 
397 void CvodeAdaptor::CvodeError(int flag, const char* msg)
398 {
399  std::stringstream err;
400  char* p_flag_name = CVodeGetReturnFlagName(flag);
401  err << msg << ": " << p_flag_name;
402  free(p_flag_name);
403  std::cerr << err.str() << std::endl
404  << std::flush;
405  EXCEPTION(err.str());
406 }
407 
409  std::vector<double>& rYValues,
410  double startTime,
411  double endTime,
412  double maxStep,
413  double timeSampling)
414 {
415  assert(endTime > startTime);
416  assert(timeSampling > 0.0);
417 
418  mStoppingEventOccurred = false;
419  if (mCheckForRoots && pOdeSystem->CalculateStoppingEvent(startTime, rYValues) == true)
420  {
421  EXCEPTION("(Solve with sampling) Stopping event is true for initial condition");
422  }
423 
424  SetupCvode(pOdeSystem, rYValues, startTime, maxStep);
425 
426  TimeStepper stepper(startTime, endTime, timeSampling);
427  N_Vector yout = N_VMake_Serial(rYValues.size(), &(rYValues[0]));
428 
429  // Set up ODE solution
430  OdeSolution solutions;
431  solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
432  solutions.rGetSolutions().push_back(rYValues);
433  solutions.rGetTimes().push_back(startTime);
434  solutions.SetOdeSystemInformation(pOdeSystem->GetSystemInformation());
435 
436  // Main time sampling loop
437  while (!stepper.IsTimeAtEnd() && !mStoppingEventOccurred)
438  {
439  // This should stop CVODE going past the end of where we wanted and interpolating back.
440  int ierr = CVodeSetStopTime(mpCvodeMem, stepper.GetNextTime());
441  assert(ierr == CV_SUCCESS);
442  UNUSED_OPT(ierr); // avoid unused var warning
443 
444  double tend;
445  ierr = CVode(mpCvodeMem, stepper.GetNextTime(), yout, &tend, CV_NORMAL);
446  if (ierr < 0)
447  {
448  FreeCvodeMemory();
449  DeleteVector(yout);
450  CvodeError(ierr, "CVODE failed to solve system");
451  }
452  // Store solution
453  solutions.rGetSolutions().push_back(rYValues);
454  solutions.rGetTimes().push_back(tend);
455  if (ierr == CV_ROOT_RETURN)
456  {
457  // Stopping event occurred
458  mStoppingEventOccurred = true;
459  mStoppingTime = tend;
460  }
461  mLastSolutionTime = tend;
462  stepper.AdvanceOneTimeStep();
463  }
464 
465  // stepper.EstimateTimeSteps may have been an overestimate...
466  solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
467 
468  int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
469  assert(ierr == CV_SUCCESS);
470  UNUSED_OPT(ierr); // avoid unused var warning
472  DeleteVector(yout);
473 
474  return solutions;
475 }
476 
478  std::vector<double>& rYValues,
479  double startTime,
480  double endTime,
481  double maxStep)
482 {
483  assert(endTime > startTime);
484 
485  mStoppingEventOccurred = false;
486  if (mCheckForRoots && pOdeSystem->CalculateStoppingEvent(startTime, rYValues) == true)
487  {
488  EXCEPTION("(Solve) Stopping event is true for initial condition");
489  }
490 
491  SetupCvode(pOdeSystem, rYValues, startTime, maxStep);
492 
493  N_Vector yout = N_VMake_Serial(rYValues.size(), &(rYValues[0]));
494 
495  // This should stop CVODE going past the end of where we wanted and interpolating back.
496  int ierr = CVodeSetStopTime(mpCvodeMem, endTime);
497  assert(ierr == CV_SUCCESS);
498  UNUSED_OPT(ierr); // avoid unused var warning
499 
500  double tend;
501  ierr = CVode(mpCvodeMem, endTime, yout, &tend, CV_NORMAL);
502  if (ierr < 0)
503  {
504  FreeCvodeMemory();
505  DeleteVector(yout);
506  CvodeError(ierr, "CVODE failed to solve system");
507  }
508  if (ierr == CV_ROOT_RETURN)
509  {
510  // Stopping event occurred
511  mStoppingEventOccurred = true;
512  mStoppingTime = tend;
513  }
514  assert(NV_DATA_S(yout) == &(rYValues[0]));
515  assert(!NV_OWN_DATA_S(yout));
516 
517  // long int steps;
518  // CVodeGetNumSteps(mpCvodeMem, &steps);
519  // std::cout << " Solved to " << endTime << " in " << steps << " steps.\n";
520 
521  ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
522  assert(ierr == CV_SUCCESS);
523  UNUSED_OPT(ierr); // avoid unused var warning
524  RecordStoppingPoint(tend, yout);
525  DeleteVector(yout);
526 }
527 
528 CvodeAdaptor::CvodeAdaptor(double relTol, double absTol)
530  mpCvodeMem(nullptr),
531  mRelTol(relTol),
532  mAbsTol(absTol),
533  mLastInternalStepSize(-0.0),
534  mMaxSteps(0),
535  mCheckForRoots(false),
536  mLastSolutionState(nullptr),
537  mLastSolutionTime(0.0),
538 #if CHASTE_SUNDIALS_VERSION >= 20400
539  mForceReset(false),
540 #else
541  mForceReset(true),
542 #endif
543  mForceMinimalReset(false)
544 #if CHASTE_SUNDIALS_VERSION >= 30000
545  ,
546  mpSundialsDenseMatrix(nullptr),
547  mpSundialsLinearSolver(nullptr)
548 #endif
549 {
550 }
551 
553 {
554  FreeCvodeMemory();
556 }
557 
558 void CvodeAdaptor::RecordStoppingPoint(double stopTime, N_Vector yEnd)
559 {
560  if (!mForceReset)
561  {
562  const unsigned size = GetVectorSize(yEnd);
564  for (unsigned i = 0; i < size; i++)
565  {
567  }
568  mLastSolutionTime = stopTime;
569  }
570 }
571 
572 void CvodeAdaptor::SetTolerances(double relTol, double absTol)
573 {
574  mRelTol = relTol;
575  mAbsTol = absTol;
576 }
577 
579 {
580  return mRelTol;
581 }
582 
584 {
585  return mAbsTol;
586 }
587 
589 {
590  return mLastInternalStepSize;
591 }
592 
594 {
595  mCheckForRoots = true;
596 }
597 
598 void CvodeAdaptor::SetMaxSteps(long int numSteps)
599 {
600  mMaxSteps = numSteps;
601 }
602 
604 {
605  return mMaxSteps;
606 }
607 
608 // Serialization for Boost >= 1.36
611 
612 #endif // CHASTE_CVODE
bool IsTimeAtEnd() const
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()
boost::shared_ptr< const AbstractOdeSystemInformation > GetSystemInformation() const
#define EXCEPTION(message)
Definition: Exception.hpp:143
std::string GetMessage() const
Definition: Exception.cpp:83
N_Vector mLastSolutionState
double GetNextTime() const
virtual void EvaluateYDerivatives(double time, const std::vector< double > &rY, std::vector< double > &rDY)=0
void FreeCvodeMemory()
unsigned GetTotalTimeStepsTaken() const
long int GetMaxSteps()
double mLastInternalStepSize
void AdvanceOneTimeStep()
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)
bool mForceMinimalReset
virtual double CalculateRootFunction(double time, const std::vector< double > &rY)
unsigned EstimateTimeSteps() const
void SetVectorComponent(VECTOR &rVec, unsigned index, double value)
std::vector< realtype > * pY
std::vector< double > & rGetTimes()
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 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)