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 #include "AbstractIvpOdeSolver.hpp"
00033 #include "TimeStepper.hpp"
00034 
00035 #include <iostream>
00036 
00037 // CVODE headers
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     // Check for no-op
00049     if (NV_DATA_S(src) == &(rDest[0])) return;
00050     // Set dest size
00051     rDest.resize(NV_LENGTH_S(src));
00052     // Copy data
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     // Get y, ydot into std::vector<>s
00064     static std::vector<realtype> ydot_vec;
00065     CopyToStdVector(y, *p_data->pY);
00066     CopyToStdVector(ydot, ydot_vec);
00067     // Call our function
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     // Copy derivative back
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     // Get y into a std::vector
00088     CopyToStdVector(y, *p_data->pY);
00089     // Call our function
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 // int CvodeDenseJacobianAdaptor(long int numberOfStateVariables, DenseMat J,
00103 //                               realtype t, N_Vector y, N_Vector fy,
00104 //                               void* pData,
00105 //                               N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
00106 // {
00107 //     AbstractOdeSystemWithAnalyticJacobian* pSystem
00108 //         = (AbstractOdeSystemWithAnalyticJacobian*) pData;
00109 //     // Get the current time step
00110 //     double dt;
00111 //     CVodeGetCurrentStep(CvodeMem, &dt);
00112 //     // Get std::vector<> for y and double** for J
00113 //     std::vector<realtype>& y_vec = *NV_DATA_STL(y);
00114 //     double** ppJ = J->data; // organised column-wise: J_{i,j} = ppJ[j][i]
00115 //     // Call our function
00116 //     try
00117 //     {
00118 //         pSystem->AnalyticJacobian(y_vec, ppJ, t, dt);
00119 //     }
00120 //     catch (Exception &e)
00121 //     {
00122 //         std::cerr << "CVODE J Exception: " << e.GetMessage() << std::endl << std::flush;
00123 //         return -1;
00124 //     }
00125 //     // Update J (if needed)
00126 //     return 0;
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 //    std::cout << " Initial values: "; N_VPrint_Serial(mInitialValues);
00153 //    std::cout << " Rtol: " << mRelTol << ", Atol: " << mAbsTol << std::endl;
00154 //    std::cout << " Start: " << startTime << " max dt=" << maxStep << std::endl << std::flush;
00155 
00156     mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
00157     if (mpCvodeMem == NULL) EXCEPTION("Failed to SetupCvode CVODE");
00158     // Set error handler
00159     CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, NULL);
00160     // Set the user data
00161     mData.pSystem = pOdeSystem;
00162     mData.pY = &rInitialY;
00163     CVodeSetFdata(mpCvodeMem, (void*)(&mData));
00164     // Setup CVODE
00165     CVodeMalloc(mpCvodeMem, CvodeRhsAdaptor, startTime, mInitialValues,
00166                 CV_SS, mRelTol, &mAbsTol);
00167     CVodeSetMaxStep(mpCvodeMem, maxStep);
00168     // Set the rootfinder function if wanted
00169     if (mCheckForRoots)
00170     {
00171         CVodeRootInit(mpCvodeMem, 1, CvodeRootAdaptor, (void*)(&mData));
00172     }
00173     // Change max steps if wanted
00174     if (mMaxSteps > 0)
00175     {
00176         CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
00177     }
00178     // Attach a linear solver for Newton iteration
00179     CVDense(mpCvodeMem, rInitialY.size());
00180 }
00181 
00182 void CvodeAdaptor::FreeCvodeMemory()
00183 {
00184 //    assert(!NV_OWN_DATA_STL(mInitialValues));
00185 //    std::vector<double>* pVec = NV_DATA_STL(mInitialValues);
00186 //    double val = (*pVec)[0];
00187     N_VDestroy_Serial(mInitialValues); mInitialValues = NULL;
00188 //    std::cout << "  a: " << val << ", b: " << (*pVec)[0] << std::endl;
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     // Set up ODE solution
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     // Main time sampling loop
00234     while (!stepper.IsTimeAtEnd() && !mStoppingEventOccurred)
00235     {
00236 //        std::cout << "Solving to time " << stepper.GetNextTime() << std::endl << std::flush;
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         // Store solution
00251         solutions.rGetSolutions().push_back(rYValues);
00252         solutions.rGetTimes().push_back(tend);
00253         if (ierr == CV_ROOT_RETURN)
00254         {
00255             // Stopping event occurred
00256             mStoppingEventOccurred = true;
00257             mStoppingTime = tend;
00258         }
00259         stepper.AdvanceOneTimeStep();
00260     }
00261 
00262     // stepper.EstimateTimeSteps may have been an overestimate...
00263     solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
00264 
00265 //    std::cout << " Solved to " << stepper.GetTime() << " in " << stepper.GetTotalTimeStepsTaken() << " samples.\n" << std::flush;
00266     int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00267     assert(ierr == CV_SUCCESS); ierr=ierr; // avoid unused var warning
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         // Stopping event occurred
00306         mStoppingEventOccurred = true;
00307         mStoppingTime = tend;
00308 //        std::cout << "CVODE Stopped at t = " << tend << std::endl;
00309     }
00310     assert(NV_DATA_S(yout) == &(rYValues[0]));
00311     assert(!NV_OWN_DATA_S(yout));
00312 
00313 //    long int steps;
00314 //    CVodeGetNumSteps(mpCvodeMem, &steps);
00315 //    std::cout << " Solved to " << endTime << " in " << steps << " steps.\n";
00316 
00317     ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00318 //    if (mStoppingEventOccurred)
00319 //    {
00320 //        std::cout << "Last internal dt was " << mLastInternalStepSize << std::endl;
00321 //    }
00322     assert(ierr == CV_SUCCESS);
00323     ierr=ierr; // avoid unused var warning
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 // Serialization for Boost >= 1.36
00374 #include "SerializationExportWrapperForCpp.hpp"
00375 CHASTE_CLASS_EXPORT(CvodeAdaptor);
00376 
00377 #endif // CHASTE_CVODE

Generated by  doxygen 1.6.2