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