Blame interpolation/steffen.c

Packit 67cb25
/* interpolation/steffen.c
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 2014 Jean-François Caron
Packit 67cb25
 * 
Packit 67cb25
 * This program is free software; you can redistribute it and/or modify
Packit 67cb25
 * it under the terms of the GNU General Public License as published by
Packit 67cb25
 * the Free Software Foundation; either version 3 of the License, or (at
Packit 67cb25
 * your option) any later version.
Packit 67cb25
 * 
Packit 67cb25
 * This program is distributed in the hope that it will be useful, but
Packit 67cb25
 * WITHOUT ANY WARRANTY; without even the implied warranty of
Packit 67cb25
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
Packit 67cb25
 * General Public License for more details.
Packit 67cb25
 * 
Packit 67cb25
 * You should have received a copy of the GNU General Public License
Packit 67cb25
 * along with this program; if not, write to the Free Software
Packit 67cb25
 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
/* Author:  J.-F. Caron
Packit 67cb25
 *
Packit 67cb25
 * This interpolation method is taken from 
Packit 67cb25
 * M.Steffen, "A simple method for monotonic interpolation in one dimension",
Packit 67cb25
 * Astron. Astrophys. 239, 443-450 (1990).
Packit 67cb25
 *
Packit 67cb25
 * This interpolation method guarantees monotonic interpolation functions between
Packit 67cb25
 * the given data points.  A consequence of this is that extremal values can only
Packit 67cb25
 * occur at the data points.  The interpolating function and its first derivative
Packit 67cb25
 * are guaranteed to be continuous, but the second derivative is not.
Packit 67cb25
 *
Packit 67cb25
 * The implementation is modelled on the existing Akima interpolation method
Packit 67cb25
 * previously included in GSL by Gerard Jungman.
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
#include <config.h>
Packit 67cb25
#include <stdlib.h>
Packit 67cb25
#include <math.h>
Packit 67cb25
#include <gsl/gsl_math.h>
Packit 67cb25
#include <gsl/gsl_errno.h>
Packit 67cb25
#include "integ_eval.h"
Packit 67cb25
#include <gsl/gsl_interp.h>
Packit 67cb25
Packit 67cb25
typedef struct
Packit 67cb25
{
Packit 67cb25
  double * a;       /* eqs 2-5 of paper */
Packit 67cb25
  double * b;
Packit 67cb25
  double * c;
Packit 67cb25
  double * d;
Packit 67cb25
Packit 67cb25
  double * y_prime; /* eq 11 of paper */
