Chaste  Release::2017.1
AbstractCvodeSystem.hpp
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 #ifndef _ABSTRACTCVODESYSTEM_HPP_
38 #define _ABSTRACTCVODESYSTEM_HPP_
39 
40 #include <vector>
41 #include <string>
42 #include <algorithm>
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 "OdeSolution.hpp"
49 #include "AbstractParameterisedSystem.hpp"
50 #include "Exception.hpp"
52 
53 // Serialiazation
54 #include "ChasteSerialization.hpp"
55 #include <boost/serialization/split_member.hpp>
56 #include <boost/serialization/vector.hpp>
57 #include "ClassIsAbstract.hpp"
58 
59 // CVODE headers
60 #include <nvector/nvector_serial.h>
61 #include <sundials/sundials_dense.h> /* definitions DlsMat DENSE_ELEM */
62 
63 // CVODE changed their dense matrix type...
64 #if CHASTE_SUNDIALS_VERSION >= 20400
65 #define CHASTE_CVODE_DENSE_MATRIX DlsMat
66 #else
67 #define CHASTE_CVODE_DENSE_MATRIX DenseMat
68 #endif
69 
120 {
121 private:
122  friend class TestAbstractCvodeSystem;
123 
124  friend class boost::serialization::access;
131  template<class Archive>
132  void save(Archive & archive, const unsigned int version) const
133  {
134  // Despite the fact that 3 of these variables actually live in our base class,
135  // we still archive them here to maintain backwards compatibility,
136  // this doesn't hurt
137  archive & mNumberOfStateVariables;
138  archive & mUseAnalyticJacobian;
139 
140  if (version >= 1u)
141  {
142  archive & mHasAnalyticJacobian;
143  }
144 
145  // Convert from N_Vector to std::vector for serialization
146  const std::vector<double> state_vars = MakeStdVec(mStateVariables);
147  archive & state_vars;
148  const std::vector<double> params = MakeStdVec(mParameters);
149  archive & params;
150  archive & rGetParameterNames();
151 
152  archive & mLastSolutionTime;
153  archive & mForceReset;
154  archive & mForceMinimalReset;
155  archive & mRelTol;
156  archive & mAbsTol;
157  archive & mMaxSteps;
158  archive & mLastInternalStepSize;
159 
160  // We don't bother archiving CVODE's internal data, because it is missing then we'll just
161  // get a new solver being initialised after a save/load.
162 
163 
164  // This is always set up by subclass constructors, and is essentially
165  // 'static' data, so shouldn't go in the archive.
166  //archive &mpSystemInfo;
167  }
174  template<class Archive>
175  void load(Archive & archive, const unsigned int version)
176  {
177  archive & mNumberOfStateVariables;
178  archive & mUseAnalyticJacobian;
179 
180  // This is pretty much what the code was saying before.
182  if (version >= 1u)
183  { // Overwrite if it has been archived though
184  archive & mHasAnalyticJacobian;
185  }
186 
187  std::vector<double> state_vars;
188  archive & state_vars;
190 
191  std::vector<double> parameters;
192  archive & parameters;
193 
194  std::vector<std::string> param_names;
195  archive & param_names;
196  archive & mLastSolutionTime;
197  archive & mForceReset;
198  archive & mForceMinimalReset;
199  archive & mRelTol;
200  archive & mAbsTol;
201  archive & mMaxSteps;
202  archive & mLastInternalStepSize;
203 
204  // We don't bother archiving CVODE's internal data, because it is missing then we'll just
205  // get a new solver being initialised after a save/load.
206 
207  // Do some checking on the parameters
208  CheckParametersOnLoad(parameters,param_names);
209  }
210  BOOST_SERIALIZATION_SPLIT_MEMBER()
211 
212 
219  void SetupCvode(N_Vector initialConditions,
220  realtype tStart,
221  realtype maxDt);
222 
227  void RecordStoppingPoint(double stopTime);
228 
230  void FreeCvodeMemory();
231 
242  void CvodeError(int flag, const char * msg, const double& rTime,
243  const double& rStartTime, const double& rEndTime);
244 
247 
250 
253 
256 
257 protected:
258 
261 
264 
266  double mRelTol;
267 
269  double mAbsTol;
270 
272  void* mpCvodeMem;
273 
278  long int mMaxSteps;
279 
282 
287  void Init();
288 
289 public:
290 
296  AbstractCvodeSystem(unsigned numberOfStateVariables);
297 
301  virtual ~AbstractCvodeSystem();
302 
310  virtual void EvaluateYDerivatives(realtype time,
311  const N_Vector y,
312  N_Vector ydot)=0;
313 
314 
331  virtual void EvaluateAnalyticJacobian(long int N, realtype time, N_Vector y, N_Vector ydot,
332  CHASTE_CVODE_DENSE_MATRIX jacobian,
333  N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
334  {
335  EXCEPTION("No analytic Jacobian has been defined for this system.");
336  }
337 
348  void SetForceReset(bool autoReset);
349 
358  void SetMinimalReset(bool minimalReset);
359 
370  void ResetSolver();
371 
391  OdeSolution Solve(realtype tStart,
392  realtype tEnd,
393  realtype maxDt,
394  realtype tSamp);
395 
411  void Solve(realtype tStart,
412  realtype tEnd,
413  realtype maxDt);
414 
421  void SetMaxSteps(long int numSteps);
422 
427  long int GetMaxSteps();
428 
436  void SetTolerances(double relTol=1e-5, double absTol=1e-7);
437 
441  double GetRelativeTolerance();
442 
446  double GetAbsoluteTolerance();
447 
451  double GetLastStepSize();
452 
453 
454  /* NB This needs making into a doxygen comment if you bring the method back in.
455  *
456  * An alternative approach to stopping events; currently only useful with CVODE.
457  * CVODE can search for roots (zeros) of this function while solving the ODE system,
458  * and home in on them to find sign transitions to high precision.
459  *
460  * The default implementation here fakes a root function using CalculateStoppingEvent.
461  *
462  * @param time the current time
463  * @param rY the current values of the state variables
464  */
465 // virtual double CalculateRootFunction(double time, const std::vector<double>& rY);
466 
470  bool GetUseAnalyticJacobian() const;
471 
475  bool HasAnalyticJacobian() const;
476 
483  void ForceUseOfNumericalJacobian(bool useNumericalJacobian = true);
484 
485 // The following method may be useful to identify problems with the Analytic Jacobians, if anything goes wrong,
486 // but #1795 seems to have got these working OK, so commented out for now.
487 
488 // /*
489 // * Compare the calculated analytic jacobian to a numerical approximation, and throw if it looks silly.
490 // *
491 // * @param time the current time
492 // * @param y the current state variables
493 // * @param jacobian the analytic jacobian matrix
494 // * @param tmp1 working memory of the correct size provided by CVODE for temporary calculations
495 // * @param tmp2 working memory of the correct size provided by CVODE for temporary calculations
496 // * @param tmp3 working memory of the correct size provided by CVODE for temporary calculations
497 // */
498 // void CheckAnalyticJacobian(realtype time, N_Vector y, N_Vector ydot,
499 // CHASTE_CVODE_DENSE_MATRIX jacobian,
500 // N_Vector tmp1, N_Vector tmp2, N_Vector tmp3);
501 };
502 
504 BOOST_CLASS_VERSION(AbstractCvodeSystem, 1u)
505 
506 #endif //_ABSTRACTCVODESYSTEM_HPP_
507 #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)
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
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