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

Generated by  doxygen 1.6.2