Packit 67cb25
} steffen_state_t;
Packit 67cb25
Packit 67cb25
static void steffen_free (void * vstate);
Packit 67cb25
static double steffen_copysign(const double x, const double y);
Packit 67cb25
Packit 67cb25
static void *
Packit 67cb25
steffen_alloc (size_t size)
Packit 67cb25
{
Packit 67cb25
  steffen_state_t *state;
Packit 67cb25
  
Packit 67cb25
  state = (steffen_state_t *) calloc (1, sizeof (steffen_state_t));
Packit 67cb25
  
Packit 67cb25
  if (state == NULL)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR_NULL("failed to allocate space for state", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->a = (double *) malloc (size * sizeof (double));
Packit 67cb25
  
Packit 67cb25
  if (state->a == NULL)
Packit 67cb25
    {
Packit 67cb25
      steffen_free(state);
Packit 67cb25
      GSL_ERROR_NULL("failed to allocate space for a", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
  
Packit 67cb25
  state->b = (double *) malloc (size * sizeof (double));
Packit 67cb25
  
Packit 67cb25
  if (state->b == NULL)
Packit 67cb25
    {
Packit 67cb25
      steffen_free(state);
Packit 67cb25
      GSL_ERROR_NULL("failed to allocate space for b", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
  
Packit 67cb25
  state->c = (double *) malloc (size * sizeof (double));
Packit 67cb25
  
Packit 67cb25
  if (state->c == NULL)
Packit 67cb25
    {
Packit 67cb25
      steffen_free(state);
Packit 67cb25
      GSL_ERROR_NULL("failed to allocate space for c", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
  
Packit 67cb25
  state->d = (double *) malloc (size * sizeof (double));
Packit 67cb25
  
Packit 67cb25
  if (state->d == NULL)
Packit 67cb25
    {
Packit 67cb25
      steffen_free(state);
Packit 67cb25
      GSL_ERROR_NULL("failed to allocate space for d", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  state->y_prime = (double *) malloc (size * sizeof (double));
Packit 67cb25
  if (state->y_prime == NULL)
Packit 67cb25
    {
Packit 67cb25
      steffen_free(state);
Packit 67cb25
      GSL_ERROR_NULL("failed to allocate space for y_prime", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  return state;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
steffen_init (void * vstate, const double x_array[],
Packit 67cb25
              const double y_array[], size_t size)
Packit 67cb25
{
Packit 67cb25
  steffen_state_t *state = (steffen_state_t *) vstate;
Packit 67cb25
  size_t i;
Packit 67cb25
  double *a = state->a;
Packit 67cb25
  double *b = state->b;
Packit 67cb25
  double *c = state->c;
Packit 67cb25
  double *d = state->d;
Packit 67cb25
  double *y_prime = state->y_prime;
Packit 67cb25
Packit 67cb25
  /*
Packit 67cb25
   * first assign the interval and slopes for the left boundary.
Packit 67cb25
   * We use the "simplest possibility" method described in the paper
Packit 67cb25
   * in section 2.2
Packit 67cb25
   */
Packit 67cb25
  double h0 = (x_array[1] - x_array[0]);
Packit 67cb25
  double s0 = (y_array[1] - y_array[0]) / h0;
Packit 67cb25
Packit 67cb25
  y_prime[0] = s0;
Packit 67cb25
Packit 67cb25
  /* Now we calculate all the necessary s, h, p, and y' variables 
Packit 67cb25
     from 1 to N-2 (0 to size - 2 inclusive) */
Packit 67cb25
  for (i = 1; i < (size - 1); i++)
Packit 67cb25
    {
Packit 67cb25
      double pi;
Packit 67cb25
Packit 67cb25
      /* equation 6 in the paper */
Packit 67cb25
      double hi = (x_array[i+1] - x_array[i]);
Packit 67cb25
      double him1 = (x_array[i] - x_array[i - 1]);
Packit 67cb25
Packit 67cb25
      /* equation 7 in the paper */
Packit 67cb25
      double si = (y_array[i+1] - y_array[i]) / hi;
Packit 67cb25
      double sim1 = (y_array[i] - y_array[i - 1]) / him1;
Packit 67cb25
Packit 67cb25
      /* equation 8 in the paper */
Packit 67cb25
      pi = (sim1*hi + si*him1) / (him1 + hi);
Packit 67cb25
Packit 67cb25
      /* This is a C equivalent of the FORTRAN statement below eqn 11 */
Packit 67cb25
      y_prime[i] = (steffen_copysign(1.0,sim1) + steffen_copysign(1.0,si)) *
Packit 67cb25
                    GSL_MIN(fabs(sim1),
Packit 67cb25
                            GSL_MIN(fabs(si), 0.5*fabs(pi))); 
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  /*
Packit 67cb25
   * we also need y' for the rightmost boundary; we use the
Packit 67cb25
   * "simplest possibility" method described in the paper in
Packit 67cb25
   * section 2.2
Packit 67cb25
   *
Packit 67cb25
   * y' = s_{n-1}
Packit 67cb25
   */
Packit 67cb25
  y_prime[size-1] = (y_array[size - 1] - y_array[size - 2]) /
Packit 67cb25
                    (x_array[size - 1] - x_array[size - 2]);
Packit 67cb25
Packit 67cb25
  /* Now we can calculate all the coefficients for the whole range. */
Packit 67cb25
  for (i = 0; i < (size - 1); i++)
Packit 67cb25
    {
Packit 67cb25
      double hi = (x_array[i+1] - x_array[i]);
Packit 67cb25
      double si = (y_array[i+1] - y_array[i]) / hi;
Packit 67cb25
Packit 67cb25
      /* These are from equations 2-5 in the paper. */
Packit 67cb25
      a[i] = (y_prime[i] + y_prime[i+1] - 2*si) / hi / hi;
Packit 67cb25
      b[i] = (3*si - 2*y_prime[i] - y_prime[i+1]) / hi;
Packit 67cb25
      c[i] = y_prime[i];
Packit 67cb25
      d[i] = y_array[i];
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static void
Packit 67cb25
steffen_free (void * vstate)
Packit 67cb25
{
Packit 67cb25
  steffen_state_t *state = (steffen_state_t *) vstate;
Packit 67cb25
Packit 67cb25
  RETURN_IF_NULL(state);
Packit 67cb25
Packit 67cb25
  if (state->a)
Packit 67cb25
    free (state->a);
Packit 67cb25
Packit 67cb25
  if (state->b)
Packit 67cb25
    free (state->b);
Packit 67cb25
Packit 67cb25
  if (state->c)
Packit 67cb25
    free (state->c);
Packit 67cb25
Packit 67cb25
  if (state->d)
Packit 67cb25
    free (state->d);
Packit 67cb25
Packit 67cb25
  if (state->y_prime)
Packit 67cb25
    free (state->y_prime);
Packit 67cb25
Packit 67cb25
  free (state);
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
steffen_eval (const void * vstate,
Packit 67cb25
              const double x_array[], const double y_array[], size_t size,
Packit 67cb25
              double x, gsl_interp_accel * a, double *y)
Packit 67cb25
{
Packit 67cb25
  const steffen_state_t *state = (const steffen_state_t *) vstate;
Packit 67cb25
Packit 67cb25
  size_t index;
Packit 67cb25
  
Packit 67cb25
  if (a != 0)
Packit 67cb25
    {
Packit 67cb25
      index = gsl_interp_accel_find (a, x_array, size, x);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      index = gsl_interp_bsearch (x_array, x, 0, size - 1);
Packit 67cb25
    }
Packit 67cb25
  
Packit 67cb25
  /* evaluate */
Packit 67cb25
  {
Packit 67cb25
    const double x_lo = x_array[index];
Packit 67cb25
    const double delx = x - x_lo;
Packit 67cb25
    const double a = state->a[index];
Packit 67cb25
    const double b = state->b[index];
Packit 67cb25
    const double c = state->c[index];
Packit 67cb25
    const double d = state->d[index];
Packit 67cb25
    /* Use Horner's scheme for efficient evaluation of polynomials. */
Packit 67cb25
    /* *y = a*delx*delx*delx + b*delx*delx + c*delx + d; */
Packit 67cb25
    *y = d + delx*(c + delx*(b + delx*a));
Packit 67cb25
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
steffen_eval_deriv (const void * vstate,
Packit 67cb25
                    const double x_array[], const double y_array[], size_t size,
Packit 67cb25
                    double x, gsl_interp_accel * a, double *dydx)
Packit 67cb25
{
Packit 67cb25
  const steffen_state_t *state = (const steffen_state_t *) vstate;
Packit 67cb25
Packit 67cb25
  size_t index;
Packit 67cb25
Packit 67cb25
  /* DISCARD_POINTER(y_array); /\* prevent warning about unused parameter *\/ */
Packit 67cb25
  
Packit 67cb25
  if (a != 0)
Packit 67cb25
    {
Packit 67cb25
      index = gsl_interp_accel_find (a, x_array, size, x);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      index = gsl_interp_bsearch (x_array, x, 0, size - 1);
Packit 67cb25
    }
Packit 67cb25
  
Packit 67cb25
  /* evaluate */
Packit 67cb25
  {
Packit 67cb25
    double x_lo = x_array[index];
Packit 67cb25
    double delx = x - x_lo;
Packit 67cb25
    double a = state->a[index];
Packit 67cb25
    double b = state->b[index];
Packit 67cb25
    double c = state->c[index];
Packit 67cb25
    /*double d = state->d[index];*/
Packit 67cb25
    /* *dydx = 3*a*delx*delx*delx + 2*b*delx + c; */
Packit 67cb25
    *dydx = c + delx*(2*b + delx*3*a);
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
steffen_eval_deriv2 (const void * vstate,
Packit 67cb25
                     const double x_array[], const double y_array[], size_t size,
Packit 67cb25
                     double x, gsl_interp_accel * a, double *y_pp)
Packit 67cb25
{
Packit 67cb25
  const steffen_state_t *state = (const steffen_state_t *) vstate;
Packit 67cb25
Packit 67cb25
  size_t index;
Packit 67cb25
Packit 67cb25
  /* DISCARD_POINTER(y_array); /\* prevent warning about unused parameter *\/ */
Packit 67cb25
Packit 67cb25
  if (a != 0)
Packit 67cb25
    {
Packit 67cb25
      index = gsl_interp_accel_find (a, x_array, size, x);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      index = gsl_interp_bsearch (x_array, x, 0, size - 1);
Packit 67cb25
    }
Packit 67cb25
  
Packit 67cb25
  /* evaluate */
Packit 67cb25
  {
Packit 67cb25
    const double x_lo = x_array[index];
Packit 67cb25
    const double delx = x - x_lo;
Packit 67cb25
    const double a = state->a[index];
Packit 67cb25
    const double b = state->b[index];
Packit 67cb25
    *y_pp = 6*a*delx + 2*b;
Packit 67cb25
    return GSL_SUCCESS;
Packit 67cb25
  }
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static int
Packit 67cb25
steffen_eval_integ (const void * vstate,
Packit 67cb25
                    const double x_array[], const double y_array[], size_t size,
Packit 67cb25
                    gsl_interp_accel * acc, double a, double b,
Packit 67cb25
                    double * result)
Packit 67cb25
{
Packit 67cb25
  /* a and b are the boundaries of the integration. */
Packit 67cb25
  
Packit 67cb25
  const steffen_state_t *state = (const steffen_state_t *) vstate;
Packit 67cb25
Packit 67cb25
  size_t i, index_a, index_b;
Packit 67cb25
Packit 67cb25
  /* Find the data points in the x_array that are nearest to the desired */
Packit 67cb25
  /* a and b integration boundaries. */
Packit 67cb25
Packit 67cb25
  if (acc != 0)
Packit 67cb25
    {
Packit 67cb25
      index_a = gsl_interp_accel_find (acc, x_array, size, a);
Packit 67cb25
      index_b = gsl_interp_accel_find (acc, x_array, size, b);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      index_a = gsl_interp_bsearch (x_array, a, 0, size - 1);
Packit 67cb25
      index_b = gsl_interp_bsearch (x_array, b, 0, size - 1);
Packit 67cb25
    }
Packit 67cb25
  
Packit 67cb25
  *result = 0.0;
Packit 67cb25
Packit 67cb25
  /* Iterate over all the segments between data points and sum the */
Packit 67cb25
  /* contributions into result. */
Packit 67cb25
  for(i=index_a; i<=index_b; i++) 
Packit 67cb25
    {
Packit 67cb25
      const double x_hi = x_array[i + 1];
Packit 67cb25
      const double x_lo = x_array[i];
Packit 67cb25
      const double dx = x_hi - x_lo;
Packit 67cb25
      if(dx != 0.0) 
Packit 67cb25
        {
Packit 67cb25
          /*
Packit 67cb25
           * check if we are at a boundary point, so take the
Packit 67cb25
           * a and b parameters instead of the data points.
Packit 67cb25
           */
Packit 67cb25
          double x1 = (i == index_a) ? a-x_lo : 0.0;
Packit 67cb25
          double x2 = (i == index_b) ? b-x_lo : x_hi-x_lo;
Packit 67cb25
Packit 67cb25
          *result += (1.0/4.0)*state->a[i]*(x2*x2*x2*x2 - x1*x1*x1*x1)
Packit 67cb25
                    +(1.0/3.0)*state->b[i]*(x2*x2*x2 - x1*x1*x1)
Packit 67cb25
                    +(1.0/2.0)*state->c[i]*(x2*x2 - x1*x1)
Packit 67cb25
                    +state->d[i]*(x2-x1);
Packit 67cb25
        }
Packit 67cb25
      else /* if the interval was zero, i.e. consecutive x values in data. */
Packit 67cb25
        {
Packit 67cb25
          *result = 0.0;
Packit 67cb25
          return GSL_EINVAL;
Packit 67cb25
        }
Packit 67cb25
    }
Packit 67cb25
  
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static double
Packit 67cb25
steffen_copysign(const double x, const double y)
Packit 67cb25
{
Packit 67cb25
  if ((x < 0 && y > 0) || (x > 0 && y < 0))
Packit 67cb25
    return -x;
Packit 67cb25
Packit 67cb25
  return x;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static const gsl_interp_type steffen_type = 
Packit 67cb25
{
Packit 67cb25
  "steffen", 
Packit 67cb25
  3,
Packit 67cb25
  &steffen_alloc,
Packit 67cb25
  &steffen_init,
Packit 67cb25
  &steffen_eval,
Packit 67cb25
  &steffen_eval_deriv,
Packit 67cb25
  &steffen_eval_deriv2,
Packit 67cb25
  &steffen_eval_integ,
Packit 67cb25
  &steffen_free
Packit 67cb25
};
Packit 67cb25
Packit 67cb25
const gsl_interp_type * gsl_interp_steffen = &steffen_type;