Chaste  Release::2018.1
AbstractCvodeSystem.hpp
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 #ifndef _ABSTRACTCVODESYSTEM_HPP_
38 #define _ABSTRACTCVODESYSTEM_HPP_
39 
40 #include <algorithm>
41 #include <string>
42 #include <vector>
43 
44 // This is only needed to prevent compilation errors on PETSc 2.2/Boost 1.33.1 combo
45 #include "UblasVectorInclude.hpp"
46 
47 // Chaste includes
48 #include "AbstractParameterisedSystem.hpp"
49 #include "Exception.hpp"
50 #include "OdeSolution.hpp"
52 
53 // Serialiazation
54 #include <boost/serialization/split_member.hpp>
55 #include <boost/serialization/vector.hpp>
56 #include "ChasteSerialization.hpp"
57 #include "ClassIsAbstract.hpp"
58 
59 // CVODE headers
60 #include <nvector/nvector_serial.h>
61 
62 #if CHASTE_SUNDIALS_VERSION >= 30000
63 #include <cvode/cvode_direct.h> /* access to CVDls interface */
64 #include <sundials/sundials_types.h> /* defs. of realtype, sunindextype */
65 #include <sunlinsol/sunlinsol_dense.h> /* access to dense SUNLinearSolver */
66 #include <sunmatrix/sunmatrix_dense.h> /* access to dense SUNMatrix */
67 #else
68 #include <sundials/sundials_dense.h> /* definitions DlsMat DENSE_ELEM */
69 #endif
70 
71 // CVODE changed their dense matrix type...
72 #if CHASTE_SUNDIALS_VERSION >= 30000
73 #define CHASTE_CVODE_DENSE_MATRIX SUNMatrix
74 #elif CHASTE_SUNDIALS_VERSION >= 20400
75 #define CHASTE_CVODE_DENSE_MATRIX DlsMat
76 #else
77 #define CHASTE_CVODE_DENSE_MATRIX DenseMat
78 #endif
79 
80 // CVODE changed their way of referencing elements of a matrix. So we will copy their notation in examples.
81 #if CHASTE_SUNDIALS_VERSION >= 30000
82 #define IJth(A, i, j) SM_ELEMENT_D(A, i, j)
83 #else
84 #define IJth(A, i, j) DENSE_ELEM(A, i, j)
85 #endif
86 
143 {
144 private:
145  friend class TestAbstractCvodeSystem;
146 
147  friend class boost::serialization::access;
154  template <class Archive>
155  void save(Archive& archive, const unsigned int version) const
156  {
157  // Despite the fact that 3 of these variables actually live in our base class,
158  // we still archive them here to maintain backwards compatibility,
159  // this doesn't hurt
160  archive& mNumberOfStateVariables;
161  archive& mUseAnalyticJacobian;
162 
163  if (version >= 1u)
164  {
165  archive& mHasAnalyticJacobian;
166  }
167 
168  // Convert from N_Vector to std::vector for serialization
169  const std::vector<double> state_vars = MakeStdVec(mStateVariables);
170  archive& state_vars;
171  const std::vector<double> params = MakeStdVec(mParameters);
172  archive& params;
173  archive& rGetParameterNames();
174 
175  archive& mLastSolutionTime;
176  archive& mForceReset;
177  archive& mForceMinimalReset;
178  archive& mRelTol;
179  archive& mAbsTol;
180  archive& mMaxSteps;
181  archive& mLastInternalStepSize;
182 
183  // We don't bother archiving CVODE's internal data, because it is missing then we'll just
184  // get a new solver being initialised after a save/load.
185 
186  // This is always set up by subclass constructors, and is essentially
187  // 'static' data, so shouldn't go in the archive.
188  //archive &mpSystemInfo;
189  }
196  template <class Archive>
197  void load(Archive& archive, const unsigned int version)
198  {
199  archive& mNumberOfStateVariables;
200  archive& mUseAnalyticJacobian;
201 
202  // This is pretty much what the code was saying before.
204  if (version >= 1u)
205  { // Overwrite if it has been archived though
206  archive& mHasAnalyticJacobian;
207  }
208 
209  std::vector<double> state_vars;
210  archive& state_vars;
211  CopyFromStdVector(state_vars, mStateVariables);
212 
213  std::vector<double> parameters;
214  archive& parameters;
215 
216  std::vector<std::string> param_names;
217  archive& param_names;
218  archive& mLastSolutionTime;
219  archive& mForceReset;
220  archive& mForceMinimalReset;
221  archive& mRelTol;
222  archive& mAbsTol;
223  archive& mMaxSteps;
224  archive& mLastInternalStepSize;
225 
226  // We don't bother archiving CVODE's internal data, because it is missing then we'll just
227  // get a new solver being initialised after a save/load.
228 
229  // Do some checking on the parameters
230  CheckParametersOnLoad(parameters, param_names);
231  }
232  BOOST_SERIALIZATION_SPLIT_MEMBER()
233 
234 
241  void SetupCvode(N_Vector initialConditions,
242  realtype tStart,
243  realtype maxDt);
244 
249  void RecordStoppingPoint(double stopTime);
250 
252  void FreeCvodeMemory();
253 
264  void CvodeError(int flag, const char* msg, const double& rTime,
265  const double& rStartTime, const double& rEndTime);
266 
269 
272 
275 
278 
279 #if CHASTE_SUNDIALS_VERSION >= 30000
280 
281  SUNMatrix mpSundialsDenseMatrix;
283  SUNLinearSolver mpSundialsLinearSolver;
284 #endif
285 
286 protected:
289 
292 
294  double mRelTol;
295 
297  double mAbsTol;
298 
300  void* mpCvodeMem;
301 
306  long int mMaxSteps;
307 
310 
315  void Init();
316 
317 public:
323  AbstractCvodeSystem(unsigned numberOfStateVariables);
324 
328  virtual ~AbstractCvodeSystem();
329 
337  virtual void EvaluateYDerivatives(realtype time,
338  const N_Vector y,
339  N_Vector ydot)
340  = 0;
341 
357  virtual void EvaluateAnalyticJacobian(realtype time, N_Vector y, N_Vector ydot,
358  CHASTE_CVODE_DENSE_MATRIX jacobian,
359  N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
360  {
361  EXCEPTION("No analytic Jacobian has been defined for this system.");
362  }
363 
374  void SetForceReset(bool autoReset);
375 
384  void SetMinimalReset(bool minimalReset);
385 
396  void ResetSolver();
397 
417  OdeSolution Solve(realtype tStart,
418  realtype tEnd,
419  realtype maxDt,
420  realtype tSamp);
421 
437  void Solve(realtype tStart,
438  realtype tEnd,
439  realtype maxDt);
440 
447  void SetMaxSteps(long int numSteps);
448 
453  long int GetMaxSteps();
454 
462  void SetTolerances(double relTol = 1e-5, double absTol = 1e-7);
463 
467  double GetRelativeTolerance();
468 
472  double GetAbsoluteTolerance();
473 
477  double GetLastStepSize();
478 
479  /* NB This needs making into a doxygen comment if you bring the method back in.
480  *
481  * An alternative approach to stopping events; currently only useful with CVODE.
482  * CVODE can search for roots (zeros) of this function while solving the ODE system,
483  * and home in on them to find sign transitions to high precision.
484  *
485  * The default implementation here fakes a root function using CalculateStoppingEvent.
486  *
487  * @param time the current time
488  * @param rY the current values of the state variables
489  */
490  // virtual double CalculateRootFunction(double time, const std::vector<double>& rY);
491 
495  bool GetUseAnalyticJacobian() const;
496 
500  bool HasAnalyticJacobian() const;
501 
508  void ForceUseOfNumericalJacobian(bool useNumericalJacobian = true);
509 
510  // The following method may be useful to identify problems with the Analytic Jacobians, if anything goes wrong,
511  // but #1795 seems to have got these working OK, so commented out for now.
512 
513  // /*
514  // * Compare the calculated analytic jacobian to a numerical approximation, and throw if it looks silly.
515  // *
516  // * @param time the current time
517  // * @param y the current state variables
518  // * @param jacobian the analytic jacobian matrix
519  // * @param tmp1 working memory of the correct size provided by CVODE for temporary calculations
520  // * @param tmp2 working memory of the correct size provided by CVODE for temporary calculations
521  // * @param tmp3 working memory of the correct size provided by CVODE for temporary calculations
522  // */
523  // void CheckAnalyticJacobian(realtype time, N_Vector y, N_Vector ydot,
524  // CHASTE_CVODE_DENSE_MATRIX jacobian,
525  // N_Vector tmp1, N_Vector tmp2, N_Vector tmp3);
526 };
527 
529 BOOST_CLASS_VERSION(AbstractCvodeSystem, 1u)
530 
531 #endif //_ABSTRACTCVODESYSTEM_HPP_
532 #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)
bool GetUseAnalyticJacobian() const
#define CLASS_IS_ABSTRACT(T)
void RecordStoppingPoint(double stopTime)
#define EXCEPTION(message)
Definition: Exception.hpp:143
void SetupCvode(N_Vector initialConditions, realtype tStart, realtype maxDt)
void load(Archive &archive, const unsigned int version)
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
void SetMaxSteps(long int numSteps)
OdeSolution Solve(realtype tStart, realtype tEnd, realtype maxDt, realtype tSamp)
void CopyFromStdVector(const std::vector< double > &rSrc, VECTOR &rDest)
void CheckParametersOnLoad(const std::vector< double > &rParameters, const std::vector< std::string > &rParameterNames)
void ForceUseOfNumericalJacobian(bool useNumericalJacobian=true)
void save(Archive &archive, const unsigned int version) const
AbstractCvodeSystem(unsigned numberOfStateVariables)
const std::vector< std::string > & rGetParameterNames() const