Blame movstat/test_qqr.c

Packit 67cb25
/* movstat/test_qqr.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_movstat.h>
Packit 67cb25
Packit 67cb25
/* compute moving QQR by explicitely constructing window and computing QQR */
Packit 67cb25
int
Packit 67cb25
slow_movqqr(const gsl_movstat_end_t etype, const double q, const gsl_vector * x,
Packit 67cb25
            gsl_vector * y, 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
  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 quant1, quant2;
Packit 67cb25
Packit 67cb25
      gsl_sort(window, 1, wsize);
Packit 67cb25
      quant1 = gsl_stats_quantile_from_sorted_data(window, 1, wsize, q);
Packit 67cb25
      quant2 = gsl_stats_quantile_from_sorted_data(window, 1, wsize, 1.0 - q);
Packit 67cb25
Packit 67cb25
      gsl_vector_set(y, i, quant2 - quant1);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  free(window);
Packit 67cb25
Packit 67cb25
  return GSL_SUCCESS;
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static double
Packit 67cb25
func_qqr(const size_t n, double x[], void * params)
Packit 67cb25
{
Packit 67cb25
  double q = *(double *) params;
Packit 67cb25
  double quant1, quant2;
Packit 67cb25
Packit 67cb25
  (void) params;
Packit 67cb25
Packit 67cb25
  gsl_sort(x, 1, n);
Packit 67cb25
  quant1 = gsl_stats_quantile_from_sorted_data(x, 1, n, q);
Packit 67cb25
  quant2 = gsl_stats_quantile_from_sorted_data(x, 1, n, 1.0 - q);
Packit 67cb25
Packit 67cb25
  return (quant2 - quant1);
Packit 67cb25
}
Packit 67cb25
Packit 67cb25
static void
Packit 67cb25
test_qqr_proc(const double tol, const double q, 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 = gsl_movstat_alloc2(H, J);
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_qqr;
Packit 67cb25
  F.params = (void *) &q;
Packit 67cb25
Packit 67cb25
  random_vector(x, rng_p);
Packit 67cb25
Packit 67cb25
  /* y = QQR(x) with slow brute force method */
Packit 67cb25
  slow_movqqr(etype, q, x, y, H, J);
Packit 67cb25
Packit 67cb25
  /* y = QQR(x) with fast method */
Packit 67cb25
  gsl_movstat_qqr(etype, x, q, z, w);
Packit 67cb25
Packit 67cb25
  /* test y = z */
Packit 67cb25
  sprintf(buf, "n=%zu H=%zu J=%zu endtype=%u QQR random", n, H, J, etype);
Packit 67cb25
  compare_vectors(tol, z, y, buf);
Packit 67cb25
Packit 67cb25
  /* z = QQR(x) in-place */
Packit 67cb25
  gsl_vector_memcpy(z, x);
Packit 67cb25
  gsl_movstat_qqr(etype, z, q, z, w);
Packit 67cb25
Packit 67cb25
  sprintf(buf, "n=%zu H=%zu J=%zu endtype=%u QQR random in-place", n, H, J, etype);
Packit 67cb25
  compare_vectors(tol, z, y, buf);
Packit 67cb25
Packit 67cb25
  /* z = QQR(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 QQR user", n, H, J, etype);
Packit 67cb25
  compare_vectors(tol, z, y, buf);
Packit 67cb25
Packit 67cb25
  gsl_movstat_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_qqr(gsl_rng * rng_p)
Packit 67cb25
{
Packit 67cb25
  const double eps = 1.0e-10;
Packit 67cb25
Packit 67cb25
  test_qqr_proc(eps, 0.1, 100, 0, 0, GSL_MOVSTAT_END_PADZERO, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.1, 1000, 3, 3, GSL_MOVSTAT_END_PADZERO, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.2, 1000, 0, 5, GSL_MOVSTAT_END_PADZERO, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.3, 1000, 5, 0, GSL_MOVSTAT_END_PADZERO, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.25, 2000, 10, 5, GSL_MOVSTAT_END_PADZERO, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.4, 2000, 5, 10, GSL_MOVSTAT_END_PADZERO, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.25, 20, 50, 50, GSL_MOVSTAT_END_PADZERO, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.25, 20, 10, 50, GSL_MOVSTAT_END_PADZERO, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.25, 20, 50, 10, GSL_MOVSTAT_END_PADZERO, rng_p);
Packit 67cb25
Packit 67cb25
  test_qqr_proc(eps, 0.1, 100, 0, 0, GSL_MOVSTAT_END_PADVALUE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.1, 1000, 3, 3, GSL_MOVSTAT_END_PADVALUE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.2, 1000, 0, 5, GSL_MOVSTAT_END_PADVALUE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.3, 1000, 5, 0, GSL_MOVSTAT_END_PADVALUE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.25, 2000, 10, 5, GSL_MOVSTAT_END_PADVALUE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.4, 2000, 5, 10, GSL_MOVSTAT_END_PADVALUE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.25, 20, 50, 50, GSL_MOVSTAT_END_PADVALUE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.25, 20, 10, 50, GSL_MOVSTAT_END_PADVALUE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.25, 20, 50, 10, GSL_MOVSTAT_END_PADVALUE, rng_p);
Packit 67cb25
Packit 67cb25
  test_qqr_proc(eps, 0.1, 100, 0, 0, GSL_MOVSTAT_END_TRUNCATE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.1, 1000, 3, 3, GSL_MOVSTAT_END_TRUNCATE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.2, 1000, 0, 5, GSL_MOVSTAT_END_TRUNCATE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.3, 1000, 5, 0, GSL_MOVSTAT_END_TRUNCATE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.25, 2000, 10, 5, GSL_MOVSTAT_END_TRUNCATE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.4, 2000, 5, 10, GSL_MOVSTAT_END_TRUNCATE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.25, 20, 50, 50, GSL_MOVSTAT_END_TRUNCATE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.25, 20, 10, 50, GSL_MOVSTAT_END_TRUNCATE, rng_p);
Packit 67cb25
  test_qqr_proc(eps, 0.25, 20, 50, 10, GSL_MOVSTAT_END_TRUNCATE, rng_p);
Packit 67cb25
}