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 "CvodeAdaptor.hpp"
00032
00033 #include "Exception.hpp"
00034 #include "TimeStepper.hpp"
00035
00036 #include <iostream>
00037
00038
00039 #include <sundials/sundials_nvector.h>
00040 #include <cvode/cvode_dense.h>
00041
00042
00049 void CopyToStdVector(N_Vector src, std::vector<realtype>& rDest)
00050 {
00051
00052 realtype* p_src = NV_DATA_S(src);
00053 if (p_src == &(rDest[0])) return;
00054
00055 long size = NV_LENGTH_S(src);
00056 rDest.resize(size);
00057
00058 for (long i=0; i<size; i++)
00059 {
00060 rDest[i] = p_src[i];
00061 }
00062 }
00063
00080 int CvodeRhsAdaptor(realtype t, N_Vector y, N_Vector ydot, void* pData)
00081 {
00082 assert(pData != NULL);
00083 CvodeData* p_data = (CvodeData*) pData;
00084
00085 static std::vector<realtype> ydot_vec;
00086 CopyToStdVector(y, *p_data->pY);
00087 CopyToStdVector(ydot, ydot_vec);
00088
00089 try
00090 {
00091 p_data->pSystem->EvaluateYDerivatives(t, *(p_data->pY), ydot_vec);
00092 }
00093 catch (Exception &e)
00094 {
00095 std::cerr << "CVODE RHS Exception: " << e.GetMessage() << std::endl << std::flush;
00096 return -1;
00097 }
00098
00099 long len = NV_LENGTH_S(ydot);
00100 realtype* p_ydot = NV_DATA_S(ydot);
00101 for (long i=0; i<len; i++)
00102 {
00103 p_ydot[i] = ydot_vec[i];
00104 }
00105 return 0;
00106 }
00107
00128 int CvodeRootAdaptor(realtype t, N_Vector y, realtype* pGOut, void* pData)
00129 {
00130 assert(pData != NULL);
00131 CvodeData* p_data = (CvodeData*) pData;
00132
00133 CopyToStdVector(y, *p_data->pY);
00134
00135 try
00136 {
00137 *pGOut = p_data->pSystem->CalculateRootFunction(t, *p_data->pY);
00138 }
00139 catch (Exception &e)
00140 {
00141 std::cerr << "CVODE Root Exception: " << e.GetMessage() << std::endl << std::flush;
00142 return -1;
00143 }
00144 return 0;
00145 }
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183
00184 void CvodeErrorHandler(int errorCode, const char *module, const char *function,
00185 char *message, void* pData)
00186 {
00187 std::stringstream err;
00188 err << "CVODE Error " << errorCode << " in module " << module
00189 << " function " << function << ": " << message;
00190 std::cerr << "*" << err.str() << std::endl << std::flush;
00191 EXCEPTION(err.str());
00192 }
00193
00194
00195
00196 void CvodeAdaptor::SetupCvode(AbstractOdeSystem* pOdeSystem,
00197 std::vector<double>& rInitialY,
00198 double startTime, double maxStep)
00199 {
00200 assert(rInitialY.size() == pOdeSystem->GetNumberOfStateVariables());
00201 assert(maxStep > 0.0);
00202
00203 mInitialValues = N_VMake_Serial(rInitialY.size(), &(rInitialY[0]));
00204 assert(NV_DATA_S(mInitialValues) == &(rInitialY[0]));
00205 assert(!NV_OWN_DATA_S(mInitialValues));
00206
00207
00208
00209
00210 mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
00211 if (mpCvodeMem == NULL) EXCEPTION("Failed to SetupCvode CVODE");
00212
00213 CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, NULL);
00214
00215 mData.pSystem = pOdeSystem;
00216 mData.pY = &rInitialY;
00217 CVodeSetFdata(mpCvodeMem, (void*)(&mData));
00218
00219 CVodeMalloc(mpCvodeMem, CvodeRhsAdaptor, startTime, mInitialValues,
00220 CV_SS, mRelTol, &mAbsTol);
00221 CVodeSetMaxStep(mpCvodeMem, maxStep);
00222
00223 if (mCheckForRoots)
00224 {
00225 CVodeRootInit(mpCvodeMem, 1, CvodeRootAdaptor, (void*)(&mData));
00226 }
00227
00228 if (mMaxSteps > 0)
00229 {
00230 CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
00231 }
00232
00233 CVDense(mpCvodeMem, rInitialY.size());
00234 }
00235
00236 void CvodeAdaptor::FreeCvodeMemory()
00237 {
00238
00239
00240
00241 N_VDestroy_Serial(mInitialValues); mInitialValues = NULL;
00242
00243
00244 CVodeFree(&mpCvodeMem);
00245 }
00246
00247
00248 #define COVERAGE_IGNORE
00249 void CvodeAdaptor::CvodeError(int flag, const char * msg)
00250 {
00251 std::stringstream err;
00252 err << msg << ": " << CVodeGetReturnFlagName(flag);
00253 std::cerr << err.str() << std::endl << std::flush;
00254 EXCEPTION(err.str());
00255 }
00256 #undef COVERAGE_IGNORE
00257
00258
00259 OdeSolution CvodeAdaptor::Solve(AbstractOdeSystem* pOdeSystem,
00260 std::vector<double>& rYValues,
00261 double startTime,
00262 double endTime,
00263 double maxStep,
00264 double timeSampling)
00265 {
00266 assert(endTime > startTime);
00267 assert(timeSampling > 0.0);
00268
00269 mStoppingEventOccurred = false;
00270 if (mCheckForRoots && pOdeSystem->CalculateStoppingEvent(startTime, rYValues) == true)
00271 {
00272 EXCEPTION("(Solve with sampling) Stopping event is true for initial condition");
00273 }
00274
00275 SetupCvode(pOdeSystem, rYValues, startTime, maxStep);
00276
00277 TimeStepper stepper(startTime, endTime, timeSampling);
00278 N_Vector yout = mInitialValues;
00279
00280
00281 OdeSolution solutions;
00282 solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
00283 solutions.rGetSolutions().push_back(rYValues);
00284 solutions.rGetTimes().push_back(startTime);
00285 solutions.SetOdeSystemInformation(pOdeSystem->GetSystemInformation());
00286
00287
00288 while (!stepper.IsTimeAtEnd() && !mStoppingEventOccurred)
00289 {
00290
00291 double tend;
00292 int ierr;
00293 try
00294 {
00295 ierr = CVode(mpCvodeMem, stepper.GetNextTime(), yout,
00296 &tend, CV_NORMAL);
00297 }
00298 catch (...)
00299 {
00300 FreeCvodeMemory();
00301 throw;
00302 }
00303 if (ierr<0) CvodeError(ierr, "CVODE failed to solve system");
00304
00305 solutions.rGetSolutions().push_back(rYValues);
00306 solutions.rGetTimes().push_back(tend);
00307 if (ierr == CV_ROOT_RETURN)
00308 {
00309
00310 mStoppingEventOccurred = true;
00311 mStoppingTime = tend;
00312 }
00313 stepper.AdvanceOneTimeStep();
00314 }
00315
00316
00317 solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
00318
00319
00320 int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00321 assert(ierr == CV_SUCCESS); ierr=ierr;
00322 FreeCvodeMemory();
00323
00324 return solutions;
00325 }
00326
00327
00328 void CvodeAdaptor::Solve(AbstractOdeSystem* pOdeSystem,
00329 std::vector<double>& rYValues,
00330 double startTime,
00331 double endTime,
00332 double maxStep)
00333 {
00334 assert(endTime > startTime);
00335
00336 mStoppingEventOccurred = false;
00337 if (mCheckForRoots && pOdeSystem->CalculateStoppingEvent(startTime, rYValues) == true)
00338 {
00339 EXCEPTION("(Solve) Stopping event is true for initial condition");
00340 }
00341
00342 SetupCvode(pOdeSystem, rYValues, startTime, maxStep);
00343
00344 N_Vector yout = mInitialValues;
00345 double tend;
00346 int ierr;
00347 try
00348 {
00349 ierr = CVode(mpCvodeMem, endTime, yout, &tend, CV_NORMAL);
00350 }
00351 catch (...)
00352 {
00353 FreeCvodeMemory();
00354 throw;
00355 }
00356 if (ierr<0) CvodeError(ierr, "CVODE failed to solve system");
00357 if (ierr == CV_ROOT_RETURN)
00358 {
00359
00360 mStoppingEventOccurred = true;
00361 mStoppingTime = tend;
00362
00363 }
00364 assert(NV_DATA_S(yout) == &(rYValues[0]));
00365 assert(!NV_OWN_DATA_S(yout));
00366
00367
00368
00369
00370
00371 ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00372
00373
00374
00375
00376 assert(ierr == CV_SUCCESS);
00377 ierr=ierr;
00378 FreeCvodeMemory();
00379 }
00380
00381 CvodeAdaptor::CvodeAdaptor(double relTol, double absTol)
00382 : AbstractIvpOdeSolver(),
00383 mpCvodeMem(NULL), mInitialValues(NULL),
00384 mRelTol(relTol), mAbsTol(absTol),
00385 mLastInternalStepSize(-0.0),
00386 mMaxSteps(0),
00387 mCheckForRoots(false)
00388 {
00389 }
00390
00391 void CvodeAdaptor::SetTolerances(double relTol, double absTol)
00392 {
00393 mRelTol = relTol;
00394 mAbsTol = absTol;
00395 }
00396
00397 double CvodeAdaptor::GetRelativeTolerance()
00398 {
00399 return mRelTol;
00400 }
00401
00402 double CvodeAdaptor::GetAbsoluteTolerance()
00403 {
00404 return mAbsTol;
00405 }
00406
00407 double CvodeAdaptor::GetLastStepSize()
00408 {
00409 return mLastInternalStepSize;
00410 }
00411
00412 void CvodeAdaptor::CheckForStoppingEvents()
00413 {
00414 mCheckForRoots = true;
00415 }
00416
00417 void CvodeAdaptor::SetMaxSteps(long int numSteps)
00418 {
00419 mMaxSteps = numSteps;
00420 }
00421
00422 long int CvodeAdaptor::GetMaxSteps()
00423 {
00424 return mMaxSteps;
00425 }
00426
00427
00428 #include "SerializationExportWrapperForCpp.hpp"
00429 CHASTE_CLASS_EXPORT(CvodeAdaptor)
00430
00431 #endif // CHASTE_CVODE