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
00030
00031
00032
00033
00034
00035
00036 #ifdef CHASTE_CVODE
00037
00038 #include <sstream>
00039 #include <cassert>
00040
00041 #include "AbstractCvodeSystem.hpp"
00042 #include "Exception.hpp"
00043 #include "VectorHelperFunctions.hpp"
00044 #include "MathsCustomFunctions.hpp"
00045 #include "TimeStepper.hpp"
00046 #include "CvodeAdaptor.hpp"
00047
00048
00049 #include <cvode/cvode.h>
00050 #include <sundials/sundials_nvector.h>
00051 #include <cvode/cvode_dense.h>
00052
00053
00054
00055
00056
00057
00058
00060
00061
00062
00063
00064
00074 int AbstractCvodeSystemRhsAdaptor(realtype t, N_Vector y, N_Vector ydot, void *pData)
00075 {
00076 assert(pData != NULL);
00077 AbstractCvodeSystem* p_ode_system = (AbstractCvodeSystem*) pData;
00078 try
00079 {
00080 p_ode_system->EvaluateYDerivatives(t, y, ydot);
00081 }
00082 catch (const Exception &e)
00083 {
00084 std::cerr << "CVODE RHS Exception: " << e.GetMessage()
00085 << std::endl << std::flush;
00086 return -1;
00087 }
00088 return 0;
00089 }
00090
00091
00092
00093
00094
00095 #if CHASTE_SUNDIALS_VERSION >= 20500
00096
00097 int AbstractCvodeSystemJacAdaptor(long int N, realtype t, N_Vector y, N_Vector ydot, DlsMat jacobian,
00098 #elif CHASTE_SUNDIALS_VERSION >= 20400
00099
00100 int AbstractCvodeSystemJacAdaptor(int N, realtype t, N_Vector y, N_Vector ydot, DlsMat jacobian,
00101 #else
00102
00103 int AbstractCvodeSystemJacAdaptor(long int N, DenseMat jacobian, realtype t, N_Vector y, N_Vector ydot,
00104 #endif
00105 void *pData, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
00106 {
00107 assert(pData != NULL);
00108 AbstractCvodeSystem* p_ode_system = (AbstractCvodeSystem*) pData;
00109 try
00110 {
00111 p_ode_system->EvaluateAnalyticJacobian((long)(N), t, y, ydot, jacobian, tmp1, tmp2, tmp3);
00112 }
00113 catch (const Exception &e)
00114 {
00115 std::cerr << "CVODE Jacobian Exception: " << e.GetMessage() << std::endl << std::flush;
00116 return -1;
00117 }
00118 return 0;
00119 }
00120
00121 AbstractCvodeSystem::AbstractCvodeSystem(unsigned numberOfStateVariables)
00122 : AbstractParameterisedSystem<N_Vector>(numberOfStateVariables),
00123 mLastSolutionState(NULL),
00124 mLastSolutionTime(0.0),
00125 #if CHASTE_SUNDIALS_VERSION >=20400
00126 mForceReset(false),
00127 #else
00128
00129
00130 mForceReset(true),
00131 #endif
00132 mForceMinimalReset(false),
00133 mHasAnalyticJacobian(false),
00134 mUseAnalyticJacobian(false),
00135 mpCvodeMem(NULL),
00136 mMaxSteps(0),
00137 mLastInternalStepSize(0)
00138 {
00139 SetTolerances();
00140 }
00141
00142 void AbstractCvodeSystem::Init()
00143 {
00144 DeleteVector(mStateVariables);
00145 mStateVariables = GetInitialConditions();
00146 DeleteVector(mParameters);
00147 mParameters = N_VNew_Serial(rGetParameterNames().size());
00148 for (int i=0; i<NV_LENGTH_S(mParameters); i++)
00149 {
00150 NV_Ith_S(mParameters, i) = 0.0;
00151 }
00152 }
00153
00154 AbstractCvodeSystem::~AbstractCvodeSystem()
00155 {
00156 FreeCvodeMemory();
00157 DeleteVector(mStateVariables);
00158 DeleteVector(mParameters);
00159 DeleteVector(mLastSolutionState);
00160 }
00161
00162
00163
00164
00165
00166
00167
00168
00169 OdeSolution AbstractCvodeSystem::Solve(realtype tStart,
00170 realtype tEnd,
00171 realtype maxDt,
00172 realtype tSamp)
00173 {
00174 assert(tEnd >= tStart);
00175 assert(tSamp > 0.0);
00176
00177 SetupCvode(mStateVariables, tStart, maxDt);
00178
00179 TimeStepper stepper(tStart, tEnd, tSamp);
00180
00181
00182 OdeSolution solutions;
00183 solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
00184 solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
00185 solutions.rGetTimes().push_back(tStart);
00186 solutions.SetOdeSystemInformation(mpSystemInfo);
00187
00188
00189 while (!stepper.IsTimeAtEnd())
00190 {
00191
00192 int ierr = CVodeSetStopTime(mpCvodeMem, stepper.GetNextTime());
00193 assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr);
00194
00195 double cvode_stopped_at;
00196 ierr = CVode(mpCvodeMem, stepper.GetNextTime(), mStateVariables,
00197 &cvode_stopped_at, CV_NORMAL);
00198 if (ierr<0)
00199 {
00200
00201 CvodeError(ierr, "CVODE failed to solve system");
00202 }
00203
00204 assert(fabs(cvode_stopped_at - stepper.GetNextTime()) < DBL_EPSILON);
00205 #ifndef NDEBUG
00206 VerifyStateVariables();
00207 #endif
00208
00209 solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
00210 solutions.rGetTimes().push_back(cvode_stopped_at);
00211 stepper.AdvanceOneTimeStep();
00212 }
00213
00214
00215 solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
00216
00217 int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00218 assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr);
00219
00220 RecordStoppingPoint(tEnd);
00221
00222 return solutions;
00223 }
00224
00225 void AbstractCvodeSystem::Solve(realtype tStart,
00226 realtype tEnd,
00227 realtype maxDt)
00228 {
00229 assert(tEnd >= tStart);
00230
00231 SetupCvode(mStateVariables, tStart, maxDt);
00232
00233
00234 int ierr = CVodeSetStopTime(mpCvodeMem, tEnd);
00235 assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr);
00236
00237 double cvode_stopped_at;
00238 ierr = CVode(mpCvodeMem, tEnd, mStateVariables, &cvode_stopped_at, CV_NORMAL);
00239 if (ierr<0)
00240 {
00241
00242 CvodeError(ierr, "CVODE failed to solve system");
00243 }
00244
00245 assert(fabs(cvode_stopped_at - tEnd) < DBL_EPSILON);
00246
00247 ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00248 assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr);
00249
00250 RecordStoppingPoint(cvode_stopped_at);
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272 #ifndef NDEBUG
00273 VerifyStateVariables();
00274 #endif
00275 }
00276
00277
00278 void AbstractCvodeSystem::SetMaxSteps(long int numSteps)
00279 {
00280 mMaxSteps = numSteps;
00281 }
00282
00283 long int AbstractCvodeSystem::GetMaxSteps()
00284 {
00285 return mMaxSteps;
00286 }
00287
00288 void AbstractCvodeSystem::SetTolerances(double relTol, double absTol)
00289 {
00290 mRelTol = relTol;
00291 mAbsTol = absTol;
00292 ResetSolver();
00293 }
00294
00295 double AbstractCvodeSystem::GetRelativeTolerance()
00296 {
00297 return mRelTol;
00298 }
00299
00300 double AbstractCvodeSystem::GetAbsoluteTolerance()
00301 {
00302 return mAbsTol;
00303 }
00304
00305 double AbstractCvodeSystem::GetLastStepSize()
00306 {
00307 return mLastInternalStepSize;
00308 }
00309
00310 void AbstractCvodeSystem::SetForceReset(bool autoReset)
00311 {
00312 mForceReset = autoReset;
00313 if (mForceReset)
00314 {
00315 ResetSolver();
00316 }
00317 }
00318
00319 void AbstractCvodeSystem::SetMinimalReset(bool minimalReset)
00320 {
00321 mForceMinimalReset = minimalReset;
00322 if (mForceMinimalReset)
00323 {
00324 SetForceReset(false);
00325 }
00326 }
00327
00328
00329 void AbstractCvodeSystem::ResetSolver()
00330 {
00331 DeleteVector(mLastSolutionState);
00332 }
00333
00334
00335 void AbstractCvodeSystem::SetupCvode(N_Vector initialConditions,
00336 realtype tStart,
00337 realtype maxDt)
00338 {
00339 assert((unsigned)NV_LENGTH_S(initialConditions) == GetNumberOfStateVariables());
00340 assert(maxDt >= 0.0);
00341
00342
00343
00344 bool reinit = !mpCvodeMem || mForceReset || !mLastSolutionState || !CompareDoubles::WithinAnyTolerance(tStart, mLastSolutionTime);
00345 if (!reinit && !mForceMinimalReset)
00346 {
00347 const unsigned size = GetNumberOfStateVariables();
00348 for (unsigned i=0; i<size; i++)
00349 {
00350 if (!CompareDoubles::WithinAnyTolerance(GetVectorComponent(mLastSolutionState, i), GetVectorComponent(mStateVariables, i)))
00351 {
00352 reinit = true;
00353 break;
00354 }
00355 }
00356 }
00357
00358 if (!mpCvodeMem)
00359 {
00360
00361 mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
00362 if (mpCvodeMem == NULL) EXCEPTION("Failed to SetupCvode CVODE");
00363
00364
00365 CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, NULL);
00366
00367 #if CHASTE_SUNDIALS_VERSION >= 20400
00368 CVodeSetUserData(mpCvodeMem, (void*)(this));
00369 #else
00370 CVodeSetFdata(mpCvodeMem, (void*)(this));
00371 #endif
00372
00373 #if CHASTE_SUNDIALS_VERSION >= 20400
00374 CVodeInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions);
00375 CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
00376 #else
00377 CVodeMalloc(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
00378 CV_SS, mRelTol, &mAbsTol);
00379 #endif
00380
00381 CVDense(mpCvodeMem, NV_LENGTH_S(initialConditions));
00382
00383 if (mUseAnalyticJacobian)
00384 {
00385 #if CHASTE_SUNDIALS_VERSION >= 20400
00386 CVDlsSetDenseJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor);
00387 #else
00388 CVDenseSetJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor, (void*)(this));
00389 #endif
00390 }
00391 }
00392 else if (reinit)
00393 {
00394
00395 #if CHASTE_SUNDIALS_VERSION >= 20400
00396 CVodeReInit(mpCvodeMem, tStart, initialConditions);
00397 CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
00398 #else
00399 CVodeReInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
00400 CV_SS, mRelTol, &mAbsTol);
00401 #endif
00402 }
00403
00404
00405 CVodeSetMaxStep(mpCvodeMem, maxDt);
00406 if (mMaxSteps > 0)
00407 {
00408 CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
00409 CVodeSetMaxErrTestFails(mpCvodeMem, 15);
00410 }
00411 }
00412
00413
00414 void AbstractCvodeSystem::RecordStoppingPoint(double stopTime)
00415 {
00416
00417
00418
00419
00420
00421 if (mForceReset) return;
00422
00423
00424
00425 const unsigned size = GetNumberOfStateVariables();
00426 CreateVectorIfEmpty(mLastSolutionState, size);
00427 for (unsigned i=0; i<size; i++)
00428 {
00429 SetVectorComponent(mLastSolutionState, i, GetVectorComponent(mStateVariables, i));
00430 }
00431 mLastSolutionTime = stopTime;
00432 }
00433
00434
00435 void AbstractCvodeSystem::FreeCvodeMemory()
00436 {
00437 if (mpCvodeMem)
00438 {
00439 CVodeFree(&mpCvodeMem);
00440 }
00441 mpCvodeMem = NULL;
00442 }
00443
00444
00445 void AbstractCvodeSystem::CvodeError(int flag, const char * msg)
00446 {
00447
00448 std::stringstream err;
00449 char* p_flag_name = CVodeGetReturnFlagName(flag);
00450 err << msg << ": " << p_flag_name;
00451 free(p_flag_name);
00452 if (flag == CV_LSETUP_FAIL)
00453 {
00454 #if CHASTE_SUNDIALS_VERSION >= 20500
00455 long int ls_flag;
00456 #else
00457 int ls_flag;
00458 #endif
00459 char* p_ls_flag_name;
00460 #if CHASTE_SUNDIALS_VERSION >= 20400
00461 CVDlsGetLastFlag(mpCvodeMem, &ls_flag);
00462 p_ls_flag_name = CVDlsGetReturnFlagName(ls_flag);
00463 #else
00464 CVDenseGetLastFlag(mpCvodeMem, &ls_flag);
00465 p_ls_flag_name = CVDenseGetReturnFlagName(ls_flag);
00466 #endif
00467 err << " (LS flag=" << ls_flag << ":" << p_ls_flag_name << ")";
00468 delete p_ls_flag_name;
00469 }
00470 FreeCvodeMemory();
00471 std::cerr << err.str() << std::endl << std::flush;
00472 EXCEPTION(err.str());
00473 }
00474
00475 bool AbstractCvodeSystem::HasAnalyticJacobian() const
00476 {
00477 return mHasAnalyticJacobian;
00478 }
00479
00480 bool AbstractCvodeSystem::GetUseAnalyticJacobian() const
00481 {
00482 return mUseAnalyticJacobian;
00483 }
00484
00485
00486 void AbstractCvodeSystem::ForceUseOfNumericalJacobian(bool useNumericalJacobian)
00487 {
00488 if (!useNumericalJacobian && !mHasAnalyticJacobian)
00489 {
00490 EXCEPTION("Analytic Jacobian requested, but this ODE system doesn't have one. You can check this with HasAnalyticJacobian().");
00491 }
00492
00493 if (mUseAnalyticJacobian == useNumericalJacobian)
00494 {
00495 mUseAnalyticJacobian = !useNumericalJacobian;
00496
00497 this->FreeCvodeMemory();
00498 }
00499 }
00500
00501
00502
00503
00504
00505
00506
00507
00508
00509
00510
00511
00512
00513
00514
00515
00516
00517
00518
00519
00520
00521
00522
00523
00524
00525
00526
00527
00528
00529
00530
00531
00532
00533
00534
00535
00536
00537
00538
00539
00540
00541
00542
00543
00544
00545
00546
00547
00548
00549
00550
00551
00552
00553 #endif // CHASTE_CVODE