|
Packit |
67cb25 |
/* interpolation/bicubic.c
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* Copyright 2012 David Zaslavsky
|
|
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 |
#include <config.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_errno.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_spline.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_vector.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_interp.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_interp2d.h>
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
#define IDX2D(i, j, w) ((j) * ((w)->xsize) + (i))
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
typedef struct
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
double * zx;
|
|
Packit |
67cb25 |
double * zy;
|
|
Packit |
67cb25 |
double * zxy;
|
|
Packit |
67cb25 |
size_t xsize;
|
|
Packit |
67cb25 |
size_t ysize;
|
|
Packit |
67cb25 |
} bicubic_state_t;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void bicubic_free (void * vstate);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void *
|
|
Packit |
67cb25 |
bicubic_alloc(size_t xsize, size_t ysize)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
bicubic_state_t *state;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state = calloc(1, sizeof (bicubic_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->zx = (double *) malloc (xsize * ysize * sizeof (double));
|
|
Packit |
67cb25 |
if (state->zx == NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
bicubic_free(state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL("failed to allocate space for zx", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->zy = (double *) malloc (xsize * ysize * sizeof (double));
|
|
Packit |
67cb25 |
if (state->zy == NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
bicubic_free(state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL("failed to allocate space for zy", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->zxy = (double *) malloc (xsize * ysize * sizeof (double));
|
|
Packit |
67cb25 |
if (state->zxy == NULL)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
bicubic_free(state);
|
|
Packit |
67cb25 |
GSL_ERROR_NULL("failed to allocate space for zxy", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
state->xsize = xsize;
|
|
Packit |
67cb25 |
state->ysize = ysize;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return state;
|
|
Packit |
67cb25 |
} /* bicubic_alloc() */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void
|
|
Packit |
67cb25 |
bicubic_free (void * vstate)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
bicubic_state_t *state = (bicubic_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
RETURN_IF_NULL(state);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->zx)
|
|
Packit |
67cb25 |
free (state->zx);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->zy)
|
|
Packit |
67cb25 |
free (state->zy);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (state->zxy)
|
|
Packit |
67cb25 |
free (state->zxy);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
free (state);
|
|
Packit |
67cb25 |
} /* bicubic_free() */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
bicubic_init(void * vstate, const double xa[], const double ya[],
|
|
Packit |
67cb25 |
const double za[], size_t xsize, size_t ysize)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
bicubic_state_t *state = (bicubic_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_interp_accel *acc = gsl_interp_accel_alloc();
|
|
Packit |
67cb25 |
gsl_spline *spline;
|
|
Packit |
67cb25 |
gsl_vector *x;
|
|
Packit |
67cb25 |
gsl_vector *y;
|
|
Packit |
67cb25 |
size_t i, j;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
x = gsl_vector_alloc(xsize);
|
|
Packit |
67cb25 |
y = gsl_vector_alloc(xsize);
|
|
Packit |
67cb25 |
spline = gsl_spline_alloc(gsl_interp_cspline, xsize);
|
|
Packit |
67cb25 |
for (j = 0; j <= ysize - 1; j++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
for (i = 0; i <= xsize - 1; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_vector_set(x, i, xa[i]);
|
|
Packit |
67cb25 |
gsl_vector_set(y, i, za[IDX2D(i, j, state)]);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
gsl_spline_init(spline, x->data, y->data, xsize);
|
|
Packit |
67cb25 |
for (i = 0; i <= xsize - 1; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
state->zx[IDX2D(i, j, state)] = gsl_spline_eval_deriv(spline, xa[i], acc);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
gsl_vector_free(x);
|
|
Packit |
67cb25 |
gsl_vector_free(y);
|
|
Packit |
67cb25 |
gsl_spline_free(spline);
|
|
Packit |
67cb25 |
gsl_interp_accel_reset(acc);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
x = gsl_vector_alloc(ysize);
|
|
Packit |
67cb25 |
y = gsl_vector_alloc(ysize);
|
|
Packit |
67cb25 |
spline = gsl_spline_alloc(gsl_interp_cspline, ysize);
|
|
Packit |
67cb25 |
for (i = 0; i <= xsize - 1; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
for (j = 0; j <= ysize - 1; j++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_vector_set(x, j, ya[j]);
|
|
Packit |
67cb25 |
gsl_vector_set(y, j, za[IDX2D(i, j, state)]);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
gsl_spline_init(spline, x->data, y->data, ysize);
|
|
Packit |
67cb25 |
for (j = 0; j <= ysize - 1; j++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
state->zy[IDX2D(i, j, state)] = gsl_spline_eval_deriv(spline, ya[j], acc);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
gsl_vector_free(x);
|
|
Packit |
67cb25 |
gsl_vector_free(y);
|
|
Packit |
67cb25 |
gsl_spline_free(spline);
|
|
Packit |
67cb25 |
gsl_interp_accel_reset(acc);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
x = gsl_vector_alloc(xsize);
|
|
Packit |
67cb25 |
y = gsl_vector_alloc(xsize);
|
|
Packit |
67cb25 |
spline = gsl_spline_alloc(gsl_interp_cspline, xsize);
|
|
Packit |
67cb25 |
for (j = 0; j <= ysize - 1; j++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
for (i = 0; i <= xsize - 1; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_vector_set(x, i, xa[i]);
|
|
Packit |
67cb25 |
gsl_vector_set(y, i, state->zy[IDX2D(i, j, state)]);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
gsl_spline_init(spline, x->data, y->data, xsize);
|
|
Packit |
67cb25 |
for (i = 0; i <= xsize - 1; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
state->zxy[IDX2D(i, j, state)] = gsl_spline_eval_deriv(spline, xa[i], acc);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
gsl_vector_free(x);
|
|
Packit |
67cb25 |
gsl_vector_free(y);
|
|
Packit |
67cb25 |
gsl_spline_free(spline);
|
|
Packit |
67cb25 |
gsl_interp_accel_free(acc);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
} /* bicubic_init() */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
bicubic_eval(const void * vstate, const double xarr[], const double yarr[],
|
|
Packit |
67cb25 |
const double zarr[], size_t xsize, size_t ysize,
|
|
Packit |
67cb25 |
double x, double y, gsl_interp_accel * xa,
|
|
Packit |
67cb25 |
gsl_interp_accel * ya, double * z)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
bicubic_state_t *state = (bicubic_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double xmin, xmax, ymin, ymax;
|
|
Packit |
67cb25 |
double zminmin, zminmax, zmaxmin, zmaxmax;
|
|
Packit |
67cb25 |
double zxminmin, zxminmax, zxmaxmin, zxmaxmax;
|
|
Packit |
67cb25 |
double zyminmin, zyminmax, zymaxmin, zymaxmax;
|
|
Packit |
67cb25 |
double zxyminmin, zxyminmax, zxymaxmin, zxymaxmax;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double dx, dy; /* size of the grid cell */
|
|
Packit |
67cb25 |
double dt, du;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
* t and u are the positions within the grid cell at which we are computing
|
|
Packit |
67cb25 |
* the interpolation, in units of grid cell size
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
double t, u;
|
|
Packit |
67cb25 |
double t0, t1, t2, t3, u0, u1, u2, u3;
|
|
Packit |
67cb25 |
double v;
|
|
Packit |
67cb25 |
size_t xi, yi;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* first compute the indices into the data arrays where we are interpolating */
|
|
Packit |
67cb25 |
if (xa != NULL)
|
|
Packit |
67cb25 |
xi = gsl_interp_accel_find(xa, xarr, xsize, x);
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
xi = gsl_interp_bsearch(xarr, x, 0, xsize - 1);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (ya != NULL)
|
|
Packit |
67cb25 |
yi = gsl_interp_accel_find(ya, yarr, ysize, y);
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
yi = gsl_interp_bsearch(yarr, y, 0, ysize - 1);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* find the minimum and maximum values on the grid cell in each dimension */
|
|
Packit |
67cb25 |
xmin = xarr[xi];
|
|
Packit |
67cb25 |
xmax = xarr[xi + 1];
|
|
Packit |
67cb25 |
ymin = yarr[yi];
|
|
Packit |
67cb25 |
ymax = yarr[yi + 1];
|
|
Packit |
67cb25 |
zminmin = zarr[IDX2D(xi, yi, state)];
|
|
Packit |
67cb25 |
zminmax = zarr[IDX2D(xi, yi + 1, state)];
|
|
Packit |
67cb25 |
zmaxmin = zarr[IDX2D(xi + 1, yi, state)];
|
|
Packit |
67cb25 |
zmaxmax = zarr[IDX2D(xi + 1, yi + 1, state)];
|
|
Packit |
67cb25 |
/* Get the width and height of the grid cell */
|
|
Packit |
67cb25 |
dx = xmax - xmin;
|
|
Packit |
67cb25 |
dy = ymax - ymin;
|
|
Packit |
67cb25 |
t = (x - xmin)/dx;
|
|
Packit |
67cb25 |
u = (y - ymin)/dy;
|
|
Packit |
67cb25 |
dt = 1./dx; /* partial t / partial x */
|
|
Packit |
67cb25 |
du = 1./dy; /* partial u / partial y */
|
|
Packit |
67cb25 |
zxminmin = state->zx[IDX2D(xi, yi, state)]/dt;
|
|
Packit |
67cb25 |
zxminmax = state->zx[IDX2D(xi, yi + 1, state)]/dt;
|
|
Packit |
67cb25 |
zxmaxmin = state->zx[IDX2D(xi + 1, yi, state)]/dt;
|
|
Packit |
67cb25 |
zxmaxmax = state->zx[IDX2D(xi + 1, yi + 1, state)]/dt;
|
|
Packit |
67cb25 |
zyminmin = state->zy[IDX2D(xi, yi, state)]/du;
|
|
Packit |
67cb25 |
zyminmax = state->zy[IDX2D(xi, yi + 1, state)]/du;
|
|
Packit |
67cb25 |
zymaxmin = state->zy[IDX2D(xi + 1, yi, state)]/du;
|
|
Packit |
67cb25 |
zymaxmax = state->zy[IDX2D(xi + 1, yi + 1, state)]/du;
|
|
Packit |
67cb25 |
zxyminmin = state->zxy[IDX2D(xi, yi, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxyminmax = state->zxy[IDX2D(xi, yi + 1, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxymaxmin = state->zxy[IDX2D(xi + 1, yi, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxymaxmax = state->zxy[IDX2D(xi + 1, yi + 1, state)]/(dt*du);
|
|
Packit |
67cb25 |
t0 = 1;
|
|
Packit |
67cb25 |
t1 = t;
|
|
Packit |
67cb25 |
t2 = t*t;
|
|
Packit |
67cb25 |
t3 = t*t2;
|
|
Packit |
67cb25 |
u0 = 1;
|
|
Packit |
67cb25 |
u1 = u;
|
|
Packit |
67cb25 |
u2 = u*u;
|
|
Packit |
67cb25 |
u3 = u*u2;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
*z = 0;
|
|
Packit |
67cb25 |
v = zminmin;
|
|
Packit |
67cb25 |
*z += v*t0*u0;
|
|
Packit |
67cb25 |
v = zyminmin;
|
|
Packit |
67cb25 |
*z += v*t0*u1;
|
|
Packit |
67cb25 |
v = -3*zminmin + 3*zminmax - 2*zyminmin - zyminmax;
|
|
Packit |
67cb25 |
*z += v*t0*u2;
|
|
Packit |
67cb25 |
v = 2*zminmin - 2*zminmax + zyminmin + zyminmax;
|
|
Packit |
67cb25 |
*z += v*t0*u3;
|
|
Packit |
67cb25 |
v = zxminmin;
|
|
Packit |
67cb25 |
*z += v*t1*u0;
|
|
Packit |
67cb25 |
v = zxyminmin;
|
|
Packit |
67cb25 |
*z += v*t1*u1;
|
|
Packit |
67cb25 |
v = -3*zxminmin + 3*zxminmax - 2*zxyminmin - zxyminmax;
|
|
Packit |
67cb25 |
*z += v*t1*u2;
|
|
Packit |
67cb25 |
v = 2*zxminmin - 2*zxminmax + zxyminmin + zxyminmax;
|
|
Packit |
67cb25 |
*z += v*t1*u3;
|
|
Packit |
67cb25 |
v = -3*zminmin + 3*zmaxmin - 2*zxminmin - zxmaxmin;
|
|
Packit |
67cb25 |
*z += v*t2*u0;
|
|
Packit |
67cb25 |
v = -3*zyminmin + 3*zymaxmin - 2*zxyminmin - zxymaxmin;
|
|
Packit |
67cb25 |
*z += v*t2*u1;
|
|
Packit |
67cb25 |
v = 9*zminmin - 9*zmaxmin + 9*zmaxmax - 9*zminmax + 6*zxminmin + 3*zxmaxmin - 3*zxmaxmax - 6*zxminmax + 6*zyminmin - 6*zymaxmin - 3*zymaxmax + 3*zyminmax + 4*zxyminmin + 2*zxymaxmin + zxymaxmax + 2*zxyminmax;
|
|
Packit |
67cb25 |
*z += v*t2*u2;
|
|
Packit |
67cb25 |
v = -6*zminmin + 6*zmaxmin - 6*zmaxmax + 6*zminmax - 4*zxminmin - 2*zxmaxmin + 2*zxmaxmax + 4*zxminmax - 3*zyminmin + 3*zymaxmin + 3*zymaxmax - 3*zyminmax - 2*zxyminmin - zxymaxmin - zxymaxmax - 2*zxyminmax;
|
|
Packit |
67cb25 |
*z += v*t2*u3;
|
|
Packit |
67cb25 |
v = 2*zminmin - 2*zmaxmin + zxminmin + zxmaxmin;
|
|
Packit |
67cb25 |
*z += v*t3*u0;
|
|
Packit |
67cb25 |
v = 2*zyminmin - 2*zymaxmin + zxyminmin + zxymaxmin;
|
|
Packit |
67cb25 |
*z += v*t3*u1;
|
|
Packit |
67cb25 |
v = -6*zminmin + 6*zmaxmin - 6*zmaxmax + 6*zminmax - 3*zxminmin - 3*zxmaxmin + 3*zxmaxmax + 3*zxminmax - 4*zyminmin + 4*zymaxmin + 2*zymaxmax - 2*zyminmax - 2*zxyminmin - 2*zxymaxmin - zxymaxmax - zxyminmax;
|
|
Packit |
67cb25 |
*z += v*t3*u2;
|
|
Packit |
67cb25 |
v = 4*zminmin - 4*zmaxmin + 4*zmaxmax - 4*zminmax + 2*zxminmin + 2*zxmaxmin - 2*zxmaxmax - 2*zxminmax + 2*zyminmin - 2*zymaxmin - 2*zymaxmax + 2*zyminmax + zxyminmin + zxymaxmin + zxymaxmax + zxyminmax;
|
|
Packit |
67cb25 |
*z += v*t3*u3;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
} /* bicubic_eval() */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
bicubic_deriv_x(const void * vstate, const double xarr[], const double yarr[],
|
|
Packit |
67cb25 |
const double zarr[], size_t xsize, size_t ysize,
|
|
Packit |
67cb25 |
double x, double y,
|
|
Packit |
67cb25 |
gsl_interp_accel * xa, gsl_interp_accel * ya, double * z_p)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
bicubic_state_t *state = (bicubic_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double xmin, xmax, ymin, ymax;
|
|
Packit |
67cb25 |
double zminmin, zminmax, zmaxmin, zmaxmax;
|
|
Packit |
67cb25 |
double zxminmin, zxminmax, zxmaxmin, zxmaxmax;
|
|
Packit |
67cb25 |
double zyminmin, zyminmax, zymaxmin, zymaxmax;
|
|
Packit |
67cb25 |
double zxyminmin, zxyminmax, zxymaxmin, zxymaxmax;
|
|
Packit |
67cb25 |
double dx, dy; /* size of the grid cell */
|
|
Packit |
67cb25 |
double dt, du;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
* t and u are the positions within the grid cell at which we are computing
|
|
Packit |
67cb25 |
* the interpolation, in units of grid cell size
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
double t, u;
|
|
Packit |
67cb25 |
double t0, t1, t2, u0, u1, u2, u3;
|
|
Packit |
67cb25 |
double v;
|
|
Packit |
67cb25 |
size_t xi, yi;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* first compute the indices into the data arrays where we are interpolating */
|
|
Packit |
67cb25 |
if (xa != NULL)
|
|
Packit |
67cb25 |
xi = gsl_interp_accel_find(xa, xarr, xsize, x);
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
xi = gsl_interp_bsearch(xarr, x, 0, xsize - 1);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (ya != NULL)
|
|
Packit |
67cb25 |
yi = gsl_interp_accel_find(ya, yarr, ysize, y);
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
yi = gsl_interp_bsearch(yarr, y, 0, ysize - 1);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* find the minimum and maximum values on the grid cell in each dimension */
|
|
Packit |
67cb25 |
xmin = xarr[xi];
|
|
Packit |
67cb25 |
xmax = xarr[xi + 1];
|
|
Packit |
67cb25 |
ymin = yarr[yi];
|
|
Packit |
67cb25 |
ymax = yarr[yi + 1];
|
|
Packit |
67cb25 |
zminmin = zarr[IDX2D(xi, yi, state)];
|
|
Packit |
67cb25 |
zminmax = zarr[IDX2D(xi, yi + 1, state)];
|
|
Packit |
67cb25 |
zmaxmin = zarr[IDX2D(xi + 1, yi, state)];
|
|
Packit |
67cb25 |
zmaxmax = zarr[IDX2D(xi + 1, yi + 1, state)];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* get the width and height of the grid cell */
|
|
Packit |
67cb25 |
dx = xmax - xmin;
|
|
Packit |
67cb25 |
dy = ymax - ymin;
|
|
Packit |
67cb25 |
t = (x - xmin)/dx;
|
|
Packit |
67cb25 |
u = (y - ymin)/dy;
|
|
Packit |
67cb25 |
dt = 1./dx; /* partial t / partial x */
|
|
Packit |
67cb25 |
du = 1./dy; /* partial u / partial y */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
zxminmin = state->zx[IDX2D(xi, yi, state)]/dt;
|
|
Packit |
67cb25 |
zxminmax = state->zx[IDX2D(xi, yi + 1, state)]/dt;
|
|
Packit |
67cb25 |
zxmaxmin = state->zx[IDX2D(xi + 1, yi, state)]/dt;
|
|
Packit |
67cb25 |
zxmaxmax = state->zx[IDX2D(xi + 1, yi + 1, state)]/dt;
|
|
Packit |
67cb25 |
zyminmin = state->zy[IDX2D(xi, yi, state)]/du;
|
|
Packit |
67cb25 |
zyminmax = state->zy[IDX2D(xi, yi + 1, state)]/du;
|
|
Packit |
67cb25 |
zymaxmin = state->zy[IDX2D(xi + 1, yi, state)]/du;
|
|
Packit |
67cb25 |
zymaxmax = state->zy[IDX2D(xi + 1, yi + 1, state)]/du;
|
|
Packit |
67cb25 |
zxyminmin = state->zxy[IDX2D(xi, yi, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxyminmax = state->zxy[IDX2D(xi, yi + 1, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxymaxmin = state->zxy[IDX2D(xi + 1, yi, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxymaxmax = state->zxy[IDX2D(xi + 1, yi + 1, state)]/(dt*du);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
t0 = 1;
|
|
Packit |
67cb25 |
t1 = t;
|
|
Packit |
67cb25 |
t2 = t*t;
|
|
Packit |
67cb25 |
u0 = 1;
|
|
Packit |
67cb25 |
u1 = u;
|
|
Packit |
67cb25 |
u2 = u*u;
|
|
Packit |
67cb25 |
u3 = u*u2;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
*z_p = 0;
|
|
Packit |
67cb25 |
v = zxminmin;
|
|
Packit |
67cb25 |
*z_p += v*t0*u0;
|
|
Packit |
67cb25 |
v = zxyminmin;
|
|
Packit |
67cb25 |
*z_p += v*t0*u1;
|
|
Packit |
67cb25 |
v = -3*zxminmin + 3*zxminmax - 2*zxyminmin - zxyminmax;
|
|
Packit |
67cb25 |
*z_p += v*t0*u2;
|
|
Packit |
67cb25 |
v = 2*zxminmin - 2*zxminmax + zxyminmin + zxyminmax;
|
|
Packit |
67cb25 |
*z_p += v*t0*u3;
|
|
Packit |
67cb25 |
v = -3*zminmin + 3*zmaxmin - 2*zxminmin - zxmaxmin;
|
|
Packit |
67cb25 |
*z_p += 2*v*t1*u0;
|
|
Packit |
67cb25 |
v = -3*zyminmin + 3*zymaxmin - 2*zxyminmin - zxymaxmin;
|
|
Packit |
67cb25 |
*z_p += 2*v*t1*u1;
|
|
Packit |
67cb25 |
v = 9*zminmin - 9*zmaxmin + 9*zmaxmax - 9*zminmax + 6*zxminmin + 3*zxmaxmin - 3*zxmaxmax - 6*zxminmax + 6*zyminmin - 6*zymaxmin - 3*zymaxmax + 3*zyminmax + 4*zxyminmin + 2*zxymaxmin + zxymaxmax + 2*zxyminmax;
|
|
Packit |
67cb25 |
*z_p += 2*v*t1*u2;
|
|
Packit |
67cb25 |
v = -6*zminmin + 6*zmaxmin - 6*zmaxmax + 6*zminmax - 4*zxminmin - 2*zxmaxmin + 2*zxmaxmax + 4*zxminmax - 3*zyminmin + 3*zymaxmin + 3*zymaxmax - 3*zyminmax - 2*zxyminmin - zxymaxmin - zxymaxmax - 2*zxyminmax;
|
|
Packit |
67cb25 |
*z_p += 2*v*t1*u3;
|
|
Packit |
67cb25 |
v = 2*zminmin - 2*zmaxmin + zxminmin + zxmaxmin;
|
|
Packit |
67cb25 |
*z_p += 3*v*t2*u0;
|
|
Packit |
67cb25 |
v = 2*zyminmin - 2*zymaxmin + zxyminmin + zxymaxmin;
|
|
Packit |
67cb25 |
*z_p += 3*v*t2*u1;
|
|
Packit |
67cb25 |
v = -6*zminmin + 6*zmaxmin - 6*zmaxmax + 6*zminmax - 3*zxminmin - 3*zxmaxmin + 3*zxmaxmax + 3*zxminmax - 4*zyminmin + 4*zymaxmin + 2*zymaxmax - 2*zyminmax - 2*zxyminmin - 2*zxymaxmin - zxymaxmax - zxyminmax;
|
|
Packit |
67cb25 |
*z_p += 3*v*t2*u2;
|
|
Packit |
67cb25 |
v = 4*zminmin - 4*zmaxmin + 4*zmaxmax - 4*zminmax + 2*zxminmin + 2*zxmaxmin - 2*zxmaxmax - 2*zxminmax + 2*zyminmin - 2*zymaxmin - 2*zymaxmax + 2*zyminmax + zxyminmin + zxymaxmin + zxymaxmax + zxyminmax;
|
|
Packit |
67cb25 |
*z_p += 3*v*t2*u3;
|
|
Packit |
67cb25 |
*z_p *= dt;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
} /* bicubic_deriv_x() */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
bicubic_deriv_y(const void * vstate, const double xarr[], const double yarr[],
|
|
Packit |
67cb25 |
const double zarr[], size_t xsize, size_t ysize,
|
|
Packit |
67cb25 |
double x, double y,
|
|
Packit |
67cb25 |
gsl_interp_accel * xa, gsl_interp_accel * ya, double * z_p)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
bicubic_state_t *state = (bicubic_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double xmin, xmax, ymin, ymax;
|
|
Packit |
67cb25 |
double zminmin, zminmax, zmaxmin, zmaxmax;
|
|
Packit |
67cb25 |
double zxminmin, zxminmax, zxmaxmin, zxmaxmax;
|
|
Packit |
67cb25 |
double zyminmin, zyminmax, zymaxmin, zymaxmax;
|
|
Packit |
67cb25 |
double zxyminmin, zxyminmax, zxymaxmin, zxymaxmax;
|
|
Packit |
67cb25 |
/* dx and dy are the size of the grid cell */
|
|
Packit |
67cb25 |
double dx, dy;
|
|
Packit |
67cb25 |
double dt, du;
|
|
Packit |
67cb25 |
/* t and u are the positions within the grid cell at which we are
|
|
Packit |
67cb25 |
* computing the interpolation, in units of grid cell size */
|
|
Packit |
67cb25 |
double t, u;
|
|
Packit |
67cb25 |
double t0, t1, t2, t3, u0, u1, u2;
|
|
Packit |
67cb25 |
double v;
|
|
Packit |
67cb25 |
size_t xi, yi;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* first compute the indices into the data arrays where we are interpolating */
|
|
Packit |
67cb25 |
if (xa != NULL)
|
|
Packit |
67cb25 |
xi = gsl_interp_accel_find(xa, xarr, xsize, x);
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
xi = gsl_interp_bsearch(xarr, x, 0, xsize - 1);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (ya != NULL)
|
|
Packit |
67cb25 |
yi = gsl_interp_accel_find(ya, yarr, ysize, y);
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
yi = gsl_interp_bsearch(yarr, y, 0, ysize - 1);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* find the minimum and maximum values on the grid cell in each dimension */
|
|
Packit |
67cb25 |
xmin = xarr[xi];
|
|
Packit |
67cb25 |
xmax = xarr[xi + 1];
|
|
Packit |
67cb25 |
ymin = yarr[yi];
|
|
Packit |
67cb25 |
ymax = yarr[yi + 1];
|
|
Packit |
67cb25 |
zminmin = zarr[IDX2D(xi, yi, state)];
|
|
Packit |
67cb25 |
zminmax = zarr[IDX2D(xi, yi + 1, state)];
|
|
Packit |
67cb25 |
zmaxmin = zarr[IDX2D(xi + 1, yi, state)];
|
|
Packit |
67cb25 |
zmaxmax = zarr[IDX2D(xi + 1, yi + 1, state)];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* get the width and height of the grid cell */
|
|
Packit |
67cb25 |
dx = xmax - xmin;
|
|
Packit |
67cb25 |
dy = ymax - ymin;
|
|
Packit |
67cb25 |
t = (x - xmin)/dx;
|
|
Packit |
67cb25 |
u = (y - ymin)/dy;
|
|
Packit |
67cb25 |
dt = 1./dx; /* partial t / partial x */
|
|
Packit |
67cb25 |
du = 1./dy; /* partial u / partial y */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
zxminmin = state->zx[IDX2D(xi, yi, state)]/dt;
|
|
Packit |
67cb25 |
zxminmax = state->zx[IDX2D(xi, yi + 1, state)]/dt;
|
|
Packit |
67cb25 |
zxmaxmin = state->zx[IDX2D(xi + 1, yi, state)]/dt;
|
|
Packit |
67cb25 |
zxmaxmax = state->zx[IDX2D(xi + 1, yi + 1, state)]/dt;
|
|
Packit |
67cb25 |
zyminmin = state->zy[IDX2D(xi, yi, state)]/du;
|
|
Packit |
67cb25 |
zyminmax = state->zy[IDX2D(xi, yi + 1, state)]/du;
|
|
Packit |
67cb25 |
zymaxmin = state->zy[IDX2D(xi + 1, yi, state)]/du;
|
|
Packit |
67cb25 |
zymaxmax = state->zy[IDX2D(xi + 1, yi + 1, state)]/du;
|
|
Packit |
67cb25 |
zxyminmin = state->zxy[IDX2D(xi, yi, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxyminmax = state->zxy[IDX2D(xi, yi + 1, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxymaxmin = state->zxy[IDX2D(xi + 1, yi, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxymaxmax = state->zxy[IDX2D(xi + 1, yi + 1, state)]/(dt*du);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
t0 = 1;
|
|
Packit |
67cb25 |
t1 = t;
|
|
Packit |
67cb25 |
t2 = t*t;
|
|
Packit |
67cb25 |
t3 = t*t2;
|
|
Packit |
67cb25 |
u0 = 1;
|
|
Packit |
67cb25 |
u1 = u;
|
|
Packit |
67cb25 |
u2 = u*u;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
*z_p = 0;
|
|
Packit |
67cb25 |
v = zyminmin;
|
|
Packit |
67cb25 |
*z_p += v*t0*u0;
|
|
Packit |
67cb25 |
v = -3*zminmin + 3*zminmax - 2*zyminmin - zyminmax;
|
|
Packit |
67cb25 |
*z_p += 2*v*t0*u1;
|
|
Packit |
67cb25 |
v = 2*zminmin-2*zminmax + zyminmin + zyminmax;
|
|
Packit |
67cb25 |
*z_p += 3*v*t0*u2;
|
|
Packit |
67cb25 |
v = zxyminmin;
|
|
Packit |
67cb25 |
*z_p += v*t1*u0;
|
|
Packit |
67cb25 |
v = -3*zxminmin + 3*zxminmax - 2*zxyminmin - zxyminmax;
|
|
Packit |
67cb25 |
*z_p += 2*v*t1*u1;
|
|
Packit |
67cb25 |
v = 2*zxminmin - 2*zxminmax + zxyminmin + zxyminmax;
|
|
Packit |
67cb25 |
*z_p += 3*v*t1*u2;
|
|
Packit |
67cb25 |
v = -3*zyminmin + 3*zymaxmin - 2*zxyminmin - zxymaxmin;
|
|
Packit |
67cb25 |
*z_p += v*t2*u0;
|
|
Packit |
67cb25 |
v = 9*zminmin - 9*zmaxmin + 9*zmaxmax - 9*zminmax + 6*zxminmin + 3*zxmaxmin - 3*zxmaxmax - 6*zxminmax + 6*zyminmin - 6*zymaxmin - 3*zymaxmax + 3*zyminmax + 4*zxyminmin + 2*zxymaxmin + zxymaxmax + 2*zxyminmax;
|
|
Packit |
67cb25 |
*z_p += 2*v*t2*u1;
|
|
Packit |
67cb25 |
v = -6*zminmin + 6*zmaxmin - 6*zmaxmax + 6*zminmax - 4*zxminmin - 2*zxmaxmin + 2*zxmaxmax + 4*zxminmax - 3*zyminmin + 3*zymaxmin + 3*zymaxmax - 3*zyminmax - 2*zxyminmin - zxymaxmin - zxymaxmax - 2*zxyminmax;
|
|
Packit |
67cb25 |
*z_p += 3*v*t2*u2;
|
|
Packit |
67cb25 |
v = 2*zyminmin - 2*zymaxmin + zxyminmin + zxymaxmin;
|
|
Packit |
67cb25 |
*z_p += v*t3*u0;
|
|
Packit |
67cb25 |
v = -6*zminmin + 6*zmaxmin - 6*zmaxmax + 6*zminmax - 3*zxminmin - 3*zxmaxmin + 3*zxmaxmax + 3*zxminmax - 4*zyminmin + 4*zymaxmin + 2*zymaxmax - 2*zyminmax - 2*zxyminmin - 2*zxymaxmin - zxymaxmax - zxyminmax;
|
|
Packit |
67cb25 |
*z_p += 2*v*t3*u1;
|
|
Packit |
67cb25 |
v = 4*zminmin - 4*zmaxmin + 4*zmaxmax - 4*zminmax + 2*zxminmin + 2*zxmaxmin - 2*zxmaxmax - 2*zxminmax + 2*zyminmin - 2*zymaxmin - 2*zymaxmax + 2*zyminmax + zxyminmin + zxymaxmin + zxymaxmax + zxyminmax;
|
|
Packit |
67cb25 |
*z_p += 3*v*t3*u2;
|
|
Packit |
67cb25 |
*z_p *= du;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
bicubic_deriv_xx(const void * vstate, const double xarr[], const double yarr[],
|
|
Packit |
67cb25 |
const double zarr[], size_t xsize, size_t ysize,
|
|
Packit |
67cb25 |
double x, double y,
|
|
Packit |
67cb25 |
gsl_interp_accel * xa, gsl_interp_accel * ya, double * z_pp)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
bicubic_state_t *state = (bicubic_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double xmin, xmax, ymin, ymax;
|
|
Packit |
67cb25 |
double zminmin, zminmax, zmaxmin, zmaxmax;
|
|
Packit |
67cb25 |
double zxminmin, zxminmax, zxmaxmin, zxmaxmax;
|
|
Packit |
67cb25 |
double zyminmin, zyminmax, zymaxmin, zymaxmax;
|
|
Packit |
67cb25 |
double zxyminmin, zxyminmax, zxymaxmin, zxymaxmax;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double dx, dy; /* size of the grid cell */
|
|
Packit |
67cb25 |
double dt, du;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
* t and u are the positions within the grid cell at which we are computing
|
|
Packit |
67cb25 |
* the interpolation, in units of grid cell size
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
double t, u;
|
|
Packit |
67cb25 |
double t0, t1, u0, u1, u2, u3;
|
|
Packit |
67cb25 |
double v;
|
|
Packit |
67cb25 |
size_t xi, yi;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* first compute the indices into the data arrays where we are interpolating */
|
|
Packit |
67cb25 |
if (xa != NULL)
|
|
Packit |
67cb25 |
xi = gsl_interp_accel_find(xa, xarr, xsize, x);
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
xi = gsl_interp_bsearch(xarr, x, 0, xsize - 1);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (ya != NULL)
|
|
Packit |
67cb25 |
yi = gsl_interp_accel_find(ya, yarr, ysize, y);
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
yi = gsl_interp_bsearch(yarr, y, 0, ysize - 1);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* find the minimum and maximum values on the grid cell in each dimension */
|
|
Packit |
67cb25 |
xmin = xarr[xi];
|
|
Packit |
67cb25 |
xmax = xarr[xi + 1];
|
|
Packit |
67cb25 |
ymin = yarr[yi];
|
|
Packit |
67cb25 |
ymax = yarr[yi + 1];
|
|
Packit |
67cb25 |
zminmin = zarr[IDX2D(xi, yi, state)];
|
|
Packit |
67cb25 |
zminmax = zarr[IDX2D(xi, yi + 1, state)];
|
|
Packit |
67cb25 |
zmaxmin = zarr[IDX2D(xi + 1, yi, state)];
|
|
Packit |
67cb25 |
zmaxmax = zarr[IDX2D(xi + 1, yi + 1, state)];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* get the width and height of the grid cell */
|
|
Packit |
67cb25 |
dx = xmax - xmin;
|
|
Packit |
67cb25 |
dy = ymax - ymin;
|
|
Packit |
67cb25 |
t = (x - xmin)/dx;
|
|
Packit |
67cb25 |
u = (y - ymin)/dy;
|
|
Packit |
67cb25 |
dt = 1./dx; /* partial t / partial x */
|
|
Packit |
67cb25 |
du = 1./dy; /* partial u / partial y */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
zxminmin = state->zx[IDX2D(xi, yi, state)]/dt;
|
|
Packit |
67cb25 |
zxminmax = state->zx[IDX2D(xi, yi + 1, state)]/dt;
|
|
Packit |
67cb25 |
zxmaxmin = state->zx[IDX2D(xi + 1, yi, state)]/dt;
|
|
Packit |
67cb25 |
zxmaxmax = state->zx[IDX2D(xi + 1, yi + 1, state)]/dt;
|
|
Packit |
67cb25 |
zyminmin = state->zy[IDX2D(xi, yi, state)]/du;
|
|
Packit |
67cb25 |
zyminmax = state->zy[IDX2D(xi, yi + 1, state)]/du;
|
|
Packit |
67cb25 |
zymaxmin = state->zy[IDX2D(xi + 1, yi, state)]/du;
|
|
Packit |
67cb25 |
zymaxmax = state->zy[IDX2D(xi + 1, yi + 1, state)]/du;
|
|
Packit |
67cb25 |
zxyminmin = state->zxy[IDX2D(xi, yi, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxyminmax = state->zxy[IDX2D(xi, yi + 1, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxymaxmin = state->zxy[IDX2D(xi + 1, yi, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxymaxmax = state->zxy[IDX2D(xi + 1, yi + 1, state)]/(dt*du);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
t0 = 1;
|
|
Packit |
67cb25 |
t1 = t;
|
|
Packit |
67cb25 |
u0 = 1;
|
|
Packit |
67cb25 |
u1 = u;
|
|
Packit |
67cb25 |
u2 = u*u;
|
|
Packit |
67cb25 |
u3 = u*u2;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
*z_pp = 0;
|
|
Packit |
67cb25 |
v = -3*zminmin + 3*zmaxmin - 2*zxminmin - zxmaxmin;
|
|
Packit |
67cb25 |
*z_pp += 2*v*t0*u0;
|
|
Packit |
67cb25 |
v = -3*zyminmin + 3*zymaxmin - 2*zxyminmin - zxymaxmin;
|
|
Packit |
67cb25 |
*z_pp += 2*v*t0*u1;
|
|
Packit |
67cb25 |
v = 9*zminmin - 9*zmaxmin + 9*zmaxmax - 9*zminmax + 6*zxminmin + 3*zxmaxmin - 3*zxmaxmax - 6*zxminmax + 6*zyminmin - 6*zymaxmin - 3*zymaxmax + 3*zyminmax + 4*zxyminmin + 2*zxymaxmin + zxymaxmax + 2*zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 2*v*t0*u2;
|
|
Packit |
67cb25 |
v = -6*zminmin + 6*zmaxmin - 6*zmaxmax + 6*zminmax - 4*zxminmin - 2*zxmaxmin + 2*zxmaxmax + 4*zxminmax - 3*zyminmin + 3*zymaxmin + 3*zymaxmax - 3*zyminmax - 2*zxyminmin - zxymaxmin - zxymaxmax - 2*zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 2*v*t0*u3;
|
|
Packit |
67cb25 |
v = 2*zminmin - 2*zmaxmin + zxminmin + zxmaxmin;
|
|
Packit |
67cb25 |
*z_pp += 6*v*t1*u0;
|
|
Packit |
67cb25 |
v = 2*zyminmin - 2*zymaxmin + zxyminmin + zxymaxmin;
|
|
Packit |
67cb25 |
*z_pp += 6*v*t1*u1;
|
|
Packit |
67cb25 |
v = -6*zminmin + 6*zmaxmin - 6*zmaxmax + 6*zminmax - 3*zxminmin - 3*zxmaxmin + 3*zxmaxmax + 3*zxminmax - 4*zyminmin + 4*zymaxmin + 2*zymaxmax - 2*zyminmax - 2*zxyminmin - 2*zxymaxmin - zxymaxmax - zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 6*v*t1*u2;
|
|
Packit |
67cb25 |
v = 4*zminmin - 4*zmaxmin + 4*zmaxmax - 4*zminmax + 2*zxminmin + 2*zxmaxmin - 2*zxmaxmax - 2*zxminmax + 2*zyminmin - 2*zymaxmin - 2*zymaxmax + 2*zyminmax + zxyminmin + zxymaxmin + zxymaxmax + zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 6*v*t1*u3;
|
|
Packit |
67cb25 |
*z_pp *= dt*dt;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
bicubic_deriv_xy(const void * vstate, const double xarr[], const double yarr[],
|
|
Packit |
67cb25 |
const double zarr[], size_t xsize, size_t ysize,
|
|
Packit |
67cb25 |
double x, double y,
|
|
Packit |
67cb25 |
gsl_interp_accel * xa, gsl_interp_accel * ya, double * z_pp)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
bicubic_state_t *state = (bicubic_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double xmin, xmax, ymin, ymax;
|
|
Packit |
67cb25 |
double zminmin, zminmax, zmaxmin, zmaxmax;
|
|
Packit |
67cb25 |
double zxminmin, zxminmax, zxmaxmin, zxmaxmax;
|
|
Packit |
67cb25 |
double zyminmin, zyminmax, zymaxmin, zymaxmax;
|
|
Packit |
67cb25 |
double zxyminmin, zxyminmax, zxymaxmin, zxymaxmax;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double dx, dy; /* size of the grid cell */
|
|
Packit |
67cb25 |
double dt, du;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
* t and u are the positions within the grid cell at which we are computing
|
|
Packit |
67cb25 |
* the interpolation, in units of grid cell size
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
double t, u;
|
|
Packit |
67cb25 |
double t0, t1, t2, u0, u1, u2;
|
|
Packit |
67cb25 |
double v;
|
|
Packit |
67cb25 |
size_t xi, yi;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* first compute the indices into the data arrays where we are interpolating */
|
|
Packit |
67cb25 |
if (xa != NULL)
|
|
Packit |
67cb25 |
xi = gsl_interp_accel_find(xa, xarr, xsize, x);
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
xi = gsl_interp_bsearch(xarr, x, 0, xsize - 1);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (ya != NULL)
|
|
Packit |
67cb25 |
yi = gsl_interp_accel_find(ya, yarr, ysize, y);
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
yi = gsl_interp_bsearch(yarr, y, 0, ysize - 1);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* find the minimum and maximum values on the grid cell in each dimension */
|
|
Packit |
67cb25 |
xmin = xarr[xi];
|
|
Packit |
67cb25 |
xmax = xarr[xi + 1];
|
|
Packit |
67cb25 |
ymin = yarr[yi];
|
|
Packit |
67cb25 |
ymax = yarr[yi + 1];
|
|
Packit |
67cb25 |
zminmin = zarr[IDX2D(xi, yi, state)];
|
|
Packit |
67cb25 |
zminmax = zarr[IDX2D(xi, yi + 1, state)];
|
|
Packit |
67cb25 |
zmaxmin = zarr[IDX2D(xi + 1, yi, state)];
|
|
Packit |
67cb25 |
zmaxmax = zarr[IDX2D(xi + 1, yi + 1, state)];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* get the width and height of the grid cell */
|
|
Packit |
67cb25 |
dx = xmax - xmin;
|
|
Packit |
67cb25 |
dy = ymax - ymin;
|
|
Packit |
67cb25 |
t = (x - xmin)/dx;
|
|
Packit |
67cb25 |
u = (y - ymin)/dy;
|
|
Packit |
67cb25 |
dt = 1./dx; /* partial t / partial x */
|
|
Packit |
67cb25 |
du = 1./dy; /* partial u / partial y */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
zxminmin = state->zx[IDX2D(xi, yi, state)]/dt;
|
|
Packit |
67cb25 |
zxminmax = state->zx[IDX2D(xi, yi + 1, state)]/dt;
|
|
Packit |
67cb25 |
zxmaxmin = state->zx[IDX2D(xi + 1, yi, state)]/dt;
|
|
Packit |
67cb25 |
zxmaxmax = state->zx[IDX2D(xi + 1, yi + 1, state)]/dt;
|
|
Packit |
67cb25 |
zyminmin = state->zy[IDX2D(xi, yi, state)]/du;
|
|
Packit |
67cb25 |
zyminmax = state->zy[IDX2D(xi, yi + 1, state)]/du;
|
|
Packit |
67cb25 |
zymaxmin = state->zy[IDX2D(xi + 1, yi, state)]/du;
|
|
Packit |
67cb25 |
zymaxmax = state->zy[IDX2D(xi + 1, yi + 1, state)]/du;
|
|
Packit |
67cb25 |
zxyminmin = state->zxy[IDX2D(xi, yi, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxyminmax = state->zxy[IDX2D(xi, yi + 1, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxymaxmin = state->zxy[IDX2D(xi + 1, yi, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxymaxmax = state->zxy[IDX2D(xi + 1, yi + 1, state)]/(dt*du);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
t0 = 1;
|
|
Packit |
67cb25 |
t1 = t;
|
|
Packit |
67cb25 |
t2 = t*t;
|
|
Packit |
67cb25 |
u0 = 1;
|
|
Packit |
67cb25 |
u1 = u;
|
|
Packit |
67cb25 |
u2 = u*u;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
*z_pp = 0;
|
|
Packit |
67cb25 |
v = zxyminmin;
|
|
Packit |
67cb25 |
*z_pp += v*t0*u0;
|
|
Packit |
67cb25 |
v = -3*zxminmin + 3*zxminmax - 2*zxyminmin - zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 2*v*t0*u1;
|
|
Packit |
67cb25 |
v = 2*zxminmin - 2*zxminmax + zxyminmin + zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 3*v*t0*u2;
|
|
Packit |
67cb25 |
v = -3*zyminmin + 3*zymaxmin - 2*zxyminmin - zxymaxmin;
|
|
Packit |
67cb25 |
*z_pp += 2*v*t1*u0;
|
|
Packit |
67cb25 |
v = 9*zminmin - 9*zmaxmin + 9*zmaxmax - 9*zminmax + 6*zxminmin + 3*zxmaxmin - 3*zxmaxmax - 6*zxminmax + 6*zyminmin - 6*zymaxmin - 3*zymaxmax + 3*zyminmax + 4*zxyminmin + 2*zxymaxmin + zxymaxmax + 2*zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 4*v*t1*u1;
|
|
Packit |
67cb25 |
v = -6*zminmin + 6*zmaxmin - 6*zmaxmax + 6*zminmax - 4*zxminmin - 2*zxmaxmin + 2*zxmaxmax + 4*zxminmax - 3*zyminmin + 3*zymaxmin + 3*zymaxmax - 3*zyminmax - 2*zxyminmin - zxymaxmin - zxymaxmax - 2*zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 6*v*t1*u2;
|
|
Packit |
67cb25 |
v = 2*zyminmin - 2*zymaxmin + zxyminmin + zxymaxmin;
|
|
Packit |
67cb25 |
*z_pp += 3*v*t2*u0;
|
|
Packit |
67cb25 |
v = -6*zminmin + 6*zmaxmin - 6*zmaxmax + 6*zminmax - 3*zxminmin - 3*zxmaxmin + 3*zxmaxmax + 3*zxminmax - 4*zyminmin + 4*zymaxmin + 2*zymaxmax - 2*zyminmax - 2*zxyminmin - 2*zxymaxmin - zxymaxmax - zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 6*v*t2*u1;
|
|
Packit |
67cb25 |
v = 4*zminmin - 4*zmaxmin + 4*zmaxmax - 4*zminmax + 2*zxminmin + 2*zxmaxmin - 2*zxmaxmax - 2*zxminmax + 2*zyminmin - 2*zymaxmin - 2*zymaxmax + 2*zyminmax + zxyminmin + zxymaxmin + zxymaxmax + zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 9*v*t2*u2;
|
|
Packit |
67cb25 |
*z_pp *= dt*du;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
bicubic_deriv_yy(const void * vstate, const double xarr[], const double yarr[],
|
|
Packit |
67cb25 |
const double zarr[], size_t xsize, size_t ysize,
|
|
Packit |
67cb25 |
double x, double y,
|
|
Packit |
67cb25 |
gsl_interp_accel * xa, gsl_interp_accel * ya, double * z_pp)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
bicubic_state_t *state = (bicubic_state_t *) vstate;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double xmin, xmax, ymin, ymax;
|
|
Packit |
67cb25 |
double zminmin, zminmax, zmaxmin, zmaxmax;
|
|
Packit |
67cb25 |
double zxminmin, zxminmax, zxmaxmin, zxmaxmax;
|
|
Packit |
67cb25 |
double zyminmin, zyminmax, zymaxmin, zymaxmax;
|
|
Packit |
67cb25 |
double zxyminmin, zxyminmax, zxymaxmin, zxymaxmax;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
double dx, dy; /* size of the grid cell */
|
|
Packit |
67cb25 |
double dt, du;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/*
|
|
Packit |
67cb25 |
* t and u are the positions within the grid cell at which we are computing
|
|
Packit |
67cb25 |
* the interpolation, in units of grid cell size
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
double t, u;
|
|
Packit |
67cb25 |
double t0, t1, t2, t3, u0, u1;
|
|
Packit |
67cb25 |
double v;
|
|
Packit |
67cb25 |
size_t xi, yi;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* first compute the indices into the data arrays where we are interpolating */
|
|
Packit |
67cb25 |
if (xa != NULL)
|
|
Packit |
67cb25 |
xi = gsl_interp_accel_find(xa, xarr, xsize, x);
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
xi = gsl_interp_bsearch(xarr, x, 0, xsize - 1);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (ya != NULL)
|
|
Packit |
67cb25 |
yi = gsl_interp_accel_find(ya, yarr, ysize, y);
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
yi = gsl_interp_bsearch(yarr, y, 0, ysize - 1);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* find the minimum and maximum values on the grid cell in each dimension */
|
|
Packit |
67cb25 |
xmin = xarr[xi];
|
|
Packit |
67cb25 |
xmax = xarr[xi + 1];
|
|
Packit |
67cb25 |
ymin = yarr[yi];
|
|
Packit |
67cb25 |
ymax = yarr[yi + 1];
|
|
Packit |
67cb25 |
zminmin = zarr[IDX2D(xi, yi, state)];
|
|
Packit |
67cb25 |
zminmax = zarr[IDX2D(xi, yi + 1, state)];
|
|
Packit |
67cb25 |
zmaxmin = zarr[IDX2D(xi + 1, yi, state)];
|
|
Packit |
67cb25 |
zmaxmax = zarr[IDX2D(xi + 1, yi + 1, state)];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* get the width and height of the grid cell */
|
|
Packit |
67cb25 |
dx = xmax - xmin;
|
|
Packit |
67cb25 |
dy = ymax - ymin;
|
|
Packit |
67cb25 |
t = (x - xmin)/dx;
|
|
Packit |
67cb25 |
u = (y - ymin)/dy;
|
|
Packit |
67cb25 |
dt = 1./dx; /* partial t / partial x */
|
|
Packit |
67cb25 |
du = 1./dy; /* partial u / partial y */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
zxminmin = state->zx[IDX2D(xi, yi, state)]/dt;
|
|
Packit |
67cb25 |
zxminmax = state->zx[IDX2D(xi, yi + 1, state)]/dt;
|
|
Packit |
67cb25 |
zxmaxmin = state->zx[IDX2D(xi + 1, yi, state)]/dt;
|
|
Packit |
67cb25 |
zxmaxmax = state->zx[IDX2D(xi + 1, yi + 1, state)]/dt;
|
|
Packit |
67cb25 |
zyminmin = state->zy[IDX2D(xi, yi, state)]/du;
|
|
Packit |
67cb25 |
zyminmax = state->zy[IDX2D(xi, yi + 1, state)]/du;
|
|
Packit |
67cb25 |
zymaxmin = state->zy[IDX2D(xi + 1, yi, state)]/du;
|
|
Packit |
67cb25 |
zymaxmax = state->zy[IDX2D(xi + 1, yi + 1, state)]/du;
|
|
Packit |
67cb25 |
zxyminmin = state->zxy[IDX2D(xi, yi, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxyminmax = state->zxy[IDX2D(xi, yi + 1, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxymaxmin = state->zxy[IDX2D(xi + 1, yi, state)]/(dt*du);
|
|
Packit |
67cb25 |
zxymaxmax = state->zxy[IDX2D(xi + 1, yi + 1, state)]/(dt*du);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
t0 = 1;
|
|
Packit |
67cb25 |
t1 = t;
|
|
Packit |
67cb25 |
t2 = t*t;
|
|
Packit |
67cb25 |
t3 = t*t2;
|
|
Packit |
67cb25 |
u0 = 1;
|
|
Packit |
67cb25 |
u1 = u;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
*z_pp = 0;
|
|
Packit |
67cb25 |
v = -3*zminmin + 3*zminmax - 2*zyminmin - zyminmax;
|
|
Packit |
67cb25 |
*z_pp += 2*v*t0*u0;
|
|
Packit |
67cb25 |
v = 2*zminmin-2*zminmax + zyminmin + zyminmax;
|
|
Packit |
67cb25 |
*z_pp += 6*v*t0*u1;
|
|
Packit |
67cb25 |
v = -3*zxminmin + 3*zxminmax - 2*zxyminmin - zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 2*v*t1*u0;
|
|
Packit |
67cb25 |
v = 2*zxminmin - 2*zxminmax + zxyminmin + zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 6*v*t1*u1;
|
|
Packit |
67cb25 |
v = 9*zminmin - 9*zmaxmin + 9*zmaxmax - 9*zminmax + 6*zxminmin + 3*zxmaxmin - 3*zxmaxmax - 6*zxminmax + 6*zyminmin - 6*zymaxmin - 3*zymaxmax + 3*zyminmax + 4*zxyminmin + 2*zxymaxmin + zxymaxmax + 2*zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 2*v*t2*u0;
|
|
Packit |
67cb25 |
v = -6*zminmin + 6*zmaxmin - 6*zmaxmax + 6*zminmax - 4*zxminmin - 2*zxmaxmin + 2*zxmaxmax + 4*zxminmax - 3*zyminmin + 3*zymaxmin + 3*zymaxmax - 3*zyminmax - 2*zxyminmin - zxymaxmin - zxymaxmax - 2*zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 6*v*t2*u1;
|
|
Packit |
67cb25 |
v = -6*zminmin + 6*zmaxmin - 6*zmaxmax + 6*zminmax - 3*zxminmin - 3*zxmaxmin + 3*zxmaxmax + 3*zxminmax - 4*zyminmin + 4*zymaxmin + 2*zymaxmax - 2*zyminmax - 2*zxyminmin - 2*zxymaxmin - zxymaxmax - zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 2*v*t3*u0;
|
|
Packit |
67cb25 |
v = 4*zminmin - 4*zmaxmin + 4*zmaxmax - 4*zminmax + 2*zxminmin + 2*zxmaxmin - 2*zxmaxmax - 2*zxminmax + 2*zyminmin - 2*zymaxmin - 2*zymaxmax + 2*zyminmax + zxyminmin + zxymaxmin + zxymaxmax + zxyminmax;
|
|
Packit |
67cb25 |
*z_pp += 6*v*t3*u1;
|
|
Packit |
67cb25 |
*z_pp *= du*du;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static const gsl_interp2d_type bicubic_type = {
|
|
Packit |
67cb25 |
"bicubic",
|
|
Packit |
67cb25 |
4,
|
|
Packit |
67cb25 |
&bicubic_alloc,
|
|
Packit |
67cb25 |
&bicubic_init,
|
|
Packit |
67cb25 |
&bicubic_eval,
|
|
Packit |
67cb25 |
&bicubic_deriv_x,
|
|
Packit |
67cb25 |
&bicubic_deriv_y,
|
|
Packit |
67cb25 |
&bicubic_deriv_xx,
|
|
Packit |
67cb25 |
&bicubic_deriv_xy,
|
|
Packit |
67cb25 |
&bicubic_deriv_yy,
|
|
Packit |
67cb25 |
&bicubic_free
|
|
Packit |
67cb25 |
};
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
const gsl_interp2d_type * gsl_interp2d_bicubic = &bicubic_type;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
#undef IDX2D
|