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 solutions.SetOdeSystemInformation(pOdeSystem->GetSystemInformation());
00232
00233
00234 while (!stepper.IsTimeAtEnd() && !mStoppingEventOccurred)
00235 {
00236
00237 double tend;
00238 int ierr;
00239 try
00240 {
00241 ierr = CVode(mpCvodeMem, stepper.GetNextTime(), yout,
00242 &tend, CV_NORMAL);
00243 }
00244 catch (...)
00245 {
00246 FreeCvodeMemory();
00247 throw;
00248 }
00249 if (ierr<0) CvodeError(ierr, "CVODE failed to solve system");
00250
00251 solutions.rGetSolutions().push_back(rYValues);
00252 solutions.rGetTimes().push_back(tend);
00253 if (ierr == CV_ROOT_RETURN)
00254 {
00255
00256 mStoppingEventOccurred = true;
00257 mStoppingTime = tend;
00258 }
00259 stepper.AdvanceOneTimeStep();
00260 }
00261
00262
00263 solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
00264
00265
00266 int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00267 assert(ierr == CV_SUCCESS); ierr=ierr;
00268 FreeCvodeMemory();
00269
00270 return solutions;
00271 }
00272
00273
00274 void CvodeAdaptor::Solve(AbstractOdeSystem* pOdeSystem,
00275 std::vector<double>& rYValues,
00276 double startTime,
00277 double endTime,
00278 double maxStep)
00279 {
00280 assert(endTime > startTime);
00281
00282 mStoppingEventOccurred = false;
00283 if (mCheckForRoots && pOdeSystem->CalculateStoppingEvent(startTime, rYValues) == true)
00284 {
00285 EXCEPTION("(Solve) Stopping event is true for initial condition");
00286 }
00287
00288 SetupCvode(pOdeSystem, rYValues, startTime, maxStep);
00289
00290 N_Vector yout = mInitialValues;
00291 double tend;
00292 int ierr;
00293 try
00294 {
00295 ierr = CVode(mpCvodeMem, endTime, yout, &tend, CV_NORMAL);
00296 }
00297 catch (...)
00298 {
00299 FreeCvodeMemory();
00300 throw;
00301 }
00302 if (ierr<0) CvodeError(ierr, "CVODE failed to solve system");
00303 if (ierr == CV_ROOT_RETURN)
00304 {
00305
00306 mStoppingEventOccurred = true;
00307 mStoppingTime = tend;
00308
00309 }
00310 assert(NV_DATA_S(yout) == &(rYValues[0]));
00311 assert(!NV_OWN_DATA_S(yout));
00312
00313
00314
00315
00316
00317 ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00318
00319
00320
00321
00322 assert(ierr == CV_SUCCESS);
00323 ierr=ierr;
00324 FreeCvodeMemory();
00325 }
00326
00327 CvodeAdaptor::CvodeAdaptor(double relTol, double absTol)
00328 : AbstractIvpOdeSolver(),
00329 mpCvodeMem(NULL), mInitialValues(NULL),
00330 mRelTol(relTol), mAbsTol(absTol),
00331 mLastInternalStepSize(-0.0),
00332 mMaxSteps(0),
00333 mCheckForRoots(false)
00334 {
00335 }
00336
00337 void CvodeAdaptor::SetTolerances(double relTol, double absTol)
00338 {
00339 mRelTol = relTol;
00340 mAbsTol = absTol;
00341 }
00342
00343 double CvodeAdaptor::GetRelativeTolerance()
00344 {
00345 return mRelTol;
00346 }
00347
00348 double CvodeAdaptor::GetAbsoluteTolerance()
00349 {
00350 return mAbsTol;
00351 }
00352
00353 double CvodeAdaptor::GetLastStepSize()
00354 {
00355 return mLastInternalStepSize;
00356 }
00357
00358 void CvodeAdaptor::CheckForStoppingEvents()
00359 {
00360 mCheckForRoots = true;
00361 }
00362
00363 void CvodeAdaptor::SetMaxSteps(long int numSteps)
00364 {
00365 mMaxSteps = numSteps;
00366 }
00367
00368 long int CvodeAdaptor::GetMaxSteps()
00369 {
00370 return mMaxSteps;
00371 }
00372
00373
00374 #include "SerializationExportWrapperForCpp.hpp"
00375 CHASTE_CLASS_EXPORT(CvodeAdaptor);
00376
00377 #endif // CHASTE_CVODE