|
Packit |
67cb25 |
/* filter/test_gaussian.c
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* Copyright (C) 2018 Patrick Alken
|
|
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 <gsl/gsl_math.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_vector.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_filter.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_movstat.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_test.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_rng.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_randist.h>
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* compute Gaussian filter by explicitely constructing window and computing weighted sum */
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
slow_gaussian(const gsl_filter_end_t etype, const double alpha, const size_t order, const gsl_vector * x,
|
|
Packit |
67cb25 |
gsl_vector * y, const size_t K)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const size_t n = x->size;
|
|
Packit |
67cb25 |
const size_t H = K / 2;
|
|
Packit |
67cb25 |
double *window = malloc(K * sizeof(double));
|
|
Packit |
67cb25 |
double *kernel = malloc(K * sizeof(double));
|
|
Packit |
67cb25 |
gsl_vector_view k = gsl_vector_view_array(kernel, K);
|
|
Packit |
67cb25 |
size_t i;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_filter_gaussian_kernel(alpha, order, 1, &k.vector);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 0; i < n; ++i)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t wsize = gsl_movstat_fill(etype, x, i, H, H, window);
|
|
Packit |
67cb25 |
double sum = 0.0;
|
|
Packit |
67cb25 |
size_t j;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (j = 0; j < wsize; ++j)
|
|
Packit |
67cb25 |
sum += window[j] * kernel[wsize - j - 1];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector_set(y, i, sum);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
free(window);
|
|
Packit |
67cb25 |
free(kernel);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void
|
|
Packit |
67cb25 |
fdiff(const gsl_vector * x, gsl_vector * dx)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const size_t N = x->size;
|
|
Packit |
67cb25 |
size_t i;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 1; i < N - 1; ++i)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
double xm1 = gsl_vector_get(x, i - 1);
|
|
Packit |
67cb25 |
double xp1 = gsl_vector_get(x, i + 1);
|
|
Packit |
67cb25 |
gsl_vector_set(dx, i, 0.5 * (xp1 - xm1));
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector_set(dx, 0, gsl_vector_get(x, 1) - gsl_vector_get(x, 0));
|
|
Packit |
67cb25 |
gsl_vector_set(dx, N - 1, gsl_vector_get(x, N - 1) - gsl_vector_get(x, N - 2));
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void
|
|
Packit |
67cb25 |
test_gaussian_kernel(const double alpha, const size_t K)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const size_t max_order = 3;
|
|
Packit |
67cb25 |
gsl_vector * kernel = gsl_vector_alloc(K);
|
|
Packit |
67cb25 |
gsl_vector * deriv = gsl_vector_alloc(K);
|
|
Packit |
67cb25 |
gsl_vector * deriv_fd = gsl_vector_alloc(K);
|
|
Packit |
67cb25 |
char buf[2048];
|
|
Packit |
67cb25 |
size_t order;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_filter_gaussian_kernel(alpha, 0, 0, kernel);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (order = 1; order <= max_order; ++order)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_filter_gaussian_kernel(alpha, order, 0, deriv);
|
|
Packit |
67cb25 |
fdiff(kernel, deriv_fd);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
sprintf(buf, "gaussian kernel order=%zu alpha=%g K=%zu", order, alpha, K);
|
|
Packit |
67cb25 |
compare_vectors(1.0e-2, deriv_fd, deriv, buf);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector_memcpy(kernel, deriv);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector_free(kernel);
|
|
Packit |
67cb25 |
gsl_vector_free(deriv);
|
|
Packit |
67cb25 |
gsl_vector_free(deriv_fd);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void
|
|
Packit |
67cb25 |
test_gaussian_proc(const double tol, const double alpha, const size_t order, const size_t n, const size_t K,
|
|
Packit |
67cb25 |
const gsl_filter_end_t etype, gsl_rng * rng_p)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_filter_gaussian_workspace * w = gsl_filter_gaussian_alloc(K);
|
|
Packit |
67cb25 |
gsl_vector * x = gsl_vector_alloc(n);
|
|
Packit |
67cb25 |
gsl_vector * y = gsl_vector_alloc(n);
|
|
Packit |
67cb25 |
gsl_vector * z = gsl_vector_alloc(n);
|
|
Packit |
67cb25 |
char buf[2048];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
random_vector(x, rng_p);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* y = filter(x) with slow brute force method */
|
|
Packit |
67cb25 |
slow_gaussian(etype, alpha, order, x, y, K);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* y = filter(x) with fast method */
|
|
Packit |
67cb25 |
gsl_filter_gaussian(etype, alpha, order, x, z, w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* test y = z */
|
|
Packit |
67cb25 |
sprintf(buf, "n=%zu K=%zu endtype=%u alpha=%g order=%zu gaussian random", n, K, etype, alpha, order);
|
|
Packit |
67cb25 |
compare_vectors(tol, z, y, buf);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* z = filter(x) in-place */
|
|
Packit |
67cb25 |
gsl_vector_memcpy(z, x);
|
|
Packit |
67cb25 |
gsl_filter_gaussian(etype, alpha, order, z, z, w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
sprintf(buf, "n=%zu K=%zu endtype=%u alpha=%g order=%zu gaussian random in-place", n, K, etype, alpha, order);
|
|
Packit |
67cb25 |
compare_vectors(tol, z, y, buf);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_filter_gaussian_free(w);
|
|
Packit |
67cb25 |
gsl_vector_free(x);
|
|
Packit |
67cb25 |
gsl_vector_free(y);
|
|
Packit |
67cb25 |
gsl_vector_free(z);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void
|
|
Packit |
67cb25 |
test_gaussian_deriv(const double alpha, const size_t n, const size_t K)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
#if 0
|
|
Packit |
67cb25 |
const double f_low = 1.0;
|
|
Packit |
67cb25 |
const double f_high = 50.0;
|
|
Packit |
67cb25 |
const double gamma = 2.0 * M_PI / (n - 1.0);
|
|
Packit |
67cb25 |
const double dt = 1.0 / (n - 1.0);
|
|
Packit |
67cb25 |
gsl_vector *x = gsl_vector_alloc(n);
|
|
Packit |
67cb25 |
gsl_vector *dx = gsl_vector_alloc(n);
|
|
Packit |
67cb25 |
gsl_vector *y1 = gsl_vector_alloc(n);
|
|
Packit |
67cb25 |
gsl_vector *y2 = gsl_vector_alloc(n);
|
|
Packit |
67cb25 |
gsl_filter_gaussian_workspace *w = gsl_filter_gaussian_alloc(K);
|
|
Packit |
67cb25 |
size_t i;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* make input signal composed of two sine waves at different frequencies */
|
|
Packit |
67cb25 |
for (i = 0; i < n; ++i)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
double xi = sin(gamma * f_low * i) + sin(gamma * f_high * i);
|
|
Packit |
67cb25 |
double dxi = gamma * f_low * cos(gamma * f_low * i) +
|
|
Packit |
67cb25 |
gamma * f_high * cos(gamma * f_high * i);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector_set(x, i, xi);
|
|
Packit |
67cb25 |
gsl_vector_set(dx, i, dxi);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* compute y1 = G * dx(t)/dt */
|
|
Packit |
67cb25 |
gsl_filter_gaussian(alpha, 0, dx, y1, w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* compute y2 = dG/dt * x(t) */
|
|
Packit |
67cb25 |
gsl_filter_gaussian(alpha, 1, x, y2, w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 0; i < n; ++i)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
printf("%zu %.12e %.12e %.12e %.12e\n",
|
|
Packit |
67cb25 |
i,
|
|
Packit |
67cb25 |
gsl_vector_get(x, i),
|
|
Packit |
67cb25 |
gsl_vector_get(dx, i),
|
|
Packit |
67cb25 |
gsl_vector_get(y1, i),
|
|
Packit |
67cb25 |
gsl_vector_get(y2, i));
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector_free(x);
|
|
Packit |
67cb25 |
gsl_vector_free(dx);
|
|
Packit |
67cb25 |
gsl_vector_free(y1);
|
|
Packit |
67cb25 |
gsl_vector_free(y2);
|
|
Packit |
67cb25 |
gsl_filter_gaussian_free(w);
|
|
Packit |
67cb25 |
#endif
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void
|
|
Packit |
67cb25 |
test_gaussian(gsl_rng * r)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const double tol = 1.0e-10;
|
|
Packit |
67cb25 |
size_t order;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
test_gaussian_kernel(3.0, 2001);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (order = 0; order <= 3; ++order)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
test_gaussian_proc(tol, 2.5, order, 1000, 21, GSL_FILTER_END_PADZERO, r);
|
|
Packit |
67cb25 |
test_gaussian_proc(tol, 3.0, order, 500, 11, GSL_FILTER_END_PADZERO, r);
|
|
Packit |
67cb25 |
test_gaussian_proc(tol, 1.0, order, 50, 101, GSL_FILTER_END_PADZERO, r);
|
|
Packit |
67cb25 |
test_gaussian_proc(tol, 2.0, order, 50, 11, GSL_FILTER_END_PADZERO, r);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
test_gaussian_proc(tol, 2.5, order, 1000, 21, GSL_FILTER_END_PADVALUE, r);
|
|
Packit |
67cb25 |
test_gaussian_proc(tol, 3.0, order, 500, 11, GSL_FILTER_END_PADVALUE, r);
|
|
Packit |
67cb25 |
test_gaussian_proc(tol, 1.0, order, 50, 101, GSL_FILTER_END_PADVALUE, r);
|
|
Packit |
67cb25 |
test_gaussian_proc(tol, 2.0, order, 50, 11, GSL_FILTER_END_PADVALUE, r);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
test_gaussian_proc(tol, 2.5, order, 1000, 21, GSL_FILTER_END_TRUNCATE, r);
|
|
Packit |
67cb25 |
test_gaussian_proc(tol, 3.0, order, 500, 11, GSL_FILTER_END_TRUNCATE, r);
|
|
Packit |
67cb25 |
test_gaussian_proc(tol, 1.0, order, 50, 101, GSL_FILTER_END_TRUNCATE, r);
|
|
Packit |
67cb25 |
test_gaussian_proc(tol, 2.0, order, 50, 11, GSL_FILTER_END_TRUNCATE, r);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|