CvodeAdaptor.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2010
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
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 // CVODE headers
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     // Check for no-op
00052     realtype* p_src = NV_DATA_S(src);
00053     if (p_src == &(rDest[0])) return;
00054     // Set dest size
00055     long size = NV_LENGTH_S(src);
00056     rDest.resize(size);
00057     // Copy data
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     // Get y, ydot into std::vector<>s
00085     static std::vector<realtype> ydot_vec;
00086     CopyToStdVector(y, *p_data->pY);
00087     CopyToStdVector(ydot, ydot_vec);
00088     // Call our function
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     // Copy derivative back
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     // Get y into a std::vector
00133     CopyToStdVector(y, *p_data->pY);
00134     // Call our function
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 //  * Jacobian computation adaptor function.
00149 //  *
00150 //  * If solving an AbstractOdeSystemWithAnalyticJacobian, this function
00151 //  * can be used to allow CVODE to compute the Jacobian analytically.
00152 //  *
00153 //  * Note to self: can test using pSystem->GetUseAnalyticJacobian().
00154 //  */
00155 // int CvodeDenseJacobianAdaptor(long int numberOfStateVariables, DenseMat J,
00156 //                               realtype t, N_Vector y, N_Vector fy,
00157 //                               void* pData,
00158 //                               N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
00159 // {
00160 //     AbstractOdeSystemWithAnalyticJacobian* pSystem
00161 //         = (AbstractOdeSystemWithAnalyticJacobian*) pData;
00162 //     // Get the current time step
00163 //     double dt;
00164 //     CVodeGetCurrentStep(CvodeMem, &dt);
00165 //     // Get std::vector<> for y and double** for J
00166 //     std::vector<realtype>& y_vec = *NV_DATA_STL(y);
00167 //     double** ppJ = J->data; // organised column-wise: J_{i,j} = ppJ[j][i]
00168 //     // Call our function
00169 //     try
00170 //     {
00171 //         pSystem->AnalyticJacobian(y_vec, ppJ, t, dt);
00172 //     }
00173 //     catch (Exception &e)
00174 //     {
00175 //         std::cerr << "CVODE J Exception: " << e.GetMessage() << std::endl << std::flush;
00176 //         return -1;
00177 //     }
00178 //     // Update J (if needed)
00179 //     return 0;
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 //    std::cout << " Initial values: "; N_VPrint_Serial(mInitialValues);
00207 //    std::cout << " Rtol: " << mRelTol << ", Atol: " << mAbsTol << std::endl;
00208 //    std::cout << " Start: " << startTime << " max dt=" << maxStep << std::endl << std::flush;
00209 
00210     mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
00211     if (mpCvodeMem == NULL) EXCEPTION("Failed to SetupCvode CVODE");
00212     // Set error handler
00213     CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, NULL);
00214     // Set the user data
00215     mData.pSystem = pOdeSystem;
00216     mData.pY = &rInitialY;
00217     CVodeSetFdata(mpCvodeMem, (void*)(&mData));
00218     // Setup CVODE
00219     CVodeMalloc(mpCvodeMem, CvodeRhsAdaptor, startTime, mInitialValues,
00220                 CV_SS, mRelTol, &mAbsTol);
00221     CVodeSetMaxStep(mpCvodeMem, maxStep);
00222     // Set the rootfinder function if wanted
00223     if (mCheckForRoots)
00224     {
00225         CVodeRootInit(mpCvodeMem, 1, CvodeRootAdaptor, (void*)(&mData));
00226     }
00227     // Change max steps if wanted
00228     if (mMaxSteps > 0)
00229     {
00230         CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
00231     }
00232     // Attach a linear solver for Newton iteration
00233     CVDense(mpCvodeMem, rInitialY.size());
00234 }
00235 
00236 void CvodeAdaptor::FreeCvodeMemory()
00237 {
00238 //    assert(!NV_OWN_DATA_STL(mInitialValues));
00239 //    std::vector<double>* pVec = NV_DATA_STL(mInitialValues);
00240 //    double val = (*pVec)[0];
00241     N_VDestroy_Serial(mInitialValues); mInitialValues = NULL;
00242 //    std::cout << "  a: " << val << ", b: " << (*pVec)[0] << std::endl;
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     // Set up ODE solution
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     // Main time sampling loop
00288     while (!stepper.IsTimeAtEnd() && !mStoppingEventOccurred)
00289     {
00290 //        std::cout << "Solving to time " << stepper.GetNextTime() << std::endl << std::flush;
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         // Store solution
00305         solutions.rGetSolutions().push_back(rYValues);
00306         solutions.rGetTimes().push_back(tend);
00307         if (ierr == CV_ROOT_RETURN)
00308         {
00309             // Stopping event occurred
00310             mStoppingEventOccurred = true;
00311             mStoppingTime = tend;
00312         }
00313         stepper.AdvanceOneTimeStep();
00314     }
00315 
00316     // stepper.EstimateTimeSteps may have been an overestimate...
00317     solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
00318 
00319 //    std::cout << " Solved to " << stepper.GetTime() << " in " << stepper.GetTotalTimeStepsTaken() << " samples.\n" << std::flush;
00320     int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00321     assert(ierr == CV_SUCCESS); ierr=ierr; // avoid unused var warning
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         // Stopping event occurred
00360         mStoppingEventOccurred = true;
00361         mStoppingTime = tend;
00362 //        std::cout << "CVODE Stopped at t = " << tend << std::endl;
00363     }
00364     assert(NV_DATA_S(yout) == &(rYValues[0]));
00365     assert(!NV_OWN_DATA_S(yout));
00366 
00367 //    long int steps;
00368 //    CVodeGetNumSteps(mpCvodeMem, &steps);
00369 //    std::cout << " Solved to " << endTime << " in " << steps << " steps.\n";
00370 
00371     ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00372 //    if (mStoppingEventOccurred)
00373 //    {
00374 //        std::cout << "Last internal dt was " << mLastInternalStepSize << std::endl;
00375 //    }
00376     assert(ierr == CV_SUCCESS);
00377     ierr=ierr; // avoid unused var warning
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 // Serialization for Boost >= 1.36
00428 #include "SerializationExportWrapperForCpp.hpp"
00429 CHASTE_CLASS_EXPORT(CvodeAdaptor)
00430 
00431 #endif // CHASTE_CVODE

Generated on Mon Nov 1 12:35:24 2010 for Chaste by  doxygen 1.5.5