37#ifndef _ABSTRACTCVODESYSTEM_HPP_
38#define _ABSTRACTCVODESYSTEM_HPP_
48#include "AbstractParameterisedSystem.hpp"
50#include "OdeSolution.hpp"
54#include <boost/serialization/split_member.hpp>
55#include <boost/serialization/vector.hpp>
60#include <nvector/nvector_serial.h>
62#if CHASTE_SUNDIALS_VERSION >= 30000
63#if CHASTE_SUNDIALS_VERSION < 70000
64#include <cvode/cvode_direct.h>
66#include <sundials/sundials_types.h>
67#include <sunlinsol/sunlinsol_dense.h>
68#include <sunmatrix/sunmatrix_dense.h>
70#include <sundials/sundials_dense.h>
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
79#define CHASTE_CVODE_DENSE_MATRIX DenseMat
83#if CHASTE_SUNDIALS_VERSION >= 30000
84#define IJth(A, i, j) SM_ELEMENT_D(A, i, j)
86#define IJth(A, i, j) DENSE_ELEM(A, i, j)
146 friend class TestAbstractCvodeSystem;
148 friend class boost::serialization::access;
155 template <
class Archive>
156 void save(Archive& archive,
const unsigned int version)
const
197 template <
class Archive>
198 void load(Archive& archive,
const unsigned int version)
210 std::vector<double> state_vars;
214 std::vector<double> parameters;
217 std::vector<std::string> param_names;
218 archive& param_names;
233 BOOST_SERIALIZATION_SPLIT_MEMBER()
265 void CvodeError(
int flag, const
char* msg, const
double& rTime,
266 const
double& rStartTime, const
double& rEndTime);
280#if CHASTE_SUNDIALS_VERSION >= 30000
282 SUNMatrix mpSundialsDenseMatrix;
284 SUNLinearSolver mpSundialsLinearSolver;
359 CHASTE_CVODE_DENSE_MATRIX jacobian,
362 EXCEPTION(
"No analytic Jacobian has been defined for this system.");
454 void Solve(realtype tStart,
479 void SetTolerances(
double relTol = 1e-5,
double absTol = 1e-7);
#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)
bool mUseAnalyticJacobian
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)
double GetRelativeTolerance()
void RecordStoppingPoint(double stopTime)
void SetMaxSteps(long int numSteps)
bool HasAnalyticJacobian() const
double GetAbsoluteTolerance()
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
bool mHasAnalyticJacobian
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 GetUseAnalyticJacobian() const
double mLastInternalStepSize
void SetMinimalReset(bool minimalReset)
bool GetForceReset()
Get whether we will force a solver reset on every call to Solve()
N_Vector mLastSolutionState
virtual ~AbstractCvodeSystem()
void CheckParametersOnLoad(const std::vector< double > &rParameters, const std::vector< std::string > &rParameterNames)
unsigned mNumberOfStateVariables
const std::vector< std::string > & rGetParameterNames() const