gsl_bspline__bspline.c 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498
  1. /* bspline/bspline.c
  2. *
  3. * Copyright (C) 2006 Patrick Alken
  4. *
  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.
  9. *
  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.
  14. *
  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.
  18. */
  19. #include "gsl__config.h"
  20. #include "gsl_errno.h"
  21. #include "gsl_bspline.h"
  22. /*
  23. * This module contains routines related to calculating B-splines.
  24. * The algorithms used are described in
  25. *
  26. * [1] Carl de Boor, "A Practical Guide to Splines", Springer
  27. * Verlag, 1978.
  28. */
  29. static int bspline_eval_all(const double x, gsl_vector *B, size_t *idx,
  30. gsl_bspline_workspace *w);
  31. static inline size_t bspline_find_interval(const double x, int *flag,
  32. gsl_bspline_workspace *w);
  33. /*
  34. gsl_bspline_alloc()
  35. Allocate space for a bspline workspace. The size of the
  36. workspace is O(5k + nbreak)
  37. Inputs: k - spline order (cubic = 4)
  38. nbreak - number of breakpoints
  39. Return: pointer to workspace
  40. */
  41. gsl_bspline_workspace *
  42. gsl_bspline_alloc(const size_t k, const size_t nbreak)
  43. {
  44. if (k == 0)
  45. {
  46. GSL_ERROR_NULL("k must be at least 1", GSL_EINVAL);
  47. }
  48. else if (nbreak < 2)
  49. {
  50. GSL_ERROR_NULL("nbreak must be at least 2", GSL_EINVAL);
  51. }
  52. else
  53. {
  54. gsl_bspline_workspace *w;
  55. w = (gsl_bspline_workspace *)
  56. calloc(1, sizeof(gsl_bspline_workspace));
  57. if (w == 0)
  58. {
  59. GSL_ERROR_NULL("failed to allocate space for workspace", GSL_ENOMEM);
  60. }
  61. w->k = k;
  62. w->km1 = k - 1;
  63. w->nbreak = nbreak;
  64. w->l = nbreak - 1;
  65. w->n = w->l + k - 1;
  66. w->knots = gsl_vector_alloc(w->n + k);
  67. if (w->knots == 0)
  68. {
  69. gsl_bspline_free(w);
  70. GSL_ERROR_NULL("failed to allocate space for knots vector", GSL_ENOMEM);
  71. }
  72. w->deltal = gsl_vector_alloc(k);
  73. w->deltar = gsl_vector_alloc(k);
  74. if (!w->deltal || !w->deltar)
  75. {
  76. gsl_bspline_free(w);
  77. GSL_ERROR_NULL("failed to allocate space for delta vectors", GSL_ENOMEM);
  78. }
  79. w->B = gsl_vector_alloc(k);
  80. if (w->B == 0)
  81. {
  82. gsl_bspline_free(w);
  83. GSL_ERROR_NULL("failed to allocate space for temporary spline vector", GSL_ENOMEM);
  84. }
  85. return (w);
  86. }
  87. } /* gsl_bspline_alloc() */
  88. /* Return number of coefficients */
  89. size_t
  90. gsl_bspline_ncoeffs (gsl_bspline_workspace * w)
  91. {
  92. return w->n;
  93. }
  94. /* Return order */
  95. size_t
  96. gsl_bspline_order (gsl_bspline_workspace * w)
  97. {
  98. return w->k;
  99. }
  100. /* Return number of breakpoints */
  101. size_t
  102. gsl_bspline_nbreak (gsl_bspline_workspace * w)
  103. {
  104. return w->nbreak;
  105. }
  106. double
  107. gsl_bspline_breakpoint (size_t i, gsl_bspline_workspace * w)
  108. {
  109. size_t j = i + w->k - 1;
  110. return gsl_vector_get(w->knots, j);
  111. }
  112. /*
  113. gsl_bspline_free()
  114. Free a bspline workspace
  115. Inputs: w - workspace to free
  116. Return: none
  117. */
  118. void
  119. gsl_bspline_free(gsl_bspline_workspace *w)
  120. {
  121. if (!w)
  122. return;
  123. if (w->knots)
  124. gsl_vector_free(w->knots);
  125. if (w->deltal)
  126. gsl_vector_free(w->deltal);
  127. if (w->deltar)
  128. gsl_vector_free(w->deltar);
  129. if (w->B)
  130. gsl_vector_free(w->B);
  131. free(w);
  132. } /* gsl_bspline_free() */
  133. /*
  134. gsl_bspline_knots()
  135. Compute the knots from the given breakpoints:
  136. knots(1:k) = breakpts(1)
  137. knots(k+1:k+l-1) = breakpts(i), i = 2 .. l
  138. knots(n+1:n+k) = breakpts(l + 1)
  139. where l is the number of polynomial pieces (l = nbreak - 1) and
  140. n = k + l - 1
  141. (using matlab syntax for the arrays)
  142. The repeated knots at the beginning and end of the interval
  143. correspond to the continuity condition there. See pg. 119
  144. of [1].
  145. Inputs: breakpts - breakpoints
  146. w - bspline workspace
  147. Return: success or error
  148. */
  149. int
  150. gsl_bspline_knots(const gsl_vector *breakpts, gsl_bspline_workspace *w)
  151. {
  152. if (breakpts->size != w->nbreak)
  153. {
  154. GSL_ERROR("breakpts vector has wrong size", GSL_EBADLEN);
  155. }
  156. else
  157. {
  158. size_t i; /* looping */
  159. for (i = 0; i < w->k; ++i)
  160. gsl_vector_set(w->knots, i, gsl_vector_get(breakpts, 0));
  161. for (i = 1; i < w->l; ++i)
  162. {
  163. gsl_vector_set(w->knots, w->k - 1 + i,
  164. gsl_vector_get(breakpts, i));
  165. }
  166. for (i = w->n; i < w->n + w->k; ++i)
  167. gsl_vector_set(w->knots, i, gsl_vector_get(breakpts, w->l));
  168. return GSL_SUCCESS;
  169. }
  170. } /* gsl_bspline_knots() */
  171. /*
  172. gsl_bspline_knots_uniform()
  173. Construct uniformly spaced knots on the interval [a,b] using
  174. the previously specified number of breakpoints. 'a' is the position
  175. of the first breakpoint and 'b' is the position of the last
  176. breakpoint.
  177. Inputs: a - left side of interval
  178. b - right side of interval
  179. w - bspline workspace
  180. Return: success or error
  181. Notes: 1) w->knots is modified to contain the uniformly spaced
  182. knots
  183. 2) The knots vector is set up as follows (using octave syntax):
  184. knots(1:k) = a
  185. knots(k+1:k+l-1) = a + i*delta, i = 1 .. l - 1
  186. knots(n+1:n+k) = b
  187. */
  188. int
  189. gsl_bspline_knots_uniform(const double a, const double b,
  190. gsl_bspline_workspace *w)
  191. {
  192. size_t i; /* looping */
  193. double delta; /* interval spacing */
  194. double x;
  195. delta = (b - a) / (double) w->l;
  196. for (i = 0; i < w->k; ++i)
  197. gsl_vector_set(w->knots, i, a);
  198. x = a + delta;
  199. for (i = 0; i < w->l - 1; ++i)
  200. {
  201. gsl_vector_set(w->knots, w->k + i, x);
  202. x += delta;
  203. }
  204. for (i = w->n; i < w->n + w->k; ++i)
  205. gsl_vector_set(w->knots, i, b);
  206. return GSL_SUCCESS;
  207. } /* gsl_bspline_knots_uniform() */
  208. /*
  209. gsl_bspline_eval()
  210. Evaluate the basis functions B_i(x) for all i. This is
  211. basically a wrapper function for bspline_eval_all()
  212. which formats the output in a nice way.
  213. Inputs: x - point for evaluation
  214. B - (output) where to store B_i(x) values
  215. the length of this vector is
  216. n = nbreak + k - 2 = l + k - 1 = w->n
  217. w - bspline workspace
  218. Return: success or error
  219. Notes: The w->knots vector must be initialized prior to calling
  220. this function (see gsl_bspline_knots())
  221. */
  222. int
  223. gsl_bspline_eval(const double x, gsl_vector *B,
  224. gsl_bspline_workspace *w)
  225. {
  226. if (B->size != w->n)
  227. {
  228. GSL_ERROR("B vector length does not match n", GSL_EBADLEN);
  229. }
  230. else
  231. {
  232. size_t i; /* looping */
  233. size_t idx = 0;
  234. size_t start; /* first non-zero spline */
  235. /* find all non-zero B_i(x) values */
  236. bspline_eval_all(x, w->B, &idx, w);
  237. /* store values in appropriate part of given vector */
  238. start = idx - w->k + 1;
  239. for (i = 0; i < start; ++i)
  240. gsl_vector_set(B, i, 0.0);
  241. for (i = start; i <= idx; ++i)
  242. gsl_vector_set(B, i, gsl_vector_get(w->B, i - start));
  243. for (i = idx + 1; i < w->n; ++i)
  244. gsl_vector_set(B, i, 0.0);
  245. return GSL_SUCCESS;
  246. }
  247. } /* gsl_bspline_eval() */
  248. /****************************************
  249. * INTERNAL ROUTINES *
  250. ****************************************/
  251. /*
  252. bspline_eval_all()
  253. Evaluate all non-zero B-splines B_i(x) using algorithm (8)
  254. of chapter X of [1]
  255. The idea is something like this. Given x \in [t_i, t_{i+1}]
  256. with t_i < t_{i+1} and the t_i are knots, the values of the
  257. B-splines not automatically zero fit into a triangular
  258. array as follows:
  259. 0
  260. 0
  261. 0 B_{i-2,3}
  262. B_{i-1,2}
  263. B_{i,1} B_{i-1,3} ...
  264. B_{i,2}
  265. 0 B_{i,3}
  266. 0
  267. 0
  268. where B_{i,k} is the ith B-spline of order k. The boundary 0s
  269. indicate that those B-splines not in the table vanish at x.
  270. To compute the non-zero B-splines of a given order k, we use
  271. Eqs. (4) and (5) of chapter X of [1]:
  272. (4) B_{i,1}(x) = { 1, t_i <= x < t_{i+1}
  273. 0, else }
  274. (5) B_{i,k}(x) = (x - t_i)
  275. ----------------- B_{i,k-1}(x)
  276. (t_{i+k-1} - t_i)
  277. t_{i+k} - x
  278. + ----------------- B_{i+1,k-1}(x)
  279. t_{i+k} - t_{i+1}
  280. So (4) gives us the first column of the table and we can use
  281. the recurrence relation (5) to get the rest of the columns.
  282. Inputs: x - point at which to evaluate splines
  283. B - (output) where to store B-spline values (length k)
  284. idx - (output) B-spline function index of last output
  285. value (B_{idx}(x) is stored in the last slot of 'B')
  286. w - bspline workspace
  287. Return: success or error
  288. Notes: 1) the w->knots vector must be initialized before calling
  289. this function
  290. 2) On output, B contains:
  291. B = [B_{i-k+1,k}, B_{i-k+2,k}, ..., B_{i-1,k}, B_{i,k}]
  292. where 'i' is stored in 'idx' on output
  293. 3) based on PPPACK bsplvb
  294. */
  295. static int
  296. bspline_eval_all(const double x, gsl_vector *B, size_t *idx,
  297. gsl_bspline_workspace *w)
  298. {
  299. if (B->size != w->k)
  300. {
  301. GSL_ERROR("B vector not of length k", GSL_EBADLEN);
  302. }
  303. else
  304. {
  305. size_t i; /* spline index */
  306. size_t j; /* looping */
  307. size_t ii; /* looping */
  308. int flag = 0; /* error flag */
  309. double saved;
  310. double term;
  311. i = bspline_find_interval(x, &flag, w);
  312. if (flag == -1)
  313. {
  314. GSL_ERROR("x outside of knot interval", GSL_EINVAL);
  315. }
  316. else if (flag == 1)
  317. {
  318. if (x <= gsl_vector_get(w->knots, i) + GSL_DBL_EPSILON)
  319. {
  320. --i;
  321. }
  322. else
  323. {
  324. GSL_ERROR("x outside of knot interval", GSL_EINVAL);
  325. }
  326. }
  327. *idx = i;
  328. gsl_vector_set(B, 0, 1.0);
  329. for (j = 0; j < w->k - 1; ++j)
  330. {
  331. gsl_vector_set(w->deltar, j,
  332. gsl_vector_get(w->knots, i + j + 1) - x);
  333. gsl_vector_set(w->deltal, j,
  334. x - gsl_vector_get(w->knots, i - j));
  335. saved = 0.0;
  336. for (ii = 0; ii <= j; ++ii)
  337. {
  338. term = gsl_vector_get(B, ii) /
  339. (gsl_vector_get(w->deltar, ii) +
  340. gsl_vector_get(w->deltal, j - ii));
  341. gsl_vector_set(B, ii,
  342. saved +
  343. gsl_vector_get(w->deltar, ii) * term);
  344. saved = gsl_vector_get(w->deltal, j - ii) * term;
  345. }
  346. gsl_vector_set(B, j + 1, saved);
  347. }
  348. return GSL_SUCCESS;
  349. }
  350. } /* bspline_eval_all() */
  351. /*
  352. bspline_find_interval()
  353. Find knot interval such that
  354. t_i <= x < t_{i + 1}
  355. where the t_i are knot values.
  356. Inputs: x - x value
  357. flag - (output) error flag
  358. w - bspline workspace
  359. Return: i (index in w->knots corresponding to left limit of interval)
  360. Notes: The error conditions are reported as follows:
  361. Condition Return value Flag
  362. --------- ------------ ----
  363. x < t_0 0 -1
  364. t_i <= x < t_{i+1} i 0
  365. t_{n+k-1} <= x l+k-1 +1
  366. */
  367. static inline size_t
  368. bspline_find_interval(const double x, int *flag,
  369. gsl_bspline_workspace *w)
  370. {
  371. size_t i;
  372. if (x < gsl_vector_get(w->knots, 0))
  373. {
  374. *flag = -1;
  375. return 0;
  376. }
  377. /* find i such that t_i <= x < t_{i+1} */
  378. for (i = w->k - 1; i < w->k + w->l - 1; ++i)
  379. {
  380. const double ti = gsl_vector_get(w->knots, i);
  381. const double tip1 = gsl_vector_get(w->knots, i + 1);
  382. if (tip1 < ti)
  383. {
  384. GSL_ERROR("knots vector is not increasing", GSL_EINVAL);
  385. }
  386. if (ti <= x && x < tip1)
  387. break;
  388. }
  389. if (i == w->k + w->l - 1)
  390. *flag = 1;
  391. else
  392. *flag = 0;
  393. return i;
  394. } /* bspline_find_interval() */