3 * Copyright (C) 2006, 2007 Brian Gough
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 3 of the License, or (at
8 * your option) any later version.
10 * This program is distributed in the hope that it will be useful, but
11 * WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 * General Public License for more details.
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
21 #include <gsl/gsl_test.h>
22 #include <gsl/gsl_errno.h>
23 #include <gsl/gsl_bspline.h>
24 #include <gsl/gsl_ieee_utils.h>
27 test_bspline (gsl_bspline_workspace * bw)
32 size_t ncoeffs = gsl_bspline_ncoeffs (bw);
33 size_t order = gsl_bspline_order (bw);
34 size_t nbreak = gsl_bspline_nbreak (bw);
35 double a = gsl_bspline_breakpoint (0, bw);
36 double b = gsl_bspline_breakpoint (nbreak - 1, bw);
38 B = gsl_vector_alloc (ncoeffs);
40 for (i = 0; i < n; i++)
42 double xi = a + (b - a) * (i / (n - 1.0));
44 gsl_bspline_eval (xi, B, bw);
46 for (j = 0; j < ncoeffs; j++)
48 double Bj = gsl_vector_get (B, j);
49 int s = (Bj < 0 || Bj > 1);
51 "basis-spline coefficient %u is in range [0,1] for x=%g",
56 gsl_test_rel (sum, 1.0, order * GSL_DBL_EPSILON,
57 "basis-spline order %u is normalized for x=%g", order,
67 main (int argc, char **argv)
69 size_t order, breakpoints, i;
71 gsl_ieee_env_setup ();
73 argc = 0; /* prevent warnings about unused parameters */
76 for (order = 1; order < 10; order++)
78 for (breakpoints = 2; breakpoints < 100; breakpoints++)
80 double a = -1.23 * order, b = 45.6 * order;
81 gsl_bspline_workspace *bw = gsl_bspline_alloc (order, breakpoints);
82 gsl_bspline_knots_uniform (a, b, bw);
84 gsl_bspline_free (bw);
89 for (order = 1; order < 10; order++)
91 for (breakpoints = 2; breakpoints < 100; breakpoints++)
93 double a = -1.23 * order, b = 45.6 * order;
94 gsl_bspline_workspace *bw = gsl_bspline_alloc (order, breakpoints);
95 gsl_vector *k = gsl_vector_alloc (breakpoints);
96 for (i = 0; i < breakpoints; i++)
99 f = sqrt (i / (breakpoints - 1.0));
100 x = (1 - f) * a + f * b;
101 gsl_vector_set (k, i, x);
103 gsl_bspline_knots (k, bw);
106 gsl_bspline_free (bw);
110 exit (gsl_test_summary ());