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 mUseAnalyticJacobian(false),
00134 mpCvodeMem(NULL),
00135 mMaxSteps(0),
00136 mLastInternalStepSize(0)
00137 {
00138 SetTolerances();
00139 }
00140
00141 void AbstractCvodeSystem::Init()
00142 {
00143 DeleteVector(mStateVariables);
00144 mStateVariables = GetInitialConditions();
00145 DeleteVector(mParameters);
00146 mParameters = N_VNew_Serial(rGetParameterNames().size());
00147 for (int i=0; i<NV_LENGTH_S(mParameters); i++)
00148 {
00149 NV_Ith_S(mParameters, i) = 0.0;
00150 }
00151 }
00152
00153 AbstractCvodeSystem::~AbstractCvodeSystem()
00154 {
00155 FreeCvodeMemory();
00156 DeleteVector(mStateVariables);
00157 DeleteVector(mParameters);
00158 DeleteVector(mLastSolutionState);
00159 }
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174
00175 OdeSolution AbstractCvodeSystem::Solve(realtype tStart,
00176 realtype tEnd,
00177 realtype maxDt,
00178 realtype tSamp)
00179 {
00180 assert(tEnd >= tStart);
00181 assert(tSamp > 0.0);
00182
00183 SetupCvode(mStateVariables, tStart, maxDt);
00184
00185 TimeStepper stepper(tStart, tEnd, tSamp);
00186
00187
00188 OdeSolution solutions;
00189 solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
00190 solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
00191 solutions.rGetTimes().push_back(tStart);
00192 solutions.SetOdeSystemInformation(mpSystemInfo);
00193
00194
00195 while (!stepper.IsTimeAtEnd())
00196 {
00197
00198 int ierr = CVodeSetStopTime(mpCvodeMem, stepper.GetNextTime());
00199 assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr);
00200
00201 double cvode_stopped_at;
00202 ierr = CVode(mpCvodeMem, stepper.GetNextTime(), mStateVariables,
00203 &cvode_stopped_at, CV_NORMAL);
00204 if (ierr<0)
00205 {
00206
00207 CvodeError(ierr, "CVODE failed to solve system");
00208 }
00209
00210 assert(fabs(cvode_stopped_at - stepper.GetNextTime()) < DBL_EPSILON);
00211 #ifndef NDEBUG
00212 VerifyStateVariables();
00213 #endif
00214
00215 solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
00216 solutions.rGetTimes().push_back(cvode_stopped_at);
00217 stepper.AdvanceOneTimeStep();
00218 }
00219
00220
00221 solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
00222
00223 int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00224 assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr);
00225
00226 RecordStoppingPoint(tEnd);
00227
00228 return solutions;
00229 }
00230
00231 void AbstractCvodeSystem::Solve(realtype tStart,
00232 realtype tEnd,
00233 realtype maxDt)
00234 {
00235 assert(tEnd >= tStart);
00236
00237 SetupCvode(mStateVariables, tStart, maxDt);
00238
00239
00240 int ierr = CVodeSetStopTime(mpCvodeMem, tEnd);
00241 assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr);
00242
00243 double cvode_stopped_at;
00244 ierr = CVode(mpCvodeMem, tEnd, mStateVariables, &cvode_stopped_at, CV_NORMAL);
00245 if (ierr<0)
00246 {
00247
00248 CvodeError(ierr, "CVODE failed to solve system");
00249 }
00250
00251 assert(fabs(cvode_stopped_at - tEnd) < DBL_EPSILON);
00252
00253 ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00254 assert(ierr == CV_SUCCESS); UNUSED_OPT(ierr);
00255
00256 RecordStoppingPoint(cvode_stopped_at);
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272
00273
00274
00275
00276
00277
00278 #ifndef NDEBUG
00279 VerifyStateVariables();
00280 #endif
00281 }
00282
00283
00284 void AbstractCvodeSystem::SetMaxSteps(long int numSteps)
00285 {
00286 mMaxSteps = numSteps;
00287 }
00288
00289 long int AbstractCvodeSystem::GetMaxSteps()
00290 {
00291 return mMaxSteps;
00292 }
00293
00294 void AbstractCvodeSystem::SetTolerances(double relTol, double absTol)
00295 {
00296 mRelTol = relTol;
00297 mAbsTol = absTol;
00298 ResetSolver();
00299 }
00300
00301 double AbstractCvodeSystem::GetRelativeTolerance()
00302 {
00303 return mRelTol;
00304 }
00305
00306 double AbstractCvodeSystem::GetAbsoluteTolerance()
00307 {
00308 return mAbsTol;
00309 }
00310
00311 double AbstractCvodeSystem::GetLastStepSize()
00312 {
00313 return mLastInternalStepSize;
00314 }
00315
00316 void AbstractCvodeSystem::SetForceReset(bool autoReset)
00317 {
00318 mForceReset = autoReset;
00319 if (mForceReset)
00320 {
00321 ResetSolver();
00322 }
00323 }
00324
00325 void AbstractCvodeSystem::SetMinimalReset(bool minimalReset)
00326 {
00327 mForceMinimalReset = minimalReset;
00328 if (mForceMinimalReset)
00329 {
00330 SetForceReset(false);
00331 }
00332 }
00333
00334
00335 void AbstractCvodeSystem::ResetSolver()
00336 {
00337 DeleteVector(mLastSolutionState);
00338 }
00339
00340
00341 void AbstractCvodeSystem::SetupCvode(N_Vector initialConditions,
00342 realtype tStart,
00343 realtype maxDt)
00344 {
00345 assert((unsigned)NV_LENGTH_S(initialConditions) == GetNumberOfStateVariables());
00346 assert(maxDt >= 0.0);
00347
00348
00349
00350 bool reinit = !mpCvodeMem || mForceReset || !mLastSolutionState || !CompareDoubles::WithinAnyTolerance(tStart, mLastSolutionTime);
00351 if (!reinit && !mForceMinimalReset)
00352 {
00353 const unsigned size = GetNumberOfStateVariables();
00354 for (unsigned i=0; i<size; i++)
00355 {
00356 if (!CompareDoubles::WithinAnyTolerance(GetVectorComponent(mLastSolutionState, i), GetVectorComponent(mStateVariables, i)))
00357 {
00358 reinit = true;
00359 break;
00360 }
00361 }
00362 }
00363
00364 if (!mpCvodeMem)
00365 {
00366
00367 mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
00368 if (mpCvodeMem == NULL) EXCEPTION("Failed to SetupCvode CVODE");
00369
00370
00371 CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, NULL);
00372
00373 #if CHASTE_SUNDIALS_VERSION >= 20400
00374 CVodeSetUserData(mpCvodeMem, (void*)(this));
00375 #else
00376 CVodeSetFdata(mpCvodeMem, (void*)(this));
00377 #endif
00378
00379 #if CHASTE_SUNDIALS_VERSION >= 20400
00380 CVodeInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions);
00381 CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
00382 #else
00383 CVodeMalloc(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
00384 CV_SS, mRelTol, &mAbsTol);
00385 #endif
00386
00387 CVDense(mpCvodeMem, NV_LENGTH_S(initialConditions));
00388
00389 if (mUseAnalyticJacobian)
00390 {
00391 #if CHASTE_SUNDIALS_VERSION >= 20400
00392 CVDlsSetDenseJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor);
00393 #else
00394 CVDenseSetJacFn(mpCvodeMem, AbstractCvodeSystemJacAdaptor, (void*)(this));
00395 #endif
00396 }
00397 }
00398 else if (reinit)
00399 {
00400
00401 #if CHASTE_SUNDIALS_VERSION >= 20400
00402 CVodeReInit(mpCvodeMem, tStart, initialConditions);
00403 CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
00404 #else
00405 CVodeReInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
00406 CV_SS, mRelTol, &mAbsTol);
00407 #endif
00408 }
00409
00410
00411 CVodeSetMaxStep(mpCvodeMem, maxDt);
00412 if (mMaxSteps > 0)
00413 {
00414 CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
00415 CVodeSetMaxErrTestFails(mpCvodeMem, 15);
00416 }
00417 }
00418
00419
00420 void AbstractCvodeSystem::RecordStoppingPoint(double stopTime)
00421 {
00422
00423
00424
00425
00426
00427 if (mForceReset) return;
00428
00429
00430
00431 const unsigned size = GetNumberOfStateVariables();
00432 CreateVectorIfEmpty(mLastSolutionState, size);
00433 for (unsigned i=0; i<size; i++)
00434 {
00435 SetVectorComponent(mLastSolutionState, i, GetVectorComponent(mStateVariables, i));
00436 }
00437 mLastSolutionTime = stopTime;
00438 }
00439
00440
00441 void AbstractCvodeSystem::FreeCvodeMemory()
00442 {
00443 if (mpCvodeMem)
00444 {
00445 CVodeFree(&mpCvodeMem);
00446 }
00447 mpCvodeMem = NULL;
00448 }
00449
00450
00451 void AbstractCvodeSystem::CvodeError(int flag, const char * msg)
00452 {
00453
00454 std::stringstream err;
00455 char* p_flag_name = CVodeGetReturnFlagName(flag);
00456 err << msg << ": " << p_flag_name;
00457 free(p_flag_name);
00458 if (flag == CV_LSETUP_FAIL)
00459 {
00460 #if CHASTE_SUNDIALS_VERSION >= 20500
00461 long int ls_flag;
00462 #else
00463 int ls_flag;
00464 #endif
00465 char* p_ls_flag_name;
00466 #if CHASTE_SUNDIALS_VERSION >= 20400
00467 CVDlsGetLastFlag(mpCvodeMem, &ls_flag);
00468 p_ls_flag_name = CVDlsGetReturnFlagName(ls_flag);
00469 #else
00470 CVDenseGetLastFlag(mpCvodeMem, &ls_flag);
00471 p_ls_flag_name = CVDenseGetReturnFlagName(ls_flag);
00472 #endif
00473 err << " (LS flag=" << ls_flag << ":" << p_ls_flag_name << ")";
00474 delete p_ls_flag_name;
00475 }
00476 FreeCvodeMemory();
00477 std::cerr << err.str() << std::endl << std::flush;
00478 EXCEPTION(err.str());
00479 }
00480
00481
00482 bool AbstractCvodeSystem::GetUseAnalyticJacobian()
00483 {
00484 return mUseAnalyticJacobian;
00485 }
00486
00487
00488 void AbstractCvodeSystem::ForceUseOfNumericalJacobian(bool useNumericalJacobian)
00489 {
00490 if (mUseAnalyticJacobian == useNumericalJacobian)
00491 {
00492 mUseAnalyticJacobian = !useNumericalJacobian;
00493
00494 this->FreeCvodeMemory();
00495 }
00496 }
00497
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 #endif // CHASTE_CVODE