Chaste Commit::675f9facbe008c5eacb9006feaeb6423206579ea
AbstractCvodeSystem.hpp
1/*
2
3Copyright (c) 2005-2025, University of Oxford.
4All rights reserved.
5
6University of Oxford means the Chancellor, Masters and Scholars of the
7University of Oxford, having an administrative office at Wellington
8Square, Oxford OX1 2JD, UK.
9
10This file is part of Chaste.
11
12Redistribution and use in source and binary forms, with or without
13modification, 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
23THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32OF 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
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>
57#include "ClassIsAbstract.hpp"
58
59// CVODE headers
60#include <nvector/nvector_serial.h>
61
62#if CHASTE_SUNDIALS_VERSION >= 30000
63#if CHASTE_SUNDIALS_VERSION < 70000
64#include <cvode/cvode_direct.h> /* access to CVDls interface */
65#endif
66#include <sundials/sundials_types.h> /* defs. of realtype, sunindextype */
67#include <sunlinsol/sunlinsol_dense.h> /* access to dense SUNLinearSolver */
68#include <sunmatrix/sunmatrix_dense.h> /* access to dense SUNMatrix */
69#else
70#include <sundials/sundials_dense.h> /* definitions DlsMat DENSE_ELEM */
71#endif
72
73// CVODE changed their dense matrix type...
74#if CHASTE_SUNDIALS_VERSION >= 30000
75#define CHASTE_CVODE_DENSE_MATRIX SUNMatrix
76#elif CHASTE_SUNDIALS_VERSION >= 20400
77#define CHASTE_CVODE_DENSE_MATRIX DlsMat
78#else
79#define CHASTE_CVODE_DENSE_MATRIX DenseMat
80#endif
81
82// CVODE changed their way of referencing elements of a matrix. So we will copy their notation in examples.
83#if CHASTE_SUNDIALS_VERSION >= 30000
84#define IJth(A, i, j) SM_ELEMENT_D(A, i, j)
85#else
86#define IJth(A, i, j) DENSE_ELEM(A, i, j)
87#endif
88
144{
145private:
146 friend class TestAbstractCvodeSystem;
147
148 friend class boost::serialization::access;
155 template <class Archive>
156 void save(Archive& archive, const unsigned int version) const
157 {
158 // Despite the fact that 3 of these variables actually live in our base class,
159 // we still archive them here to maintain backwards compatibility,
160 // this doesn't hurt
162 archive& mUseAnalyticJacobian;
163
164 if (version >= 1u)
165 {
166 archive& mHasAnalyticJacobian;
167 }
168
169 // Convert from N_Vector to std::vector for serialization
170 const std::vector<double> state_vars = MakeStdVec(mStateVariables);
171 archive& state_vars;
172 const std::vector<double> params = MakeStdVec(mParameters);
173 archive& params;
174 archive& rGetParameterNames();
175
176 archive& mLastSolutionTime;
177 archive& mForceReset;
178 archive& mForceMinimalReset;
179 archive& mRelTol;
180 archive& mAbsTol;
181 archive& mMaxSteps;
182 archive& mLastInternalStepSize;
183
184 // We don't bother archiving CVODE's internal data, because it is missing then we'll just
185 // get a new solver being initialised after a save/load.
186
187 // This is always set up by subclass constructors, and is essentially
188 // 'static' data, so shouldn't go in the archive.
189 //archive &mpSystemInfo;
190 }
197 template <class Archive>
198 void load(Archive& archive, const unsigned int version)
199 {
201 archive& mUseAnalyticJacobian;
202
203 // This is pretty much what the code was saying before.
205 if (version >= 1u)
206 { // Overwrite if it has been archived though
207 archive& mHasAnalyticJacobian;
208 }
209
210 std::vector<double> state_vars;
211 archive& state_vars;
213
214 std::vector<double> parameters;
215 archive& parameters;
216
217 std::vector<std::string> param_names;
218 archive& param_names;
219 archive& mLastSolutionTime;
220 archive& mForceReset;
221 archive& mForceMinimalReset;
222 archive& mRelTol;
223 archive& mAbsTol;
224 archive& mMaxSteps;
225 archive& mLastInternalStepSize;
226
227 // We don't bother archiving CVODE's internal data, because it is missing then we'll just
228 // get a new solver being initialised after a save/load.
229
230 // Do some checking on the parameters
231 CheckParametersOnLoad(parameters, param_names);
232 }
233 BOOST_SERIALIZATION_SPLIT_MEMBER()
234
235
242 void SetupCvode(N_Vector initialConditions,
243 realtype tStart,
244 realtype maxDt);
245
250 void RecordStoppingPoint(double stopTime);
251
253 void FreeCvodeMemory();
254
265 void CvodeError(int flag, const char* msg, const double& rTime,
266 const double& rStartTime, const double& rEndTime);
267
270
273
276
279
280#if CHASTE_SUNDIALS_VERSION >= 30000
282 SUNMatrix mpSundialsDenseMatrix;
284 SUNLinearSolver mpSundialsLinearSolver;
285#endif
286
287protected:
290
293
295 double mRelTol;
296
298 double mAbsTol;
299
302
307 long int mMaxSteps;
308
311
316 void Init();
317
318public:
324 AbstractCvodeSystem(unsigned numberOfStateVariables);
325
329 virtual ~AbstractCvodeSystem();
330
338 virtual void EvaluateYDerivatives(realtype time,
339 const N_Vector y,
340 N_Vector ydot)
341 = 0;
342
358 virtual void EvaluateAnalyticJacobian(realtype time, N_Vector y, N_Vector ydot,
359 CHASTE_CVODE_DENSE_MATRIX jacobian,
360 N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
361 {
362 EXCEPTION("No analytic Jacobian has been defined for this system.");
363 }
364
375 void SetForceReset(bool autoReset);
376
385 void SetMinimalReset(bool minimalReset);
386
393 bool GetMinimalReset();
394
401 bool GetForceReset();
402
413 void ResetSolver();
414
434 OdeSolution Solve(realtype tStart,
435 realtype tEnd,
436 realtype maxDt,
437 realtype tSamp);
438
454 void Solve(realtype tStart,
455 realtype tEnd,
456 realtype maxDt);
457
464 void SetMaxSteps(long int numSteps);
465
470 long int GetMaxSteps();
471
479 void SetTolerances(double relTol = 1e-5, double absTol = 1e-7);
480
484 double GetRelativeTolerance();
485
489 double GetAbsoluteTolerance();
490
494 double GetLastStepSize();
495
496 /* NB This needs making into a doxygen comment if you bring the method back in.
497 *
498 * An alternative approach to stopping events; currently only useful with CVODE.
499 * CVODE can search for roots (zeros) of this function while solving the ODE system,
500 * and home in on them to find sign transitions to high precision.
501 *
502 * The default implementation here fakes a root function using CalculateStoppingEvent.
503 *
504 * @param time the current time
505 * @param rY the current values of the state variables
506 */
507 // virtual double CalculateRootFunction(double time, const std::vector<double>& rY);
508
512 bool GetUseAnalyticJacobian() const;
513
517 bool HasAnalyticJacobian() const;
518
525 void ForceUseOfNumericalJacobian(bool useNumericalJacobian = true);
526
527 // The following method may be useful to identify problems with the Analytic Jacobians, if anything goes wrong,
528 // but #1795 seems to have got these working OK, so commented out for now.
529
530 // /*
531 // * Compare the calculated analytic jacobian to a numerical approximation, and throw if it looks silly.
532 // *
533 // * @param time the current time
534 // * @param y the current state variables
535 // * @param jacobian the analytic jacobian matrix
536 // * @param tmp1 working memory of the correct size provided by CVODE for temporary calculations
537 // * @param tmp2 working memory of the correct size provided by CVODE for temporary calculations
538 // * @param tmp3 working memory of the correct size provided by CVODE for temporary calculations
539 // */
540 // void CheckAnalyticJacobian(realtype time, N_Vector y, N_Vector ydot,
541 // CHASTE_CVODE_DENSE_MATRIX jacobian,
542 // N_Vector tmp1, N_Vector tmp2, N_Vector tmp3);
543};
544
546BOOST_CLASS_VERSION(AbstractCvodeSystem, 1u)
547
548#endif //_ABSTRACTCVODESYSTEM_HPP_
549#endif // CHASTE_CVODE
#define CLASS_IS_ABSTRACT(T)
#define EXCEPTION(message)
std::vector< double > MakeStdVec(N_Vector v)
void CopyFromStdVector(const std::vector< double > &rSrc, VECTOR &rDest)
void SetTolerances(double relTol=1e-5, double absTol=1e-7)
bool GetMinimalReset()
Get whether we want to run with minimal reset or not (no reinitialisation of the solver if variables ...
void SetForceReset(bool autoReset)
void load(Archive &archive, const unsigned int version)
void SetupCvode(N_Vector initialConditions, realtype tStart, realtype maxDt)
void CvodeError(int flag, const char *msg, const double &rTime, const double &rStartTime, const double &rEndTime)
void RecordStoppingPoint(double stopTime)
void SetMaxSteps(long int numSteps)
OdeSolution Solve(realtype tStart, realtype tEnd, realtype maxDt, realtype tSamp)
void ForceUseOfNumericalJacobian(bool useNumericalJacobian=true)
virtual void EvaluateYDerivatives(realtype time, const N_Vector y, N_Vector ydot)=0
void save(Archive &archive, const unsigned int version) const
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)
void SetMinimalReset(bool minimalReset)
bool GetForceReset()
Get whether we will force a solver reset on every call to Solve()
void CheckParametersOnLoad(const std::vector< double > &rParameters, const std::vector< std::string > &rParameterNames)
const std::vector< std::string > & rGetParameterNames() const