|
Packit |
67cb25 |
/* filter/test_median.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_test.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_rng.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_filter.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_movstat.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_statistics.h>
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* compute filtered data by explicitely constructing window, sorting it and finding median */
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
slow_movmedian(const gsl_filter_end_t etype, const gsl_vector * x, gsl_vector * y, const int K)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const size_t H = K / 2;
|
|
Packit |
67cb25 |
const size_t n = x->size;
|
|
Packit |
67cb25 |
double *window = malloc(K * sizeof(double));
|
|
Packit |
67cb25 |
size_t i;
|
|
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 yi = gsl_stats_median(window, 1, wsize);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector_set(y, i, yi);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
free(window);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void
|
|
Packit |
67cb25 |
test_median_proc(const double tol, 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_median_workspace *w = gsl_filter_median_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 |
/* test moving median with random input */
|
|
Packit |
67cb25 |
random_vector(x, rng_p);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* y = median(x) with slow brute force algorithm */
|
|
Packit |
67cb25 |
slow_movmedian(etype, x, y, K);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* z = median(x) */
|
|
Packit |
67cb25 |
gsl_filter_median(etype, x, z, w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* test y = z */
|
|
Packit |
67cb25 |
sprintf(buf, "n=%zu K=%zu endtype=%u median random", n, K, etype);
|
|
Packit |
67cb25 |
compare_vectors(tol, z, y, buf);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* z = median(x) in-place */
|
|
Packit |
67cb25 |
gsl_vector_memcpy(z, x);
|
|
Packit |
67cb25 |
gsl_filter_median(etype, z, z, w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
sprintf(buf, "n=%zu K=%zu endtype=%u median random in-place", n, K, etype);
|
|
Packit |
67cb25 |
compare_vectors(tol, z, y, buf);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_vector_free(x);
|
|
Packit |
67cb25 |
gsl_vector_free(y);
|
|
Packit |
67cb25 |
gsl_vector_free(z);
|
|
Packit |
67cb25 |
gsl_filter_median_free(w);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void
|
|
Packit |
67cb25 |
test_median(gsl_rng * rng_p)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 1000, 1, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 1000, 3, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 1000, 5, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 100, 301, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 5000, 17, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 5000, 9, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 50, 901, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 50, 201, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 1000, 1, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 1000, 3, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 1000, 5, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 100, 301, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 5000, 17, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 5000, 9, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 50, 901, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 50, 201, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 1000, 1, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 1000, 3, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 1000, 5, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 100, 301, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 5000, 17, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 5000, 9, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 50, 901, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_median_proc(GSL_DBL_EPSILON, 50, 201, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
}
|