00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #ifdef CHASTE_CVODE
00030
00031 #include <sstream>
00032 #include <iostream>
00033 #include <cmath>
00034
00035 #include "AbstractCvodeCell.hpp"
00036 #include "CvodeAdaptor.hpp"
00037 #include "TimeStepper.hpp"
00038 #include "Exception.hpp"
00039 #include "HeartConfig.hpp"
00040 #include "VectorHelperFunctions.hpp"
00041
00042
00043 #include <cvode/cvode.h>
00044 #include <sundials/sundials_nvector.h>
00045 #include <cvode/cvode_dense.h>
00046
00056 int AbstractCvodeCellRhsAdaptor(realtype t, N_Vector y, N_Vector ydot, void *pData)
00057 {
00058 assert(pData != NULL);
00059 AbstractCvodeCell* pCell = (AbstractCvodeCell*) pData;
00060 try
00061 {
00062 pCell->EvaluateRhs(t, y, ydot);
00063 }
00064 catch (const Exception &e)
00065 {
00066 std::cerr << "CVODE RHS Exception: " << e.GetMessage()
00067 << std::endl << std::flush;
00068 return -1;
00069 }
00070 return 0;
00071 }
00072
00073
00074 AbstractCvodeCell::AbstractCvodeCell(boost::shared_ptr<AbstractIvpOdeSolver> ,
00075 unsigned numberOfStateVariables,
00076 unsigned voltageIndex,
00077 boost::shared_ptr<AbstractStimulusFunction> pIntracellularStimulus)
00078 : AbstractCardiacCellInterface(boost::shared_ptr<AbstractIvpOdeSolver>(), voltageIndex, pIntracellularStimulus),
00079 AbstractParameterisedSystem<N_Vector>(numberOfStateVariables),
00080 mpCvodeMem(NULL),
00081 mMaxSteps(0),
00082 mMaxDt(DOUBLE_UNSET)
00083 {
00084 SetTolerances();
00085 }
00086
00087
00088 AbstractCvodeCell::~AbstractCvodeCell()
00089 {
00090 DeleteVector(mStateVariables);
00091 DeleteVector(mParameters);
00092 }
00093
00094
00095 void AbstractCvodeCell::Init()
00096 {
00097 DeleteVector(mStateVariables);
00098 mStateVariables = GetInitialConditions();
00099 DeleteVector(mParameters);
00100 mParameters = N_VNew_Serial(rGetParameterNames().size());
00101 for (int i=0; i<NV_LENGTH_S(mParameters); i++)
00102 {
00103 NV_Ith_S(mParameters, i) = 0.0;
00104 }
00105 }
00106
00107
00108 double AbstractCvodeCell::GetVoltage()
00109 {
00110 assert(mStateVariables);
00111 return GetAnyVariable(mVoltageIndex);
00112 }
00113
00114 void AbstractCvodeCell::SetVoltage(double voltage)
00115 {
00116 assert(mStateVariables);
00117 SetAnyVariable(mVoltageIndex, voltage);
00118 }
00119
00120 void AbstractCvodeCell::SetVoltageDerivativeToZero(bool clamp)
00121 {
00122 mSetVoltageDerivativeToZero = clamp;
00123 }
00124
00125
00126 void AbstractCvodeCell::SetTimestep(double maxDt)
00127 {
00128 mMaxDt = maxDt;
00129 }
00130
00131
00132 void AbstractCvodeCell::SolveAndUpdateState(double tStart, double tEnd)
00133 {
00134 if (mMaxDt == DOUBLE_UNSET)
00135 {
00136 SetTimestep(HeartConfig::Instance()->GetPrintingTimeStep());
00137 }
00138 Solve(tStart, tEnd, mMaxDt);
00139 }
00140
00141 OdeSolution AbstractCvodeCell::Compute(double tStart, double tEnd, double tSamp)
00142 {
00143 if (tSamp == 0.0)
00144 {
00145 tSamp = HeartConfig::Instance()->GetPrintingTimeStep();
00146 }
00147 if (mMaxDt == DOUBLE_UNSET)
00148 {
00149 SetTimestep(tSamp);
00150 }
00151 return Solve(tStart, tEnd, mMaxDt, tSamp);
00152 }
00153
00154
00155 void AbstractCvodeCell::ComputeExceptVoltage(double tStart, double tEnd)
00156 {
00157 EXCEPTION("This method is not yet implemented for CVODE cells.");
00158 }
00159
00160
00161 OdeSolution AbstractCvodeCell::Solve(realtype tStart,
00162 realtype tEnd,
00163 realtype maxDt,
00164 realtype tSamp)
00165 {
00166 assert(tEnd >= tStart);
00167 assert(tSamp > 0.0);
00168
00169 SetupCvode(mStateVariables, tStart, maxDt);
00170
00171 TimeStepper stepper(tStart, tEnd, tSamp);
00172
00173
00174 OdeSolution solutions;
00175 solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
00176 solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
00177 solutions.rGetTimes().push_back(tStart);
00178 solutions.SetOdeSystemInformation(mpSystemInfo);
00179
00180
00181 while (!stepper.IsTimeAtEnd())
00182 {
00183 double cvode_stopped_at;
00184 int ierr = CVode(mpCvodeMem, stepper.GetNextTime(), mStateVariables,
00185 &cvode_stopped_at, CV_NORMAL);
00186 if (ierr<0)
00187 {
00188 FreeCvodeMemory();
00189 CvodeError(ierr, "CVODE failed to solve system");
00190 }
00191
00192 assert(fabs(cvode_stopped_at - stepper.GetNextTime()) < DBL_EPSILON);
00193
00194 VerifyStateVariables();
00195
00196
00197 solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
00198 solutions.rGetTimes().push_back(cvode_stopped_at);
00199 stepper.AdvanceOneTimeStep();
00200 }
00201
00202
00203 solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
00204
00205 int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00206 assert(ierr == CV_SUCCESS); ierr=ierr;
00207 FreeCvodeMemory();
00208
00209 return solutions;
00210 }
00211
00212 void AbstractCvodeCell::Solve(realtype tStart,
00213 realtype tEnd,
00214 realtype maxDt)
00215 {
00216 assert(tEnd >= tStart);
00217
00218 SetupCvode(mStateVariables, tStart, maxDt);
00219
00220 double cvode_stopped_at;
00221 int ierr = CVode(mpCvodeMem, tEnd, mStateVariables, &cvode_stopped_at, CV_NORMAL);
00222 if (ierr<0)
00223 {
00224 FreeCvodeMemory();
00225 CvodeError(ierr, "CVODE failed to solve system");
00226 }
00227
00228 assert(fabs(cvode_stopped_at - tEnd) < DBL_EPSILON);
00229
00230 ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00231 assert(ierr == CV_SUCCESS); ierr=ierr;
00232 FreeCvodeMemory();
00233
00234 VerifyStateVariables();
00235 }
00236
00237
00238 void AbstractCvodeCell::SetMaxSteps(long int numSteps)
00239 {
00240 mMaxSteps = numSteps;
00241 }
00242
00243 long int AbstractCvodeCell::GetMaxSteps()
00244 {
00245 return mMaxSteps;
00246 }
00247
00248 void AbstractCvodeCell::SetTolerances(double relTol, double absTol)
00249 {
00250 mRelTol = relTol;
00251 mAbsTol = absTol;
00252 }
00253
00254 double AbstractCvodeCell::GetRelativeTolerance()
00255 {
00256 return mRelTol;
00257 }
00258
00259 double AbstractCvodeCell::GetAbsoluteTolerance()
00260 {
00261 return mAbsTol;
00262 }
00263
00264 double AbstractCvodeCell::GetLastStepSize()
00265 {
00266 return mLastInternalStepSize;
00267 }
00268
00269
00270 void AbstractCvodeCell::SetupCvode(N_Vector initialConditions,
00271 realtype tStart,
00272 realtype maxDt)
00273 {
00274 assert((unsigned)NV_LENGTH_S(initialConditions) == GetNumberOfStateVariables());
00275 assert(maxDt >= 0.0);
00276
00277 mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
00278 if (mpCvodeMem == NULL) EXCEPTION("Failed to SetupCvode CVODE");
00279
00280 CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, NULL);
00281
00282 #if CHASTE_SUNDIALS_VERSION >= 20400
00283 CVodeSetUserData(mpCvodeMem, (void*)(this));
00284 CVodeInit(mpCvodeMem, AbstractCvodeCellRhsAdaptor, tStart, initialConditions);
00285 CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
00286 #else
00287 CVodeSetFdata(mpCvodeMem, (void*)(this));
00288 CVodeMalloc(mpCvodeMem, AbstractCvodeCellRhsAdaptor, tStart, initialConditions,
00289 CV_SS, mRelTol, &mAbsTol);
00290 #endif
00291 CVodeSetMaxStep(mpCvodeMem, maxDt);
00292
00293 if (mMaxSteps > 0)
00294 {
00295 CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
00296 }
00297
00298 CVDense(mpCvodeMem, NV_LENGTH_S(initialConditions));
00299 }
00300
00301
00302 void AbstractCvodeCell::FreeCvodeMemory()
00303 {
00304 CVodeFree(&mpCvodeMem);
00305 }
00306
00307
00308 void AbstractCvodeCell::CvodeError(int flag, const char * msg)
00309 {
00310
00311 std::stringstream err;
00312 char* p_flag_name = CVodeGetReturnFlagName(flag);
00313 err << msg << ": " << p_flag_name;
00314 free(p_flag_name);
00315 std::cerr << err.str() << std::endl << std::flush;
00316 EXCEPTION(err.str());
00317 }
00318
00319
00320 std::vector<double> AbstractCvodeCell::MakeStdVec(N_Vector v)
00321 {
00322 std::vector<double> sv;
00323 CopyToStdVector(v, sv);
00324 return sv;
00325 }
00326
00327
00328 std::string AbstractCvodeCell::DumpState(const std::string& message,
00329 N_Vector Y)
00330 {
00331 std::stringstream res;
00332 res << message << "\nState:\n";
00333 if (Y == NULL)
00334 {
00335 Y = mStateVariables;
00336 }
00337 for (int i=0; i<NV_LENGTH_S(Y); i++)
00338 {
00339 res << "\t" << rGetStateVariableNames()[i] << ":" << NV_Ith_S(Y, i) << "\n";
00340 }
00341 return res.str();
00342 }
00343
00344
00345 #endif // CHASTE_CVODE