-
Notifications
You must be signed in to change notification settings - Fork 214
/
Copy pathtest1d_resample.c
99 lines (81 loc) · 2.55 KB
/
test1d_resample.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
/* histogram/test1d_resample.c
*
* Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Brian Gough
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 3 of the License, or (at
* your option) any later version.
*
* This program is distributed in the hope that it will be useful, but
* WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, write to the Free Software
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
*/
#include <config.h>
#include <stdlib.h>
#include <math.h>
#include <gsl/gsl_histogram.h>
#include <gsl/gsl_test.h>
#include <gsl/gsl_ieee_utils.h>
#include "urand.c"
void
test1d_resample (void)
{
size_t i;
int status = 0;
gsl_histogram *h;
gsl_ieee_env_setup ();
h = gsl_histogram_calloc_uniform (10, 0.0, 1.0);
gsl_histogram_increment (h, 0.1);
gsl_histogram_increment (h, 0.2);
gsl_histogram_increment (h, 0.2);
gsl_histogram_increment (h, 0.3);
{
gsl_histogram_pdf *p = gsl_histogram_pdf_alloc (10);
gsl_histogram *hh = gsl_histogram_calloc_uniform (100, 0.0, 1.0);
gsl_histogram_pdf_init (p, h);
for (i = 0; i < 100000; i++)
{
double u = urand();
double x = gsl_histogram_pdf_sample (p, u);
gsl_histogram_increment (hh, x);
}
for (i = 0; i < 100; i++)
{
double y = gsl_histogram_get (hh, i) / 2500;
double x, xmax;
size_t k;
double ya;
gsl_histogram_get_range (hh, i, &x, &xmax);
gsl_histogram_find (h, x, &k);
ya = gsl_histogram_get (h, k);
if (ya == 0)
{
if (y != 0)
{
printf ("%d: %g vs %g\n", (int) i, y, ya);
status = 1;
}
}
else
{
double err = 1 / sqrt (gsl_histogram_get (hh, i));
double sigma = fabs ((y - ya) / (ya * err));
if (sigma > 3)
{
status = 1;
printf ("%g vs %g err=%g sigma=%g\n", y, ya, err, sigma);
}
}
}
gsl_histogram_pdf_free (p) ;
gsl_histogram_free (hh);
gsl_test (status, "gsl_histogram_pdf_sample within statistical errors");
}
gsl_histogram_free (h);
}