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

Generated on Mon Apr 18 11:35:28 2011 for Chaste by  doxygen 1.5.5