Chaste Release::3.1
CvodeAdaptor.cpp
00001 /*
00002 
00003 Copyright (c) 2005-2012, University of Oxford.
00004 All rights reserved.
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Redistribution and use in source and binary forms, with or without
00013 modification, are permitted provided that the following conditions are met:
00014  * Redistributions of source code must retain the above copyright notice,
00015    this list of conditions and the following disclaimer.
00016  * Redistributions in binary form must reproduce the above copyright notice,
00017    this list of conditions and the following disclaimer in the documentation
00018    and/or other materials provided with the distribution.
00019  * Neither the name of the University of Oxford nor the names of its
00020    contributors may be used to endorse or promote products derived from this
00021    software without specific prior written permission.
00022 
00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00033 
00034 */
00035 
00036 #ifdef CHASTE_CVODE
00037 
00038 #include "CvodeAdaptor.hpp"
00039 
00040 #include "Exception.hpp"
00041 #include "TimeStepper.hpp"
00042 #include "VectorHelperFunctions.hpp"
00043 #include "MathsCustomFunctions.hpp" // For tolerance check.
00044 
00045 #include <iostream>
00046 #include <sstream>
00047 
00048 // CVODE headers
00049 #include <sundials/sundials_nvector.h>
00050 #include <cvode/cvode_dense.h>
00051 
00052 
00069 int CvodeRhsAdaptor(realtype t, N_Vector y, N_Vector ydot, void* pData)
00070 {
00071     assert(pData != NULL);
00072     CvodeData* p_data = (CvodeData*) pData;
00073     // Get y, ydot into std::vector<>s
00074     static std::vector<realtype> ydot_vec;
00075     CopyToStdVector(y, *p_data->pY);
00076     CopyToStdVector(ydot, ydot_vec);
00077     // Call our function
00078     try
00079     {
00080         p_data->pSystem->EvaluateYDerivatives(t, *(p_data->pY), ydot_vec);
00081     }
00082     catch (const Exception &e)
00083     {
00084         std::cerr << "CVODE RHS Exception: " << e.GetMessage() << std::endl << std::flush;
00085         return -1;
00086     }
00087     // Copy derivative back
00088     CopyFromStdVector(ydot_vec, ydot);
00089     return 0;
00090 }
00091 
00112 int CvodeRootAdaptor(realtype t, N_Vector y, realtype* pGOut, void* pData)
00113 {
00114     assert(pData != NULL);
00115     CvodeData* p_data = (CvodeData*) pData;
00116     // Get y into a std::vector
00117     CopyToStdVector(y, *p_data->pY);
00118     // Call our function
00119     try
00120     {
00121         *pGOut = p_data->pSystem->CalculateRootFunction(t, *p_data->pY);
00122     }
00123     catch (const Exception &e)
00124     {
00125         std::cerr << "CVODE Root Exception: " << e.GetMessage() << std::endl << std::flush;
00126         return -1;
00127     }
00128     return 0;
00129 }
00130 
00131 // /**
00132 //  * Jacobian computation adaptor function.
00133 //  *
00134 //  * If solving an AbstractOdeSystemWithAnalyticJacobian, this function
00135 //  * can be used to allow CVODE to compute the Jacobian analytically.
00136 //  *
00137 //  * Note to self: can test using pSystem->GetUseAnalyticJacobian().
00138 //  */
00139 // int CvodeDenseJacobianAdaptor(long int numberOfStateVariables, DenseMat J,
00140 //                               realtype t, N_Vector y, N_Vector fy,
00141 //                               void* pData,
00142 //                               N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
00143 // {
00144 //     AbstractOdeSystemWithAnalyticJacobian* pSystem
00145 //         = (AbstractOdeSystemWithAnalyticJacobian*) pData;
00146 //     // Get the current time step
00147 //     double dt;
00148 //     CVodeGetCurrentStep(CvodeMem, &dt);
00149 //     // Get std::vector<> for y and double** for J
00150 //     std::vector<realtype>& y_vec = *NV_DATA_STL(y);
00151 //     double** ppJ = J->data; // organised column-wise: J_{i,j} = ppJ[j][i]
00152 //     // Call our function
00153 //     try
00154 //     {
00155 //         pSystem->AnalyticJacobian(y_vec, ppJ, t, dt);
00156 //     }
00157 //     catch (const Exception &e)
00158 //     {
00159 //         std::cerr << "CVODE J Exception: " << e.GetMessage() << std::endl << std::flush;
00160 //         return -1;
00161 //     }
00162 //     // Update J (if needed)
00163 //     return 0;
00164 // }
00165 
00166 
00167 
00168 void CvodeErrorHandler(int errorCode, const char *module, const char *function,
00169                        char *message, void* pData)
00170 {
00171     std::stringstream err;
00172     err << "CVODE Error " << errorCode << " in module " << module
00173         << " function " << function << ": " << message;
00174     std::cerr << "*" << err.str() << std::endl << std::flush;
00175     // Throwing the exception here causes termination on Maverick (g++ 4.4)
00176     //EXCEPTION(err.str());
00177 }
00178 
00179 
00180 
00181 void CvodeAdaptor::SetupCvode(AbstractOdeSystem* pOdeSystem,
00182                               std::vector<double>& rInitialY,
00183                               double startTime,
00184                               double maxStep)
00185 {
00186     assert(rInitialY.size() == pOdeSystem->GetNumberOfStateVariables());
00187     assert(maxStep > 0.0);
00188 
00190     N_Vector initial_values = N_VMake_Serial(rInitialY.size(), &(rInitialY[0]));
00191     assert(NV_DATA_S(initial_values) == &(rInitialY[0]));
00192     assert(!NV_OWN_DATA_S(initial_values));
00193 //    std::cout << " Initial values: "; N_VPrint_Serial(initial_values);
00194 //    std::cout << " Rtol: " << mRelTol << ", Atol: " << mAbsTol << std::endl;
00195 //    std::cout << " Start: " << startTime << " max dt=" << maxStep << std::endl << std::flush;
00196 
00197     // Find out if we need to (re-)initialise
00198     bool reinit = !mpCvodeMem || mAutoReset || !mLastSolutionState || !CompareDoubles::WithinAnyTolerance(startTime, mLastSolutionTime);
00199     if (!reinit && !mForceMinimalReset)
00200     {
00201         const unsigned size = GetVectorSize(rInitialY);
00202         for (unsigned i=0; i<size; i++)
00203         {
00204             if (!CompareDoubles::WithinAnyTolerance(GetVectorComponent(mLastSolutionState, i), GetVectorComponent(rInitialY, i)))
00205             {
00206                 reinit = true;
00207                 break;
00208             }
00209         }
00210     }
00211 
00212     if (!mpCvodeMem) // First run of this solver, set up CVODE memory
00213     {
00214         // Set up CVODE's memory.
00215         mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
00216         if (mpCvodeMem == NULL) EXCEPTION("Failed to SetupCvode CVODE"); // In one line to avoid coverage problem!
00217         // Set error handler
00218         CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, NULL);
00219         // Set the user data
00220         mData.pSystem = pOdeSystem;
00221         mData.pY = &rInitialY;
00222         #if CHASTE_SUNDIALS_VERSION >= 20400
00223             CVodeSetUserData(mpCvodeMem, (void*)(&mData));
00224         #else
00225             CVodeSetFdata(mpCvodeMem, (void*)(&mData));
00226         #endif
00227 
00228         // Setup CVODE
00229         #if CHASTE_SUNDIALS_VERSION >= 20400
00230             CVodeInit(mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values);
00231             CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
00232         #else
00233             CVodeMalloc(mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values,
00234                         CV_SS, mRelTol, &mAbsTol);
00235         #endif
00236 
00237         // Set the rootfinder function if wanted
00238         if (mCheckForRoots)
00239         {
00240             #if CHASTE_SUNDIALS_VERSION >= 20400
00241                 CVodeRootInit(mpCvodeMem, 1, CvodeRootAdaptor);
00242             #else
00243                 CVodeRootInit(mpCvodeMem, 1, CvodeRootAdaptor, (void*)(&mData));
00244             #endif
00245         }
00246         // Attach a linear solver for Newton iteration
00247         CVDense(mpCvodeMem, rInitialY.size());
00248     }
00249     else if (reinit) // Could be new ODE system, or new Y values
00250     {
00251         // Set the user data
00252         mData.pSystem = pOdeSystem; // stays the same on a re-initialize
00253         mData.pY = &rInitialY; // changes on a re-initialize
00254         #if CHASTE_SUNDIALS_VERSION >= 20400
00255             CVodeSetUserData(mpCvodeMem, (void*)(&mData));
00256         #else
00257             CVodeSetFdata(mpCvodeMem, (void*)(&mData));
00258         #endif
00259 
00260         #if CHASTE_SUNDIALS_VERSION >= 20400
00261             CVodeReInit(mpCvodeMem, startTime, initial_values);
00262             CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
00263         #else
00264             CVodeReInit(mpCvodeMem, CvodeRhsAdaptor, startTime, initial_values,
00265                         CV_SS, mRelTol, &mAbsTol);
00266         #endif
00267 
00268         // Attach a linear solver for Newton iteration
00269         CVDense(mpCvodeMem, rInitialY.size());
00270     }
00271 
00272     CVodeSetMaxStep(mpCvodeMem, maxStep);
00273     // Change max steps if wanted
00274     if (mMaxSteps > 0)
00275     {
00276         CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
00277         CVodeSetMaxErrTestFails(mpCvodeMem, 15);
00278     }
00279     DeleteVector(initial_values);
00280 }
00281 
00282 void CvodeAdaptor::FreeCvodeMemory()
00283 {
00284     if (mpCvodeMem)
00285     {
00286         CVodeFree(&mpCvodeMem);
00287     }
00288     mpCvodeMem = NULL;
00289 }
00290 
00291 
00292 void CvodeAdaptor::SetAutoReset(bool autoReset)
00293 {
00294     mAutoReset = autoReset;
00295     if (mAutoReset)
00296     {
00297         SetMinimalReset(false);
00298         ResetSolver();
00299     }
00300 }
00301 
00302 void CvodeAdaptor::SetMinimalReset(bool minimalReset)
00303 {
00304     mForceMinimalReset = minimalReset;
00305     if (mForceMinimalReset)
00306     {
00307         SetAutoReset(false);
00308     }
00309 }
00310 
00311 void CvodeAdaptor::ResetSolver()
00312 {
00313     // Doing this makes the SetupCvode think it needs to reset things.
00314     DeleteVector(mLastSolutionState);
00315 }
00316 
00317 void CvodeAdaptor::CvodeError(int flag, const char * msg)
00318 {
00319     std::stringstream err;
00320     char* p_flag_name = CVodeGetReturnFlagName(flag);
00321     err << msg << ": " << p_flag_name;
00322     free(p_flag_name);
00323     std::cerr << err.str() << std::endl << std::flush;
00324     EXCEPTION(err.str());
00325 }
00326 
00327 
00328 OdeSolution CvodeAdaptor::Solve(AbstractOdeSystem* pOdeSystem,
00329                                 std::vector<double>& rYValues,
00330                                 double startTime,
00331                                 double endTime,
00332                                 double maxStep,
00333                                 double timeSampling)
00334 {
00335     assert(endTime > startTime);
00336     assert(timeSampling > 0.0);
00337 
00338     mStoppingEventOccurred = false;
00339     if (mCheckForRoots && pOdeSystem->CalculateStoppingEvent(startTime, rYValues) == true)
00340     {
00341         EXCEPTION("(Solve with sampling) Stopping event is true for initial condition");
00342     }
00343 
00344     SetupCvode(pOdeSystem, rYValues, startTime, maxStep);
00345 
00346     TimeStepper stepper(startTime, endTime, timeSampling);
00347     N_Vector yout = N_VMake_Serial(rYValues.size(), &(rYValues[0]));
00348 
00349     // Set up ODE solution
00350     OdeSolution solutions;
00351     solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
00352     solutions.rGetSolutions().push_back(rYValues);
00353     solutions.rGetTimes().push_back(startTime);
00354     solutions.SetOdeSystemInformation(pOdeSystem->GetSystemInformation());
00355 
00356     // Main time sampling loop
00357     while (!stepper.IsTimeAtEnd() && !mStoppingEventOccurred)
00358     {
00359         double tend;
00360         int ierr = CVode(mpCvodeMem, stepper.GetNextTime(), yout, &tend, CV_NORMAL);
00361         if (ierr<0)
00362         {
00363             FreeCvodeMemory();
00364             DeleteVector(yout);
00365             CvodeError(ierr, "CVODE failed to solve system");
00366         }
00367         // Store solution
00368         solutions.rGetSolutions().push_back(rYValues);
00369         solutions.rGetTimes().push_back(tend);
00370         if (ierr == CV_ROOT_RETURN)
00371         {
00372             // Stopping event occurred
00373             mStoppingEventOccurred = true;
00374             mStoppingTime = tend;
00375         }
00376         stepper.AdvanceOneTimeStep();
00377     }
00378 
00379     // stepper.EstimateTimeSteps may have been an overestimate...
00380     solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
00381 
00382     int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00383     assert(ierr == CV_SUCCESS); ierr=ierr; // avoid unused var warning
00384     RecordStoppingPoint(endTime, yout);
00385     DeleteVector(yout);
00386 
00387     return solutions;
00388 }
00389 
00390 
00391 void CvodeAdaptor::Solve(AbstractOdeSystem* pOdeSystem,
00392                          std::vector<double>& rYValues,
00393                          double startTime,
00394                          double endTime,
00395                          double maxStep)
00396 {
00397     assert(endTime > startTime);
00398 
00399     mStoppingEventOccurred = false;
00400     if (mCheckForRoots && pOdeSystem->CalculateStoppingEvent(startTime, rYValues) == true)
00401     {
00402         EXCEPTION("(Solve) Stopping event is true for initial condition");
00403     }
00404 
00405     SetupCvode(pOdeSystem, rYValues, startTime, maxStep);
00406 
00407     N_Vector yout = N_VMake_Serial(rYValues.size(), &(rYValues[0]));
00408     double tend;
00409     int ierr = CVode(mpCvodeMem, endTime, yout, &tend, CV_NORMAL);
00410     if (ierr<0)
00411     {
00412         FreeCvodeMemory();
00413         DeleteVector(yout);
00414         CvodeError(ierr, "CVODE failed to solve system");
00415     }
00416     if (ierr == CV_ROOT_RETURN)
00417     {
00418         // Stopping event occurred
00419         mStoppingEventOccurred = true;
00420         mStoppingTime = tend;
00421     }
00422     assert(NV_DATA_S(yout) == &(rYValues[0]));
00423     assert(!NV_OWN_DATA_S(yout));
00424 
00425 //    long int steps;
00426 //    CVodeGetNumSteps(mpCvodeMem, &steps);
00427 //    std::cout << " Solved to " << endTime << " in " << steps << " steps.\n";
00428 
00429     ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00430     assert(ierr == CV_SUCCESS); ierr=ierr; // avoid unused var warning
00431     RecordStoppingPoint(tend, yout);
00432     DeleteVector(yout);
00433 }
00434 
00435 CvodeAdaptor::CvodeAdaptor(double relTol, double absTol)
00436     : AbstractIvpOdeSolver(),
00437       mpCvodeMem(NULL),
00438       mRelTol(relTol),
00439       mAbsTol(absTol),
00440       mLastInternalStepSize(-0.0),
00441       mMaxSteps(0),
00442       mCheckForRoots(false),
00443       mLastSolutionState(NULL),
00444       mLastSolutionTime(0.0),
00445       mAutoReset(true),
00446       mForceMinimalReset(false)
00447 {
00448 }
00449 
00450 CvodeAdaptor::~CvodeAdaptor()
00451 {
00452     FreeCvodeMemory();
00453     DeleteVector(mLastSolutionState);
00454 }
00455 
00456 void CvodeAdaptor::RecordStoppingPoint(double stopTime, N_Vector yEnd)
00457 {
00458     if (!mAutoReset)
00459     {
00460         const unsigned size = GetVectorSize(yEnd);
00461         CreateVectorIfEmpty(mLastSolutionState, size);
00462         for (unsigned i=0; i<size; i++)
00463         {
00464             SetVectorComponent(mLastSolutionState, i, GetVectorComponent(yEnd, i));
00465         }
00466         mLastSolutionTime = stopTime;
00467     }
00468 }
00469 
00470 void CvodeAdaptor::SetTolerances(double relTol, double absTol)
00471 {
00472     mRelTol = relTol;
00473     mAbsTol = absTol;
00474 }
00475 
00476 double CvodeAdaptor::GetRelativeTolerance()
00477 {
00478     return mRelTol;
00479 }
00480 
00481 double CvodeAdaptor::GetAbsoluteTolerance()
00482 {
00483     return mAbsTol;
00484 }
00485 
00486 double CvodeAdaptor::GetLastStepSize()
00487 {
00488     return mLastInternalStepSize;
00489 }
00490 
00491 void CvodeAdaptor::CheckForStoppingEvents()
00492 {
00493     mCheckForRoots = true;
00494 }
00495 
00496 void CvodeAdaptor::SetMaxSteps(long int numSteps)
00497 {
00498     mMaxSteps = numSteps;
00499 }
00500 
00501 long int CvodeAdaptor::GetMaxSteps()
00502 {
00503     return mMaxSteps;
00504 }
00505 
00506 // Serialization for Boost >= 1.36
00507 #include "SerializationExportWrapperForCpp.hpp"
00508 CHASTE_CLASS_EXPORT(CvodeAdaptor)
00509 
00510 #endif // CHASTE_CVODE