Ralston's Second Order Method

Description

Ralston's second order method is a Runge-Kutta method for approximating the solution of the initial value problem y'(x) = f(x,y);  y(x0) = y0 which evaluates the integrand,f(x,y), twice for step. For step i+1,
yi+1 = yi + 1/4 ( k1 + 3 k2 ),
where
k1 = h f(xi, yi),
k2 = h f(xi + 2 h / 3, yi + 2 k1 / 3 ),
and xi = x0 + i h.
Ralston's second order method is a second order procedure for which Richardson extrapolation can be used.

Function List

C Source

////////////////////////////////////////////////////////////////////////////////
// File: runge_kutta_ralston_2.c                                              //
// Routines:                                                                  //
//    Runge_Kutta_Ralston_2_Method                                            //
//    Runge_Kutta_Ralston_2_Richardson                                        //
//    Runge_Kutta_Ralston_2_Integral_Curve                                    //
//    Runge_Kutta_Ralston_2_Richardson_Integral_Curve                         //
////////////////////////////////////////////////////////////////////////////////
 
////////////////////////////////////////////////////////////////////////////////
//                                                                            //
//  Description:                                                              //
//     The second order Ralston's method is a second order Runge-Kutta method //
//     for approximating the solution of the differential equation            //
//     y'(x) = f(x,y) with initial condition y(x0) = c.  The second order     //
//     Ralston's method evaluates f(x,y) twice per step. For step i+1,        //
//     y[i+1] = y[i] + 0.25 * ( k1 + 3k2 ),  where k1 = h * f(x[i],y[i]),     //
//     k2 = h * f(x[i]+2/3 h, y[i]+2/3 k1), and x[i] = x0 + i * h.            //
//                                                                            //
//     Ralston's method is a second order procedure for which Richardson      //
//     extrapolation can be used.                                             //
//                                                                            //
////////////////////////////////////////////////////////////////////////////////

static const double two_thirds = 2.0 / 3.0;

////////////////////////////////////////////////////////////////////////////////
//  double Runge_Kutta_Ralston_2_Method( double (*f)(double, double),         //
//                    double y0, double x0, double h, int number_of_steps );  //
//                                                                            //
//  Description:                                                              //
//     This routine uses the second order Ralston method to approximate the   //
//     solution at x = x0 + h * number_of_steps of the initial value problem  //
//     y'=f(x,y), y(x0) = y0.                                                 //
//                                                                            //
//  Arguments:                                                                //
//     double *f                                                              //
//            Pointer to the function which returns the slope at (x,y) of the //
//            integral curve of the differential equation y' = f(x,y) which   //
//            passes through the point (x0,y0).                               //
//     double y0                                                              //
//            The initial value of y at x = x0.                               //
//     double x0                                                              //
//            The initial value of x.                                         //
//     double h                                                               //
//            The step size.                                                  //
//     int    number_of_steps                                                 //
//            The number of steps. Must be a nonnegative integer.             //
//                                                                            //
//  Return Values:                                                            //
//     The solution of the initial value problem y' = f(x,y), y(x0) = y0 at   //
//     x = x0 + number_of_steps * h.                                          //
//                                                                            //
////////////////////////////////////////////////////////////////////////////////
//                                                                            //
double Runge_Kutta_Ralston_2_Method( double (*f)(double, double), double y0,
                                   double x0, double h, int number_of_steps ) {

   double k1;
   double h3 = two_thirds * h;
   double h4 = h / 4.0;

   while ( --number_of_steps >= 0 ) {
      k1 = (*f)(x0,y0);
      y0 += h4 * ( k1 + 3.0 * (*f)(x0 + h3, y0 + h3 * k1) ); 
      x0 += h;
   }

   return y0;
}


