RungeKuttaFehlbergIvpOdeSolver.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2009
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 
00033 #include "RungeKuttaFehlbergIvpOdeSolver.hpp"
00034 #include "AbstractIvpOdeSolver.hpp"
00035 #include "AbstractOdeSystem.hpp"
00036 #include "OdeSolution.hpp"
00037 
00038 #include <cmath>
00039 #include <cfloat>
00040 //#include <iostream>
00041 #include <vector>
00042 
00043 const double smidge=1e-10;
00044 
00045 /*
00046  * PROTECTED FUNCTIONS =========================================================
00047  */
00048 
00049 void RungeKuttaFehlbergIvpOdeSolver::InternalSolve(OdeSolution& rSolution,
00050                                                 AbstractOdeSystem* pOdeSystem,
00051                                                 std::vector<double>& rYValues,
00052                                                 std::vector<double>& rWorkingMemory,
00053                                                 double startTime,
00054                                                 double endTime,
00055                                                 double maxTimeStep,
00056                                                 double minTimeStep,
00057                                                 double tolerance,
00058                                                 bool outputSolution)
00059 {
00060     const unsigned number_of_variables = pOdeSystem->GetNumberOfStateVariables();
00061     mError.reserve(number_of_variables);
00062     mk1.reserve(number_of_variables);
00063     mk2.reserve(number_of_variables);
00064     mk3.reserve(number_of_variables);
00065     mk4.reserve(number_of_variables);
00066     mk5.reserve(number_of_variables);
00067     mk6.reserve(number_of_variables);
00068     myk2.reserve(number_of_variables);
00069     myk3.reserve(number_of_variables);
00070     myk4.reserve(number_of_variables);
00071     myk5.reserve(number_of_variables);
00072     myk6.reserve(number_of_variables);
00073 
00074     double current_time = startTime;
00075     double time_step = maxTimeStep;
00076     bool got_to_end = false;
00077     bool accurate_enough = false;
00078     unsigned number_of_time_steps = 0;
00079 
00080     if (outputSolution)
00081     {   // Write out ICs
00082         rSolution.rGetTimes().push_back(current_time);
00083         rSolution.rGetSolutions().push_back(rYValues);
00084     }
00085 
00086     // should never get here if this bool has been set to true;
00087     assert(!mStoppingEventOccurred);
00088     while ( !got_to_end )
00089     {
00090         //std::cout << "New timestep\n" << std::flush;
00091         while (!accurate_enough)
00092         {
00093             accurate_enough = true; // assume it is OK until we check and find otherwise
00094 
00095             // Function that calls the appropriate one-step solver
00096             CalculateNextYValue(pOdeSystem,
00097                                 time_step,
00098                                 current_time,
00099                                 rYValues,
00100                                 rWorkingMemory);
00101 
00102             // Find the maximum error in this vector.
00103             double max_error = -DBL_MAX;
00104             for (unsigned i=0; i<number_of_variables; i++)
00105             {
00106                 if (mError[i] > max_error)
00107                 {
00108                     max_error = mError[i];
00109                 }
00110             }
00111 
00112             if (max_error > tolerance)
00113             {// Reject the step-size and do it again.
00114                 accurate_enough = false;
00115                 //std::cout << "Approximation rejected\n" << std::flush;
00116             }
00117             else
00118             {
00119                 // step forward the time since step has now been made
00120                 current_time = current_time + time_step;
00121                 //std::cout << "Approximation accepted with time step = "<< time_step << "\n" << std::flush;
00122                 //std::cout << "max_error = " << max_error << " tolerance = " << tolerance << "\n" << std::flush;
00123                 if (outputSolution)
00124                 {   // Write out ICs
00125                     //std::cout << "In solver Time = " << current_time << " y = " << rWorkingMemory[0] << "\n" << std::flush;
00126                     rSolution.rGetTimes().push_back(current_time);
00127                     rSolution.rGetSolutions().push_back(rWorkingMemory);
00128                     number_of_time_steps++;
00129                 }
00130             }
00131 
00132             // Set a new step size based on the accuracy here
00133             AdjustStepSize(time_step, max_error, tolerance, maxTimeStep, minTimeStep);
00134         }
00135 
00136         // For the next timestep check the step doesn't go past the end...
00137         if (current_time + time_step > endTime)
00138         {   // Allow a smaller timestep for the final step.
00139             time_step = endTime - current_time;
00140         }
00141 
00142         if ( pOdeSystem->CalculateStoppingEvent(current_time,
00143                                                 rWorkingMemory) == true )
00144         {
00145             mStoppingTime = current_time;
00146             mStoppingEventOccurred = true;
00147         }
00148 
00149         if (mStoppingEventOccurred || current_time>=endTime)
00150         {
00151             got_to_end = true;
00152         }
00153 
00154         // Approximation accepted - put it in rYValues
00155         rYValues.assign(rWorkingMemory.begin(), rWorkingMemory.end());
00156         accurate_enough = false; // for next loop.
00157         //std::cout << "Finished Time Step\n" << std::flush;
00158     }
00159     rSolution.SetNumberOfTimeSteps(number_of_time_steps);
00160 }
00161 
00162 
00163 void RungeKuttaFehlbergIvpOdeSolver::CalculateNextYValue(AbstractOdeSystem* pAbstractOdeSystem,
00164                                                   double timeStep,
00165                                                   double time,
00166                                                   std::vector<double>& currentYValues,
00167                                                   std::vector<double>& nextYValues)
00168 {
00169     const unsigned num_equations = pAbstractOdeSystem->GetNumberOfStateVariables();
00170 
00171 
00172     std::vector<double>& dy = nextYValues; // re-use memory (not that it makes much difference here!)
00173 
00174     pAbstractOdeSystem->EvaluateYDerivatives(time, currentYValues, dy);
00175 
00176     for (unsigned i=0;i<num_equations; i++)
00177     {
00178         mk1[i] = timeStep*dy[i];
00179         myk2[i] = currentYValues[i] + 0.25*mk1[i];
00180     }
00181 
00182     pAbstractOdeSystem->EvaluateYDerivatives(time + 0.25*timeStep, myk2, dy);
00183     for (unsigned i=0;i<num_equations; i++)
00184     {
00185         mk2[i] = timeStep*dy[i];
00186         myk3[i] = currentYValues[i] + 0.09375*mk1[i] + 0.28125*mk2[i];
00187     }
00188 
00189     pAbstractOdeSystem->EvaluateYDerivatives(time + 0.375*timeStep, myk3, dy);
00190     for (unsigned i=0;i<num_equations; i++)
00191     {
00192         mk3[i] = timeStep*dy[i];
00193         myk4[i] = currentYValues[i] + m1932o2197*mk1[i] - m7200o2197*mk2[i]
00194                     + m7296o2197*mk3[i];
00195     }
00196 
00197     pAbstractOdeSystem->EvaluateYDerivatives(time+m12o13*timeStep, myk4, dy);
00198     for (unsigned i=0;i<num_equations; i++)
00199     {
00200         mk4[i] = timeStep*dy[i];
00201         myk5[i] = currentYValues[i] + m439o216*mk1[i] - 8*mk2[i]
00202                     + m3680o513*mk3[i]- m845o4104*mk4[i];
00203     }
00204 
00205     pAbstractOdeSystem->EvaluateYDerivatives(time+timeStep, myk5, dy);
00206     for (unsigned i=0;i<num_equations; i++)
00207     {
00208         mk5[i] = timeStep*dy[i];
00209         myk6[i] = currentYValues[i] - m8o27*mk1[i] + 2*mk2[i] - m3544o2565*mk3[i]
00210                         + m1859o4104*mk4[i] - 0.275*mk5[i];
00211     }
00212 
00213     pAbstractOdeSystem->EvaluateYDerivatives(time+0.5*timeStep, myk6, dy);
00214     for (unsigned i=0;i<num_equations; i++)
00215     {
00216         mk6[i]=timeStep*dy[i];
00217         mError[i] = (1/timeStep)*fabs(m1o360*mk1[i] - m128o4275*mk3[i]
00218                     - m2197o75240*mk4[i] + 0.02*mk5[i]+ m2o55*mk6[i]);
00219         nextYValues[i] = currentYValues[i] + m25o216*mk1[i] + m1408o2565*mk3[i]
00220                         + m2197o4104*mk4[i] - 0.2*mk5[i];
00221     }
00222 }
00223 
00224 void RungeKuttaFehlbergIvpOdeSolver::AdjustStepSize(double& rCurrentStepSize,
00225                                 const double& rError,
00226                                 const double& rTolerance,
00227                                 const double& rMaxTimeStep,
00228                                 const double& rMinTimeStep)
00229 {
00230     // Work out scaling factor delta for the step size
00231     double delta = pow(rTolerance/(2.0*rError), 0.25);
00232 
00233     // Maximum adjustment is *0.1 or *4
00234     if (delta <= 0.1)
00235     {
00236         rCurrentStepSize *= 0.1;
00237     }
00238     else if (delta >= 4.0)
00239     {
00240         rCurrentStepSize *= 4.0;
00241     }
00242     else
00243     {
00244         rCurrentStepSize *= delta;
00245     }
00246 
00247     if (rCurrentStepSize > rMaxTimeStep)
00248     {
00249         rCurrentStepSize = rMaxTimeStep;
00250     }
00251 
00252     if (rCurrentStepSize < rMinTimeStep)
00253     {
00254         std::cout << "rCurrentStepSize = " << rCurrentStepSize << "\n" << std::flush;
00255         std::cout << "rMinTimeStep = " << rMinTimeStep << "\n" << std::flush;
00256 
00257         EXCEPTION("RKF45 Solver: Ode needs a smaller timestep than the set minimum\n");
00258     }
00259 
00260 }
00261 
00262 
00263 /*
00264  * PUBLIC FUNCTIONS=============================================================
00265  */
00266 
00267 OdeSolution RungeKuttaFehlbergIvpOdeSolver::Solve(AbstractOdeSystem* pOdeSystem,
00268                                                std::vector<double>& rYValues,
00269                                                double startTime,
00270                                                double endTime,
00271                                                double timeStep,
00272                                                double tolerance)
00273 {
00274     assert(rYValues.size()==pOdeSystem->GetNumberOfStateVariables());
00275     assert(endTime > startTime);
00276     assert(timeStep > 0.0);
00277 
00278     mStoppingEventOccurred = false;
00279     if ( pOdeSystem->CalculateStoppingEvent(startTime, rYValues) == true )
00280     {
00281         EXCEPTION("(Solve with sampling) Stopping event is true for initial condition");
00282     }
00283     // Allocate working memory
00284     std::vector<double> working_memory(rYValues.size());
00285     // And solve...
00286     OdeSolution solutions;
00287     //solutions.SetNumberOfTimeSteps((unsigned)(10.0*(startTime-endTime)/timeStep));
00288     bool return_solution = true;
00289     InternalSolve(solutions, pOdeSystem, rYValues, working_memory, startTime, endTime, timeStep, 1e-5, tolerance, return_solution);
00290     return solutions;
00291 }
00292 
00293 void RungeKuttaFehlbergIvpOdeSolver::Solve(AbstractOdeSystem* pOdeSystem,
00294                                         std::vector<double>& rYValues,
00295                                         double startTime,
00296                                         double endTime,
00297                                         double timeStep)
00298 {
00299     assert(rYValues.size()==pOdeSystem->GetNumberOfStateVariables());
00300     assert(endTime > startTime);
00301     assert(timeStep > 0.0);
00302 
00303     mStoppingEventOccurred = false;
00304     if ( pOdeSystem->CalculateStoppingEvent(startTime, rYValues) == true )
00305     {
00306         EXCEPTION("(Solve without sampling) Stopping event is true for initial condition");
00307     }
00308     // Allocate working memory
00309     std::vector<double> working_memory(rYValues.size());
00310     // And solve...
00311     OdeSolution not_required_solution;
00312     bool return_solution = false;
00313     InternalSolve(not_required_solution, pOdeSystem, rYValues, working_memory, startTime, endTime, timeStep, 1e-4, 1e-5, return_solution);
00314 }
00315 
00316 

Generated on Wed Mar 18 12:51:56 2009 for Chaste by  doxygen 1.5.5