AbstractCvodeSystem.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2011
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 <sstream>
00032 #include <cassert>
00033 
00034 #include "AbstractCvodeSystem.hpp"
00035 #include "Exception.hpp"
00036 #include "VectorHelperFunctions.hpp"
00037 #include "TimeStepper.hpp"
00038 #include "CvodeAdaptor.hpp" // For CvodeErrorHandler
00039 
00040 // CVODE headers
00041 #include <cvode/cvode.h>
00042 #include <sundials/sundials_nvector.h>
00043 #include <cvode/cvode_dense.h>
00044 
00045 
00055 int AbstractCvodeSystemRhsAdaptor(realtype t, N_Vector y, N_Vector ydot, void *pData)
00056 {
00057     assert(pData != NULL);
00058     AbstractCvodeSystem* p_ode_system = (AbstractCvodeSystem*) pData;
00059     try
00060     {
00061         p_ode_system->EvaluateYDerivatives(t, y, ydot);
00062     }
00063     catch (const Exception &e)
00064     {
00065         std::cerr << "CVODE RHS Exception: " << e.GetMessage()
00066                   << std::endl << std::flush;
00067         return -1;
00068     }
00069     return 0;
00070 }
00071 
00072 AbstractCvodeSystem::AbstractCvodeSystem(unsigned numberOfStateVariables)
00073     : AbstractParameterisedSystem<N_Vector>(numberOfStateVariables),
00074       mLastSolutionState(NULL),
00075       mLastSolutionTime(0.0),
00076       mAutoReset(true),
00077       mUseAnalyticJacobian(false),
00078       mpCvodeMem(NULL),
00079       mMaxSteps(0)
00080 {
00081     SetTolerances();
00082 }
00083 
00084 void AbstractCvodeSystem::Init()
00085 {
00086     DeleteVector(mStateVariables);
00087     mStateVariables = GetInitialConditions();
00088     DeleteVector(mParameters);
00089     mParameters = N_VNew_Serial(rGetParameterNames().size());
00090     for (int i=0; i<NV_LENGTH_S(mParameters); i++)
00091     {
00092         NV_Ith_S(mParameters, i) = 0.0;
00093     }
00094 }
00095 
00096 AbstractCvodeSystem::~AbstractCvodeSystem()
00097 {
00098     FreeCvodeMemory();
00099     DeleteVector(mStateVariables);
00100     DeleteVector(mParameters);
00101     DeleteVector(mLastSolutionState);
00102 }
00103 
00104 //
00105 //double AbstractCvodeSystem::CalculateRootFunction(double time, const std::vector<double>& rY)
00106 //{
00107 //    bool stop = CalculateStoppingEvent(time, rY);
00108 //    return stop ? 0.0 : 1.0;
00109 //}
00110 //
00111 //bool AbstractCvodeSystem::GetUseAnalyticJacobian()
00112 //{
00113 //    return mUseAnalyticJacobian;
00114 //}
00115 
00116 
00117 
00118 OdeSolution AbstractCvodeSystem::Solve(realtype tStart,
00119                                        realtype tEnd,
00120                                        realtype maxDt,
00121                                        realtype tSamp)
00122 {
00123     assert(tEnd >= tStart);
00124     assert(tSamp > 0.0);
00125 
00126     SetupCvode(mStateVariables, tStart, maxDt);
00127 
00128     TimeStepper stepper(tStart, tEnd, tSamp);
00129 
00130     // Set up ODE solution
00131     OdeSolution solutions;
00132     solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
00133     solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
00134     solutions.rGetTimes().push_back(tStart);
00135     solutions.SetOdeSystemInformation(mpSystemInfo);
00136 
00137     // Main time sampling loop
00138     while (!stepper.IsTimeAtEnd())
00139     {
00140         double cvode_stopped_at;
00141         int ierr = CVode(mpCvodeMem, stepper.GetNextTime(), mStateVariables,
00142                          &cvode_stopped_at, CV_NORMAL);
00143         if (ierr<0)
00144         {
00145             FreeCvodeMemory();
00146             CvodeError(ierr, "CVODE failed to solve system");
00147         }
00148         // Not root finding, so should have reached requested time
00149         assert(fabs(cvode_stopped_at - stepper.GetNextTime()) < DBL_EPSILON);
00150 
00151         VerifyStateVariables();
00152 
00153         // Store solution
00154         solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
00155         solutions.rGetTimes().push_back(cvode_stopped_at);
00156         stepper.AdvanceOneTimeStep();
00157     }
00158 
00159     // stepper.EstimateTimeSteps may have been an overestimate...
00160     solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
00161 
00162     int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00163     assert(ierr == CV_SUCCESS); ierr=ierr; // avoid unused var warning
00164 
00165     RecordStoppingPoint(tEnd);
00166 
00167     return solutions;
00168 }
00169 
00170 void AbstractCvodeSystem::Solve(realtype tStart,
00171                                 realtype tEnd,
00172                                 realtype maxDt)
00173 {
00174     assert(tEnd >= tStart);
00175 
00176     SetupCvode(mStateVariables, tStart, maxDt);
00177 
00178     double cvode_stopped_at;
00179     int ierr = CVode(mpCvodeMem, tEnd, mStateVariables, &cvode_stopped_at, CV_NORMAL);
00180     if (ierr<0)
00181     {
00182         FreeCvodeMemory();
00183         CvodeError(ierr, "CVODE failed to solve system");
00184     }
00185     // Not root finding, so should have reached requested time
00186     assert(fabs(cvode_stopped_at - tEnd) < DBL_EPSILON);
00187 
00188     ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00189     assert(ierr == CV_SUCCESS); ierr=ierr; // avoid unused var warning
00190 
00191     RecordStoppingPoint(tEnd);
00192 
00193     VerifyStateVariables();
00194 }
00195 
00196 
00197 void AbstractCvodeSystem::SetMaxSteps(long int numSteps)
00198 {
00199     mMaxSteps = numSteps;
00200 }
00201 
00202 long int AbstractCvodeSystem::GetMaxSteps()
00203 {
00204     return mMaxSteps;
00205 }
00206 
00207 void AbstractCvodeSystem::SetTolerances(double relTol, double absTol)
00208 {
00209     mRelTol = relTol;
00210     mAbsTol = absTol;
00211     ResetSolver();
00212 }
00213 
00214 double AbstractCvodeSystem::GetRelativeTolerance()
00215 {
00216     return mRelTol;
00217 }
00218 
00219 double AbstractCvodeSystem::GetAbsoluteTolerance()
00220 {
00221     return mAbsTol;
00222 }
00223 
00224 double AbstractCvodeSystem::GetLastStepSize()
00225 {
00226     return mLastInternalStepSize;
00227 }
00228 
00229 
00230 bool Differs(double v1, double v2)
00231 {
00232     return fabs(v1 - v2) > DBL_EPSILON;
00233 }
00234 
00235 
00236 void AbstractCvodeSystem::SetAutoReset(bool autoReset)
00237 {
00238     mAutoReset = autoReset;
00239     ResetSolver();
00240 }
00241 
00242 
00243 void AbstractCvodeSystem::ResetSolver()
00244 {
00245     DeleteVector(mLastSolutionState);
00246 }
00247 
00248 
00249 void AbstractCvodeSystem::SetupCvode(N_Vector initialConditions,
00250                                      realtype tStart,
00251                                      realtype maxDt)
00252 {
00253     assert((unsigned)NV_LENGTH_S(initialConditions) == GetNumberOfStateVariables());
00254     assert(maxDt >= 0.0);
00255 
00256     // Find out if we need to (re-)initialise
00257     bool reinit = !mpCvodeMem || mAutoReset || !mLastSolutionState || Differs(tStart, mLastSolutionTime);
00258     if (!reinit)
00259     {
00260         const unsigned size = GetNumberOfStateVariables();
00261         for (unsigned i=0; i<size; i++)
00262         {
00263             if (Differs(GetVectorComponent(mLastSolutionState, i), GetVectorComponent(mStateVariables, i)))
00264             {
00265                 reinit = true;
00266                 break;
00267             }
00268         }
00269     }
00270 
00271     if (!mpCvodeMem)
00272     {
00273         mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
00274         if (mpCvodeMem == NULL) EXCEPTION("Failed to SetupCvode CVODE");
00275         // Set error handler
00276         CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, NULL);
00277         // Set the user data
00278 #if CHASTE_SUNDIALS_VERSION >= 20400
00279         CVodeSetUserData(mpCvodeMem, (void*)(this));
00280 #else
00281         CVodeSetFdata(mpCvodeMem, (void*)(this));
00282 #endif
00283         // Setup CVODE
00284 #if CHASTE_SUNDIALS_VERSION >= 20400
00285         CVodeInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions);
00286         CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
00287 #else
00288         CVodeMalloc(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
00289                     CV_SS, mRelTol, &mAbsTol);
00290 #endif
00291         // Attach a linear solver for Newton iteration
00292         CVDense(mpCvodeMem, NV_LENGTH_S(initialConditions));
00293     }
00294     else if (reinit)
00295     {
00296 #if CHASTE_SUNDIALS_VERSION >= 20400
00297         CVodeReInit(mpCvodeMem, tStart, initialConditions);
00298         CVodeSStolerances(mpCvodeMem, mRelTol, mAbsTol);
00299 #else
00300         CVodeReInit(mpCvodeMem, AbstractCvodeSystemRhsAdaptor, tStart, initialConditions,
00301                     CV_SS, mRelTol, &mAbsTol);
00302 #endif
00303     }
00304     // Set max dt and change max steps if wanted
00305     CVodeSetMaxStep(mpCvodeMem, maxDt);
00306     if (mMaxSteps > 0)
00307     {
00308         CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
00309         CVodeSetMaxErrTestFails(mpCvodeMem, 15);
00310     }
00311 }
00312 
00313 
00314 void AbstractCvodeSystem::RecordStoppingPoint(double stopTime)
00315 {
00316     if (!mAutoReset)
00317     {
00318         const unsigned size = GetNumberOfStateVariables();
00319         CreateVectorIfEmpty(mLastSolutionState, size);
00320         for (unsigned i=0; i<size; i++)
00321         {
00322             SetVectorComponent(mLastSolutionState, i, GetVectorComponent(mStateVariables, i));
00323         }
00324         mLastSolutionTime = stopTime;
00325     }
00326 }
00327 
00328 
00329 void AbstractCvodeSystem::FreeCvodeMemory()
00330 {
00331     if (mpCvodeMem)
00332     {
00333         CVodeFree(&mpCvodeMem);
00334     }
00335     mpCvodeMem = NULL;
00336 }
00337 
00338 
00339 void AbstractCvodeSystem::CvodeError(int flag, const char * msg)
00340 {
00341 
00342     std::stringstream err;
00343     char* p_flag_name = CVodeGetReturnFlagName(flag);
00344     err << msg << ": " << p_flag_name;
00345     free(p_flag_name);
00346     std::cerr << err.str() << std::endl << std::flush;
00347     EXCEPTION(err.str());
00348 }
00349 
00350 
00351 #endif // CHASTE_CVODE
Generated on Thu Dec 22 13:00:16 2011 for Chaste by  doxygen 1.6.3