////////////////////////////////////////////////////////////////////////////////
//  double Runge_Kutta_Ralston_2_Richardson( double (*f)(double, double),     //
//                       double y0, double x0, double h, int number_of_steps, //
//                                                   int richardson_columns)  //
//                                                                            //
//  Description:                                                              //
//     This routine uses the second order Ralston method together with        //
//     Richardson extrapolation to approximate the solution at                //
//     x = x0 + h * number_of_steps of the initial value problem              //
//     y'=f(x,y), y(x0) = y0.                                                 //
//                                                                            //
//  Arguments:                                                                //
//     double *f                                                              //
//            Pointer to the function which returns the slope at (x,y) of the //
//            integral curve of the differential equation y' = f(x,y) which   //
//            passes through the point (x0,y0).                               //
//     double y0                                                              //
//            The initial value of y at x = x0.                               //
//     double x0                                                              //
//            The initial value of x.                                         //
//     double h                                                               //
//            The step size.                                                  //
//     int    number_of_steps                                                 //
//            The number of steps. Must be nonnegative.                       //
//     int    richardson_columns                                              //
//            The maximum number of columns to use in the Richardson          //
//            extrapolation to the limit.                                     //
//                                                                            //
//  Return Values:                                                            //
//     The solution of the initial value problem y' = f(x,y), y(x0) = y0 at   //
//     x = x0 + number_of_steps * h.                                          //
//                                                                            //
////////////////////////////////////////////////////////////////////////////////
//                                                                            //
static const double richardson[] = {  
  1.0 / 3.0, 1.0 / 7.0, 1.0 / 15.0, 1.0 / 31.0, 1.0 / 63.0
};

#define MAX_COLUMNS 1+sizeof(richardson)/sizeof(richardson[0])

#define max(x,y) ( (x) < (y) ? (y) : (x) )
#define min(x,y) ( (x) < (y) ? (x) : (y) )

//                                                                            //
////////////////////////////////////////////////////////////////////////////////
//                                                                            //
double Runge_Kutta_Ralston_2_Richardson( double (*f)(double, double), double y0,
           double x0, double h, int number_of_steps, int richardson_columns ) {

   double dt[MAX_COLUMNS];         // dt[i] is the last element in column i.
   double integral, delta, h_used;
   int j,k, number_sub_intervals;

   richardson_columns = max(1, min(MAX_COLUMNS, richardson_columns));
   while ( --number_of_steps >= 0 ) {
      h_used = h;
      number_sub_intervals = 1;
      for (j = 0; j < richardson_columns; j++) {
         integral = Runge_Kutta_Ralston_2_Method( f, y0, x0, h_used, 
                                                        number_sub_intervals);
         for ( k = 0; k < j; k++) {
            delta = integral - dt[k];
            dt[k] = integral;
            integral += richardson[k] * delta;
         }
         dt[j] = integral;
         h_used *= 0.5;
         number_sub_intervals += number_sub_intervals;
      }
      y0 = integral;
      x0 += h;
   }
   return y0;
}


////////////////////////////////////////////////////////////////////////////////
//  void Runge_Kutta_Ralston_2_Integral_Curve( double (*f)(double, double),   //
//        double y[], double x0, double h, int number_of_steps_per_interval,  //
//                                               int number_of_intervals );   //
//                                                                            //
//  Description:                                                              //
//     This routine uses the second order Ralston method to approximate the   //
//     solution of the differential equation y'=f(x,y) with the initial       //
//     condition y = y[0] at x = x0.  The values are returned in y[n] which   //
//     is the value of y evaluated at x = x0 + n * m * h, where m is the      //
//     number of steps per interval and n is the interval number,             //
//     0 <= n <= number_of_intervals.                                         //
//                                                                            //
//  Arguments:                                                                //
//     double *f                                                              //
//            Pointer to the function which returns the slope at (x,y) of the //
//            integral curve of the differential equation y' = f(x,y) which   //
//            which passes through the point (x0,y[0]).                       //
//     double y[]                                                             //
//            On input y[0] is the initial value of y at x = x0. On output    //
//            for n >= 1,  y[n] is the appproximation of the solution y(x) of //
//            the initial value problem y' = f(x,y), y(x0) = y[0],  at        //
//            x = x0 + nmh, where m is the number of steps per interval and n //
//            is the interval number.                                         //
//     double x0                                                              //
//            Initial value of x.                                             //
//     double h                                                               //
//            Step size                                                       //
//     int    number_of_steps_per_interval                                    //
//            The number of steps of length h used to calculate y[i+1]        //
//            starting from y[i].                                             //
//     int    number_of_intervals                                             //
//            The number of intervals,  y[] should be dimensioned at least    //
//            number_of_intervals + 1.                                        //
//                                                                            //
//  Return Values:                                                            //
//     This routine is of type void and hence does not return a value.        //
//     The solution of y' = f(x,y) from x = x0 to x = x0 + n * m * h,         //
//     where n is the number of intervals and m is the number of steps per    //
//     interval, is stored in the input array y[].                            //
//                                                                            //
////////////////////////////////////////////////////////////////////////////////
//                                                                            //
void Runge_Kutta_Ralston_2_Integral_Curve( double (*f)(double, double),
            double y[], double x0, double h, int number_of_steps_per_interval,
                                                    int number_of_intervals ) {

   double k1;
   double h3 = two_thirds * h;
   double h4 = h / 4.0;
   int i;

   while ( --number_of_intervals >= 0 ) {
      *(y+1) = *y;
      y++;
      for (i = 0; i < number_of_steps_per_interval; i++) {
         k1 = (*f)(x0,*y);
         *y += h4 * ( k1 + 3.0 * (*f)(x0+h3, *y + h3 * k1) ); 
         x0 += h;
      }
   }
}


