|
Packit |
67cb25 |
/* movstat/test_Qn.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_statistics.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_sort.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_test.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_movstat.h>
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* calculate Q_n statistic for input vector using slow/naive algorithm */
|
|
Packit |
67cb25 |
static int
|
|
Packit |
67cb25 |
slow_movQn(const gsl_movstat_end_t etype, const gsl_vector * x, gsl_vector * y,
|
|
Packit |
67cb25 |
const int H, const int J)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const size_t n = x->size;
|
|
Packit |
67cb25 |
const int K = H + J + 1;
|
|
Packit |
67cb25 |
double *window = malloc(K * sizeof(double));
|
|
Packit |
67cb25 |
double *work = malloc(3 * K * sizeof(double));
|
|
Packit |
67cb25 |
int *work_int = malloc(5 * K * sizeof(int));
|
|
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, J, window);
|
|
Packit |
67cb25 |
double Qn;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_sort(window, 1, wsize);
|
|
Packit |
67cb25 |
Qn = gsl_stats_Qn_from_sorted_data(window, 1, wsize, work, work_int);
|
|
Packit |
67cb25 |
gsl_vector_set(y, i, Qn);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
free(window);
|
|
Packit |
67cb25 |
free(work);
|
|
Packit |
67cb25 |
free(work_int);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static double
|
|
Packit |
67cb25 |
func_Qn(const size_t n, double x[], void * params)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
double *work = malloc(3 * n * sizeof(double));
|
|
Packit |
67cb25 |
int *work_int = malloc(5 * n * sizeof(int));
|
|
Packit |
67cb25 |
double Qn;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
(void) params;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gsl_sort(x, 1, n);
|
|
Packit |
67cb25 |
Qn = gsl_stats_Qn_from_sorted_data(x, 1, n, work, work_int);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
free(work);
|
|
Packit |
67cb25 |
free(work_int);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return Qn;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void
|
|
Packit |
67cb25 |
test_Qn_proc(const double tol, const size_t n, const size_t H, const size_t J,
|
|
Packit |
67cb25 |
const gsl_movstat_end_t etype, gsl_rng *rng_p)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
gsl_movstat_workspace *w;
|
|
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 |
gsl_movstat_function F;
|
|
Packit |
67cb25 |
char buf[2048];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
F.function = func_Qn;
|
|
Packit |
67cb25 |
F.params = NULL;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (H == J)
|
|
Packit |
67cb25 |
w = gsl_movstat_alloc(2*H + 1);
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
w = gsl_movstat_alloc2(H, J);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* test moving median with random input */
|
|
Packit |
67cb25 |
random_vector(x, rng_p);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* y = Q_n(x) with slow brute force algorithm */
|
|
Packit |
67cb25 |
slow_movQn(etype, x, y, H, J);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* z = Q_n(x) */
|
|
Packit |
67cb25 |
gsl_movstat_Qn(etype, x, z, w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* test y = z */
|
|
Packit |
67cb25 |
sprintf(buf, "n=%zu H=%zu J=%zu endtype=%u Qn random", n, H, J, etype);
|
|
Packit |
67cb25 |
compare_vectors(tol, z, y, buf);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* z = Q_n(x) in-place */
|
|
Packit |
67cb25 |
gsl_vector_memcpy(z, x);
|
|
Packit |
67cb25 |
gsl_movstat_Qn(etype, z, z, w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
sprintf(buf, "n=%zu H=%zu J=%zu endtype=%u Qn random in-place", n, H, J, etype);
|
|
Packit |
67cb25 |
compare_vectors(tol, z, y, buf);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* z = Q_n(x) with user-defined function */
|
|
Packit |
67cb25 |
gsl_movstat_apply(etype, &F, x, z, w);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
sprintf(buf, "n=%zu H=%zu J=%zu endtype=%u Qn user", n, H, J, 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_movstat_free(w);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
static void
|
|
Packit |
67cb25 |
test_Qn(gsl_rng * rng_p)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 1000, 0, 0, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 1000, 5, 5, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 1000, 5, 2, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 1000, 2, 5, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 2000, 50, 0, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 2000, 0, 50, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 20, 50, 50, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 20, 1, 50, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 20, 50, 1, GSL_MOVSTAT_END_PADZERO, rng_p);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 1000, 0, 0, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 1000, 5, 5, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 1000, 5, 2, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 1000, 2, 5, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 2000, 50, 0, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 2000, 0, 50, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 20, 50, 50, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 20, 1, 50, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 20, 50, 1, GSL_MOVSTAT_END_PADVALUE, rng_p);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 1000, 0, 0, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 1000, 5, 5, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 1000, 5, 2, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 1000, 2, 5, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 2000, 50, 0, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 2000, 0, 50, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 20, 50, 50, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 20, 1, 50, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
test_Qn_proc(GSL_DBL_EPSILON, 20, 50, 1, GSL_MOVSTAT_END_TRUNCATE, rng_p);
|
|
Packit |
67cb25 |
}
|