AbstractCvodeCell.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 <sstream>
00032 #include <iostream>
00033 #include <cmath>
00034 
00035 #include "AbstractCvodeCell.hpp"
00036 #include "CvodeAdaptor.hpp" // For CvodeErrorHandler
00037 #include "TimeStepper.hpp"
00038 #include "Exception.hpp"
00039 #include "HeartConfig.hpp"
00040 
00041 // CVODE headers
00042 #include <cvode/cvode.h>
00043 #include <sundials/sundials_nvector.h>
00044 #include <cvode/cvode_dense.h>
00045 
00055 int AbstractCvodeCellRhsAdaptor(realtype t, N_Vector y, N_Vector ydot, void *pData)
00056 {
00057     assert(pData != NULL);
00058     AbstractCvodeCell* pCell = (AbstractCvodeCell*) pData;
00059     try
00060     {
00061         pCell->EvaluateRhs(t, y, ydot);
00062     }
00063     catch (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 
00073 AbstractCvodeCell::AbstractCvodeCell(boost::shared_ptr<AbstractIvpOdeSolver> pSolver /* unused */,
00074                                      unsigned numberOfStateVariables,
00075                                      unsigned voltageIndex,
00076                                      boost::shared_ptr<AbstractStimulusFunction> pIntracellularStimulus)
00077     : AbstractCardiacCellInterface(pSolver, voltageIndex, pIntracellularStimulus),
00078       AbstractParameterisedSystem<N_Vector>(numberOfStateVariables),
00079       mpCvodeMem(NULL),
00080       mMaxSteps(0)
00081 {
00082     SetTolerances();
00083 }
00084 
00085 
00086 AbstractCvodeCell::~AbstractCvodeCell()
00087 {
00088     if (mStateVariables)
00089     {
00090         mStateVariables->ops->nvdestroy(mStateVariables);
00091         mStateVariables = NULL;
00092     }
00093     if (mParameters)
00094     {
00095         mParameters->ops->nvdestroy(mParameters);
00096         mParameters = NULL;
00097     }
00098 }
00099 
00100 
00101 void AbstractCvodeCell::Init()
00102 {
00103     SetStateVariables(GetInitialConditions());
00104     if (mParameters)
00105     {
00106         mParameters->ops->nvdestroy(mParameters);
00107         mParameters = NULL;
00108     }
00109     mParameters = N_VNew_Serial(rGetParameterNames().size());
00110     for (unsigned i=0; i<NV_LENGTH_S(mParameters); i++)
00111     {
00112         NV_Ith_S(mParameters, i) = 0.0;
00113     }
00114 }
00115 
00116 
00117 double AbstractCvodeCell::GetVoltage()
00118 {
00119     assert(mStateVariables);
00120     return NV_Ith_S(mStateVariables, mVoltageIndex);
00121 }
00122 
00123 void AbstractCvodeCell::SetVoltage(double voltage)
00124 {
00125     assert(mStateVariables);
00126     NV_Ith_S(mStateVariables, mVoltageIndex) = voltage;
00127 }
00128 
00129 
00130 void AbstractCvodeCell::ResetToInitialConditions()
00131 {
00132     SetStateVariables(GetInitialConditions());
00133 }
00134 
00135 N_Vector AbstractCvodeCell::GetInitialConditions()
00136 {
00137     assert(mpSystemInfo);
00138     std::vector<double> inits = mpSystemInfo->GetInitialConditions();
00139     N_Vector v = N_VNew_Serial(inits.size());
00140     for (unsigned i=0; i<inits.size(); i++)
00141     {
00142         NV_Ith_S(v, i) = inits[i];
00143     }
00144     return v;
00145 }
00146 
00147 N_Vector AbstractCvodeCell::GetStateVariables()
00148 {
00149     assert(mpSystemInfo);
00150     assert(mStateVariables);
00151     return CopyVector(mStateVariables);
00152 }
00153 
00154 void AbstractCvodeCell::SetStateVariables(N_Vector stateVars)
00155 {
00156     if (mStateVariables and stateVars != mStateVariables)
00157     {
00158         mStateVariables->ops->nvdestroy(mStateVariables);
00160     }
00161     mStateVariables = stateVars;
00162 }
00163 
00164 void AbstractCvodeCell::SetStateVariablesUsingACopyOfThisVector(N_Vector stateVars)
00165 {
00166     // Cope if stateVars == mStateVariables
00167     N_Vector temp = CopyVector(stateVars);
00168     if (mStateVariables)
00169     {
00170         mStateVariables->ops->nvdestroy(mStateVariables);
00171     }
00172     mStateVariables = temp;
00173 }
00174 
00175 
00176 void AbstractCvodeCell::SetVoltageDerivativeToZero(bool clamp)
00177 {
00178     mSetVoltageDerivativeToZero = clamp;
00179 }
00180 
00181 
00182 OdeSolution AbstractCvodeCell::Compute(double tStart, double tEnd, double tSamp)
00183 {
00184     if (tSamp == 0.0)
00185     {
00186         tSamp = HeartConfig::Instance()->GetPrintingTimeStep();
00187     }
00188     return Solve(tStart, tEnd, tSamp, tSamp);
00189 }
00190 
00191 
00192 void AbstractCvodeCell::ComputeExceptVoltage(double tStart, double tEnd)
00193 {
00194     EXCEPTION("This method is not yet implemented for CVODE cells.");
00195 }
00196 
00197 
00198 OdeSolution AbstractCvodeCell::Solve(realtype tStart,
00199                                      realtype tEnd,
00200                                      realtype maxDt,
00201                                      realtype tSamp)
00202 {
00203     assert(tEnd >= tStart);
00204     assert(tSamp > 0.0);
00205 
00206     SetupCvode(mStateVariables, tStart, maxDt);
00207 
00208     TimeStepper stepper(tStart, tEnd, tSamp);
00209 
00210     // Set up ODE solution
00211     OdeSolution solutions;
00212     solutions.SetNumberOfTimeSteps(stepper.EstimateTimeSteps());
00213     solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
00214     solutions.rGetTimes().push_back(tStart);
00215     solutions.SetOdeSystemInformation(mpSystemInfo);
00216 
00217     // Main time sampling loop
00218     while (!stepper.IsTimeAtEnd())
00219     {
00220 //        std::cout << "Solving to time " << stepper.GetNextTime() << std::endl << std::flush;
00221         double cvode_stopped_at;
00222         int ierr;
00223         try
00224         {
00225             ierr = CVode(mpCvodeMem, stepper.GetNextTime(), mStateVariables,
00226                          &cvode_stopped_at, CV_NORMAL);
00227         }
00228         catch (...)
00229         {
00230             FreeCvodeMemory();
00231             throw;
00232         }
00233         if (ierr<0) CvodeError(ierr, "CVODE failed to solve system");
00234         // Not root finding, so should have reached requested time
00235         assert(fabs(cvode_stopped_at - stepper.GetNextTime()) < DBL_EPSILON);
00236 
00237         VerifyStateVariables();
00238 
00239         // Store solution
00240         solutions.rGetSolutions().push_back(MakeStdVec(mStateVariables));
00241         solutions.rGetTimes().push_back(cvode_stopped_at);
00242         stepper.AdvanceOneTimeStep();
00243     }
00244 
00245     // stepper.EstimateTimeSteps may have been an overestimate...
00246     solutions.SetNumberOfTimeSteps(stepper.GetTotalTimeStepsTaken());
00247 
00248     int ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00249     assert(ierr == CV_SUCCESS); ierr=ierr; // avoid unused var warning
00250     FreeCvodeMemory();
00251 
00252     return solutions;
00253 }
00254 
00255 void AbstractCvodeCell::Solve(realtype tStart,
00256                               realtype tEnd,
00257                               realtype maxDt)
00258 {
00259     assert(tEnd >= tStart);
00260 
00261     SetupCvode(mStateVariables, tStart, maxDt);
00262 
00263     double cvode_stopped_at;
00264     int ierr;
00265     try
00266     {
00267         ierr = CVode(mpCvodeMem, tEnd, mStateVariables, &cvode_stopped_at, CV_NORMAL);
00268     }
00269     catch (...)
00270     {
00271         FreeCvodeMemory();
00272         throw;
00273     }
00274     if (ierr<0) CvodeError(ierr, "CVODE failed to solve system");
00275     // Not root finding, so should have reached requested time
00276     assert(fabs(cvode_stopped_at - tEnd) < DBL_EPSILON);
00277 
00278     ierr = CVodeGetLastStep(mpCvodeMem, &mLastInternalStepSize);
00279     assert(ierr == CV_SUCCESS); ierr=ierr; // avoid unused var warning
00280     FreeCvodeMemory();
00281 
00282     VerifyStateVariables();
00283 }
00284 
00285 
00286 void AbstractCvodeCell::SetMaxSteps(long int numSteps)
00287 {
00288     mMaxSteps = numSteps;
00289 }
00290 
00291 long int AbstractCvodeCell::GetMaxSteps()
00292 {
00293     return mMaxSteps;
00294 }
00295 
00296 void AbstractCvodeCell::SetTolerances(double relTol, double absTol)
00297 {
00298     mRelTol = relTol;
00299     mAbsTol = absTol;
00300 }
00301 
00302 double AbstractCvodeCell::GetRelativeTolerance()
00303 {
00304     return mRelTol;
00305 }
00306 
00307 double AbstractCvodeCell::GetAbsoluteTolerance()
00308 {
00309     return mAbsTol;
00310 }
00311 
00312 double AbstractCvodeCell::GetLastStepSize()
00313 {
00314     return mLastInternalStepSize;
00315 }
00316 
00317 
00318 void AbstractCvodeCell::SetupCvode(N_Vector initialConditions,
00319                                    realtype tStart,
00320                                    realtype maxDt)
00321 {
00322     assert(NV_LENGTH_S(initialConditions) == GetNumberOfStateVariables());
00323     assert(maxDt >= 0.0);
00324 
00325     mpCvodeMem = CVodeCreate(CV_BDF, CV_NEWTON);
00326     if (mpCvodeMem == NULL) EXCEPTION("Failed to SetupCvode CVODE");
00327     // Set error handler
00328     CVodeSetErrHandlerFn(mpCvodeMem, CvodeErrorHandler, NULL);
00329     // Set the user data
00330     CVodeSetFdata(mpCvodeMem, (void*)(this));
00331     // Setup CVODE
00332     CVodeMalloc(mpCvodeMem, AbstractCvodeCellRhsAdaptor, tStart, initialConditions,
00333                 CV_SS, mRelTol, &mAbsTol);
00334     CVodeSetMaxStep(mpCvodeMem, maxDt);
00335     // Change max steps if wanted
00336     if (mMaxSteps > 0)
00337     {
00338         CVodeSetMaxNumSteps(mpCvodeMem, mMaxSteps);
00339     }
00340     // Attach a linear solver for Newton iteration
00341     CVDense(mpCvodeMem, NV_LENGTH_S(initialConditions));
00342 }
00343 
00344 
00345 void AbstractCvodeCell::FreeCvodeMemory()
00346 {
00347     CVodeFree(&mpCvodeMem);
00348 }
00349 
00350 
00351 // Errors get caught before they hit this, but we keep it just in case...
00352 #define COVERAGE_IGNORE
00353 void AbstractCvodeCell::CvodeError(int flag, const char * msg)
00354 {
00355     std::stringstream err;
00356     err << msg << ": " << CVodeGetReturnFlagName(flag);
00357     std::cerr << err.str() << std::endl << std::flush;
00358     EXCEPTION(err.str());
00359 }
00360 #undef COVERAGE_IGNORE
00361 
00362 
00363 std::vector<double> AbstractCvodeCell::MakeStdVec(N_Vector v)
00364 {
00365     unsigned size = NV_LENGTH_S(v);
00366     std::vector<double> sv(size);
00367     for (unsigned i=0; i<size; i++)
00368     {
00369         sv[i] = NV_Ith_S(v, i);
00370     }
00371     return sv;
00372 }
00373 
00374 
00375 std::string AbstractCvodeCell::DumpState(const std::string& message,
00376                                          N_Vector Y)
00377 {
00378     std::stringstream res;
00379     res << message << "\nState:\n";
00380     if (Y == NULL)
00381     {
00382         Y = mStateVariables;
00383     }
00384     for (unsigned i=0; i<NV_LENGTH_S(Y); i++)
00385     {
00386         res << "\t" << rGetStateVariableNames()[i] << ":" << NV_Ith_S(Y, i) << "\n";
00387     }
00388     return res.str();
00389 }
00390 
00391 N_Vector AbstractCvodeCell::CopyVector(N_Vector originalVec)
00392 {
00393     unsigned size = NV_LENGTH_S(originalVec);
00394     N_Vector v = N_VClone(originalVec);
00395     for (unsigned i=0; i<size; i++)
00396     {
00397         NV_Ith_S(v, i) = NV_Ith_S(originalVec, i);
00398     }
00399     return v;
00400 }
00401 
00402 
00403 #endif // CHASTE_CVODE

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