////////////////////////////////////////////////////////////////////////////////
//  void Runge_Kutta_Ralston_2_Richardson_Integral_Curve(                     //
//              double (*f)(double, double), double y[], double x0, double h, //
//                 int number_of_steps_per_interval, int number_of_intervals, //
//                                                  int richardson_columns )  //
//                                                                            //
//  Description:                                                              //
//     This routine uses the second order Ralston method together with        //
//     Richardson extrapolation to the limit (as h -> 0) to approximate the   //
//     solution of the differential equation y'=f(x,y) with the initial       //
//     condition y = y[0] at x = x0.                                          //
//     The values are returned in y[], in which y[n] is the value of y        //
//     evaluated at x = x0 + n * m * h, where m is the number of steps per    //
//     interval and n is the interval number, 0 <= n <= number_of_intervals.  //
//                                                                            //
//  Arguments:                                                                //
//     double *f                                                              //
//            Pointer to the function which returns the slope at (x,y) of the //
//            integral curve of the differential equation y' = f(x,y) which   //
//            which passes through the point (x0,y[0]).                       //
//     double y[]                                                             //
//            On input y[0] is the initial value of y at x = x0. On output    //
//            for n >= 1,  y[n] is the appproximation of the solution y(x) of //
//            the initial value problem y' = f(x,y), y(x0) = y[0],  at        //
//            x = x0 + nmh, where m is the number of steps per interval and n //
//            is the interval number.                                         //
//     double x0                                                              //
//            Initial value of x.                                             //
//     double h                                                               //
//            Step size                                                       //
//     int    number_of_steps_per_interval                                    //
//            The number of steps of length h used to calculate y[i+1]        //
//            starting from y[i].                                             //
//     int    number_of_intervals                                             //
//            The number of intervals,  y[] should be dimensioned at least    //
//            number_of_intervals + 1.                                        //
//     int    richardson_columns                                              //
//            The maximum number of columns to use in the Richardson          //
//            extrapolation to the limit.                                     //
//                                                                            //
//  Return Values:                                                            //
//     This routine is of type void and hence does not return a value.        //
//     The solution of y' = f(x,y) from x = x0 to x = x0 + n * m * h,         //
//     where n is the number of intervals and m is the number of steps per    //
//     interval, is stored in the input array y[].                            //
//                                                                            //
////////////////////////////////////////////////////////////////////////////////
//                                                                            //
void Runge_Kutta_Ralston_2_Richardson_Integral_Curve(
                  double (*f)(double, double), double y[], double x0, double h,
                  int number_of_steps_per_interval, int number_of_intervals,
                                                     int richardson_columns ) {

   double mh = (double) number_of_steps_per_interval * h;

   while ( --number_of_intervals >= 0 ) {
      *(++y) = Runge_Kutta_Ralston_2_Richardson( f, *y, x0, h,
                            number_of_steps_per_interval, richardson_columns );
      x0 += mh;
   }
   return;
}