Page MenuHomec4science

runge_kutta_step.inl
No OneTemporary

File Metadata

Created
Sat, Nov 2, 06:42

runge_kutta_step.inl

/*-------------------------------------------------------------------------------
Copyright (c) 2014,2015 F. Georget <fabieng@princeton.edu>, Princeton University
All rights reserved.
Redistribution and use in source and binary forms, with or without modification,
are permitted provided that the following conditions are met:
1. Redistributions of source code must retain the above copyright notice, this
list of conditions and the following disclaimer.
2. Redistributions in binary form must reproduce the above copyright notice,
this list of conditions and the following disclaimer in the documentation and/or
other materials provided with the distribution.
3. Neither the name of the copyright holder nor the names of its contributors
may be used to endorse or promote products derived from this software without
specific prior written permission.
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR
ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
(INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
(INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
-----------------------------------------------------------------------------*/
#include "runge_kutta_step.hpp" // syntaxic coloration only
namespace specmicp {
namespace odeint {
template <int dim, class butcher_tableau_t>
void EmbeddedRungeKuttaStep45<dim, butcher_tableau_t>::run_step(const vec_t& y,
const vec_t& dydx,
scalar_t x,
scalar_t timestep,
vec_t& y_out,
vec_t& y_err)
{
vec_t ak2, ak3, ak4, ak5, ak6;
vec_t ytemp;
ytemp = y + butcher().b(2, 1)*timestep*dydx;
get_rhs(x + butcher().a(2)*timestep, ytemp, ak2);
ytemp = y + timestep*(butcher().b(3,1)*dydx + butcher().b(3, 2)*ak2);
get_rhs(x + butcher().a(3)*timestep, ytemp, ak3);
ytemp = y + timestep*(butcher().b(4,1)*dydx + butcher().b(4, 2)*ak2 + butcher().b(4, 3)*ak3);
get_rhs(x + butcher().a(4)*timestep, ytemp, ak4);
ytemp = y + timestep*(butcher().b(5,1)*dydx + butcher().b(5, 2)*ak2
+ butcher().b(5, 3)*ak3 + butcher().b(5, 4)*ak4);
get_rhs(x + butcher().a(5)*timestep, ytemp, ak5);
ytemp = y + timestep*(butcher().b(6,1)*dydx + butcher().b(6, 2)*ak2
+ butcher().b(6, 3)*ak3 + butcher().b(6, 4)*ak4
+ butcher().b(6, 5)*ak5);
get_rhs(x + butcher().a(6)*timestep, ytemp, ak6);
if (butcher_tableau_t::RK_evaluations == 6)
{
y_out = y + timestep*(butcher().c(1)*dydx + butcher().c(2)*ak2 + butcher().c(3)*ak3
+ butcher().c(4)*ak4 + butcher().c(5)*ak5 + butcher().c(6)*ak6);
y_err = timestep * ( (butcher().c(1) - butcher().cs(1))*dydx + (butcher().c(2) - butcher().cs(2))*ak2 +
(butcher().c(3) - butcher().cs(3))*ak3 + (butcher().c(4) - butcher().cs(4))*ak4 +
(butcher().c(5) - butcher().cs(5))*ak5 + (butcher().c(6) - butcher().cs(6))*ak6
);
}
else
{
vec_t ak7;
ytemp = y + timestep*(butcher().b(7,1)*dydx + butcher().b(7, 2)*ak2
+ butcher().b(7, 3)*ak3 + butcher().b(7, 4)*ak4
+ butcher().b(7, 5)*ak5 + butcher().b(7, 6)*ak6);
get_rhs(x + butcher().a(7)*timestep, ytemp, ak7);
y_out = y + timestep*(butcher().c(1)*dydx + butcher().c(2)*ak2 + butcher().c(3)*ak3
+ butcher().c(4)*ak4 + butcher().c(5)*ak5 + butcher().c(6)*ak6
+ butcher().c(7)*ak7);
y_err = timestep * ( (butcher().c(1) - butcher().cs(1))*dydx + (butcher().c(2) - butcher().cs(2))*ak2 +
(butcher().c(3) - butcher().cs(3))*ak3 + (butcher().c(4) - butcher().cs(4))*ak4 +
(butcher().c(5) - butcher().cs(5))*ak5 + (butcher().c(6) - butcher().cs(6))*ak6 +
(butcher().c(7) - butcher().cs(7))*ak7
);
}
if (get_options().non_negativity)
set_non_negative<dim>(y_out);
}
template <int dim, class butcher_tableau_t>
void EmbeddedRungeKuttaStep45<dim, butcher_tableau_t>::rk_step(vector_t<dim>& y,
const vector_t<dim>& dydx,
double& x,
StepLength& stepl)
{
vector_t<dim> y_err;
vector_t<dim> y_temp;
double h = stepl.test;
double xnew;
scalar_t errmax;
for (int i=0; i<get_options().max_try; ++i)
{
run_step(y, dydx, x, h, y_temp, y_err);
errmax = max_coeff<dim>(y_err) / get_options().eps;
if (errmax < 1) break;
double htemp = get_options().safety*h*std::pow(errmax, get_options().p_shrink);
h = (h >= 0.0 ? std::max(htemp, 0.1*h) : std::min(htemp, 0.1*h));
xnew = x+h;
if (xnew == x) {throw std::runtime_error("stepsize underflow");}
}
stepl.did = h;
x += h;
if (errmax > get_options().errcon) {stepl.next = get_options().safety*h*std::pow(errmax, get_options().p_grow);}
else {stepl.next = 5.0*h;}
y = y_temp;
}
} // end namespace odeint
} // end namespace specmicp

Event Timeline