CvodeAdaptor.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2009
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 <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 // CVODE headers
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     // Check for no-op
00054     if (NV_DATA_S(src) == &(rDest[0])) return;
00055     // Set dest size
00056     rDest.resize(NV_LENGTH_S(src));
00057     // Copy data
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     // Get y, ydot into std::vector<>s
00069     static std::vector<realtype> ydot_vec;
00070     CopyToStdVector(y, * p_data->pY);
00071     CopyToStdVector(ydot, ydot_vec);
00072     // Call our function
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     // Copy derivative back
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     // Get y into a std::vector
00093     CopyToStdVector(y, * p_data->pY);
00094     // Call our function
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 // int CvodeDenseJacobianAdaptor(long int numberOfStateVariables, DenseMat J,
00108 //                               realtype t, N_Vector y, N_Vector fy,
00109 //                               void *pData,
00110 //                               N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
00111 // {
00112 //     AbstractOdeSystemWithAnalyticJacobian* pSystem
00113 //         = (AbstractOdeSystemWithAnalyticJacobian*) pData;
00114 //     // Get the current time step
00115 //     double dt;
00116 //     CVodeGetCurrentStep(CvodeMem, &dt);
00117 //     // Get std::vector<> for y and double** for J
00118 //     std::vector<realtype>& y_vec = *NV_DATA_STL(y);
00119 //     double** ppJ = J->data; // organised column-wise: J_{i,j} = ppJ[j][i]
00120 //     // Call our function
00121 //     try
00122 //     {
00123 //         pSystem->AnalyticJacobian(y_vec, ppJ, t, dt);
00124 //     }
00125 //     catch (Exception &e)
00126 //     {
00127 //         std::cerr << "CVODE J Exception: " << e.GetMessage() << std::endl << std::flush;
00128 //         return -1;
00129 //     }
00130 //     // Update J (if needed)
00131 //     return 0;
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 //    std::cout << " Initial values: "; N_VPrint_Serial(mInitialValues);
00158 //    std::cout << " Rtol: " << mRelTol << ", Atol: " << mAbsTol << std::endl;
00159 //    std::cout << " Start: " << startTime << " max dt=" << maxStep << std::endl << std::flush;
00160 
00161     mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
00162     if (mpCvodeMem == NULL) EXCEPTION("Failed to SetupCvode CVODE");
00163     // Set error handler
00164     CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, NULL);
00165     // Set the user data
00166     mData.pSystem = pOdeSystem;
00167     mData.pY = &rInitialY;
00168     CVodeSetFdata(mpCvodeMem, (void*)(&mData));
00169     // Setup CVODE
00170     CVodeMalloc(mpCvodeMem, CvodeRhsAdaptor, startTime, mInitialValues,
00171                 CV_SS, mRelTol, &mAbsTol);
00172     CVodeSetMaxStep(mpCvodeMem, maxStep);
00173     // Set the rootfinder function if wanted
00174     if (mCheckForRoots)
00175     {
00176         CVodeRootInit(mpCvodeMem, 1, CvodeRootAdaptor, (void*)(&mData));
00177     }
00178     // Change max steps if wanted
00179     if (mMaxSteps > 0)
00180     {
00181         CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
00182     }
00183     // Attach a linear solver for Newton iteration
00184     CVDense(mpCvodeMem, rInitialY.size());
00185 }
00186 
00187 void CvodeAdaptor::FreeCvodeMemory()
00188 {
00189 //    assert(!NV_OWN_DATA_STL(mInitialValues));
00190 //    std::vector<double>* pVec = NV_DATA_STL(mInitialValues);
00191 //    double val = (*pVec)[0];
00192     N_VDestroy_Serial(mInitialValues); mInitialValues = NULL;
00193 //    std::cout << "  a: " << val << ", b: " << (*pVec)[0] << std::endl;
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     // Set up ode solution
00232     OdeSolution solutions;
00233     solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
00234     solutions.rGetSolutions().push_back(rYValues);
00235     solutions.rGetTimes().push_back(startTime);
00236 
00237     // Main time sampling loop
00238     while (!stepper.IsTimeAtEnd() && !mStoppingEventOccurred)
00239     {
00240 //        std::cout << "Solving to time " << stepper.GetNextTime() << std::endl << std::flush;
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         // Store solution
00255         solutions.rGetSolutions().push_back(rYValues);
00256         solutions.rGetTimes().push_back(tend);
00257         if (ierr == CV_ROOT_RETURN)
00258         {
00259             // Stopping event occurred
00260             mStoppingEventOccurred = true;
00261             mStoppingTime = tend;
00262         }
00263         stepper.AdvanceOneTimeStep();
00264     }
00265 
00266     // stepper.EstimateTimeSteps may have been an overestimate...
00267     solutions.SetNumberOfTimeSteps(stepper.GetTimeStepsElapsed());
00268 
00269 //    std::cout << " Solved to " << stepper.GetTime() << " in " << stepper.GetTimeStepsElapsed() << " samples.\n" << std::flush;
00270     int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00271     assert(ierr == CV_SUCCESS); ierr=ierr; // avoid unused var warning
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         // Stopping event occurred
00310         mStoppingEventOccurred = true;
00311         mStoppingTime = tend;
00312 //        std::cout << "CVODE Stopped at t = " << tend << std::endl;
00313     }
00314     assert(NV_DATA_S(yout) == &(rYValues[0]));
00315     assert(!NV_OWN_DATA_S(yout));
00316 
00317 //    long int steps;
00318 //    CVodeGetNumSteps(mpCvodeMem, &steps);
00319 //    std::cout << " Solved to " << endTime << " in " << steps << " steps.\n";
00320 
00321     ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00322 //    if (mStoppingEventOccurred)
00323 //    {
00324 //        std::cout << "Last internal dt was " << mLastInternalStepSize << std::endl;
00325 //    }
00326     assert(ierr == CV_SUCCESS); ierr=ierr; // avoid unused var warning
00327     FreeCvodeMemory();
00328 }
00329 
00330 #endif // CHASTE_CVODE

Generated on Wed Mar 18 12:51:56 2009 for Chaste by  doxygen 1.5.5