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