gsl_wavelet__dwt.c 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391
  1. /* wavelet/dwt.c
  2. *
  3. * Copyright (C) 2004 Ivo Alxneit
  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. /* function dwt_step is based on the public domain function pwt.c
  20. * available from http://www.numerical-recipes.com
  21. */
  22. #include "gsl__config.h"
  23. #include "gsl_errno.h"
  24. #include "gsl_wavelet.h"
  25. #include "gsl_wavelet2d.h"
  26. #define ELEMENT(a,stride,i) ((a)[(stride)*(i)])
  27. static int binary_logn (const size_t n);
  28. static void dwt_step (const gsl_wavelet * w, double *a, size_t stride, size_t n, gsl_wavelet_direction dir, gsl_wavelet_workspace * work);
  29. static int
  30. binary_logn (const size_t n)
  31. {
  32. size_t ntest;
  33. size_t logn = 0;
  34. size_t k = 1;
  35. while (k < n)
  36. {
  37. k *= 2;
  38. logn++;
  39. }
  40. ntest = (1 << logn);
  41. if (n != ntest)
  42. {
  43. return -1; /* n is not a power of 2 */
  44. }
  45. return logn;
  46. }
  47. static void
  48. dwt_step (const gsl_wavelet * w, double *a, size_t stride, size_t n,
  49. gsl_wavelet_direction dir, gsl_wavelet_workspace * work)
  50. {
  51. double ai, ai1;
  52. size_t i, ii;
  53. size_t jf;
  54. size_t k;
  55. size_t n1, ni, nh, nmod;
  56. for (i = 0; i < work->n; i++)
  57. {
  58. work->scratch[i] = 0.0;
  59. }
  60. nmod = w->nc * n;
  61. nmod -= w->offset; /* center support */
  62. n1 = n - 1;
  63. nh = n >> 1;
  64. if (dir == gsl_wavelet_forward)
  65. {
  66. for (ii = 0, i = 0; i < n; i += 2, ii++)
  67. {
  68. ni = i + nmod;
  69. for (k = 0; k < w->nc; k++)
  70. {
  71. jf = n1 & (ni + k);
  72. work->scratch[ii] += w->h1[k] * ELEMENT (a, stride, jf);
  73. work->scratch[ii + nh] += w->g1[k] * ELEMENT (a, stride, jf);
  74. }
  75. }
  76. }
  77. else
  78. {
  79. for (ii = 0, i = 0; i < n; i += 2, ii++)
  80. {
  81. ai = ELEMENT (a, stride, ii);
  82. ai1 = ELEMENT (a, stride, ii + nh);
  83. ni = i + nmod;
  84. for (k = 0; k < w->nc; k++)
  85. {
  86. jf = (n1 & (ni + k));
  87. work->scratch[jf] += (w->h2[k] * ai + w->g2[k] * ai1);
  88. }
  89. }
  90. }
  91. for (i = 0; i < n; i++)
  92. {
  93. ELEMENT (a, stride, i) = work->scratch[i];
  94. }
  95. }
  96. int
  97. gsl_wavelet_transform (const gsl_wavelet * w,
  98. double *data, size_t stride, size_t n,
  99. gsl_wavelet_direction dir,
  100. gsl_wavelet_workspace * work)
  101. {
  102. size_t i;
  103. if (work->n < n)
  104. {
  105. GSL_ERROR ("not enough workspace provided", GSL_EINVAL);
  106. }
  107. if (binary_logn (n) == -1)
  108. {
  109. GSL_ERROR ("n is not a power of 2", GSL_EINVAL);
  110. }
  111. if (n < 2)
  112. {
  113. return GSL_SUCCESS;
  114. }
  115. if (dir == gsl_wavelet_forward)
  116. {
  117. for (i = n; i >= 2; i >>= 1)
  118. {
  119. dwt_step (w, data, stride, i, dir, work);
  120. }
  121. }
  122. else
  123. {
  124. for (i = 2; i <= n; i <<= 1)
  125. {
  126. dwt_step (w, data, stride, i, dir, work);
  127. }
  128. }
  129. return GSL_SUCCESS;
  130. }
  131. int
  132. gsl_wavelet_transform_forward (const gsl_wavelet * w,
  133. double *data, size_t stride, size_t n,
  134. gsl_wavelet_workspace * work)
  135. {
  136. return gsl_wavelet_transform (w, data, stride, n, gsl_wavelet_forward, work);
  137. }
  138. int
  139. gsl_wavelet_transform_inverse (const gsl_wavelet * w,
  140. double *data, size_t stride, size_t n,
  141. gsl_wavelet_workspace * work)
  142. {
  143. return gsl_wavelet_transform (w, data, stride, n, gsl_wavelet_backward, work);
  144. }
  145. /* Leaving this out for now BJG */
  146. #if 0
  147. int
  148. gsl_dwt_vector (const gsl_wavelet * w, gsl_vector *v, gsl_wavelet_direction
  149. dir, gsl_wavelet_workspace * work)
  150. {
  151. return gsl_dwt (w, v->data, v->stride, v->size, dir, work);
  152. }
  153. #endif
  154. int
  155. gsl_wavelet2d_transform (const gsl_wavelet * w,
  156. double *data, size_t tda, size_t size1,
  157. size_t size2, gsl_wavelet_direction dir,
  158. gsl_wavelet_workspace * work)
  159. {
  160. size_t i;
  161. if (size1 != size2)
  162. {
  163. GSL_ERROR ("2d dwt works only with square matrix", GSL_EINVAL);
  164. }
  165. if (work->n < size1)
  166. {
  167. GSL_ERROR ("not enough workspace provided", GSL_EINVAL);
  168. }
  169. if (binary_logn (size1) == -1)
  170. {
  171. GSL_ERROR ("n is not a power of 2", GSL_EINVAL);
  172. }
  173. if (size1 < 2)
  174. {
  175. return GSL_SUCCESS;
  176. }
  177. if (dir == gsl_wavelet_forward)
  178. {
  179. for (i = 0; i < size1; i++) /* for every row j */
  180. {
  181. gsl_wavelet_transform (w, &ELEMENT(data, tda, i), 1, size1, dir, work);
  182. }
  183. for (i = 0; i < size2; i++) /* for every column j */
  184. {
  185. gsl_wavelet_transform (w, &ELEMENT(data, 1, i), tda, size2, dir, work);
  186. }
  187. }
  188. else
  189. {
  190. for (i = 0; i < size2; i++) /* for every column j */
  191. {
  192. gsl_wavelet_transform (w, &ELEMENT(data, 1, i), tda, size2, dir, work);
  193. }
  194. for (i = 0; i < size1; i++) /* for every row j */
  195. {
  196. gsl_wavelet_transform (w, &ELEMENT(data, tda, i), 1, size1, dir, work);
  197. }
  198. }
  199. return GSL_SUCCESS;
  200. }
  201. int
  202. gsl_wavelet2d_nstransform (const gsl_wavelet * w,
  203. double *data, size_t tda, size_t size1,
  204. size_t size2, gsl_wavelet_direction dir,
  205. gsl_wavelet_workspace * work)
  206. {
  207. size_t i, j;
  208. if (size1 != size2)
  209. {
  210. GSL_ERROR ("2d dwt works only with square matrix", GSL_EINVAL);
  211. }
  212. if (work->n < size1)
  213. {
  214. GSL_ERROR ("not enough workspace provided", GSL_EINVAL);
  215. }
  216. if (binary_logn (size1) == -1)
  217. {
  218. GSL_ERROR ("n is not a power of 2", GSL_EINVAL);
  219. }
  220. if (size1 < 2)
  221. {
  222. return GSL_SUCCESS;
  223. }
  224. if (dir == gsl_wavelet_forward)
  225. {
  226. for (i = size1; i >= 2; i >>= 1)
  227. {
  228. for (j = 0; j < i; j++) /* for every row j */
  229. {
  230. dwt_step (w, &ELEMENT(data, tda, j), 1, i, dir, work);
  231. }
  232. for (j = 0; j < i; j++) /* for every column j */
  233. {
  234. dwt_step (w, &ELEMENT(data, 1, j), tda, i, dir, work);
  235. }
  236. }
  237. }
  238. else
  239. {
  240. for (i = 2; i <= size1; i <<= 1)
  241. {
  242. for (j = 0; j < i; j++) /* for every column j */
  243. {
  244. dwt_step (w, &ELEMENT(data, 1, j), tda, i, dir, work);
  245. }
  246. for (j = 0; j < i; j++) /* for every row j */
  247. {
  248. dwt_step (w, &ELEMENT(data, tda, j), 1, i, dir, work);
  249. }
  250. }
  251. }
  252. return GSL_SUCCESS;
  253. }
  254. int
  255. gsl_wavelet2d_transform_forward (const gsl_wavelet * w,
  256. double *data, size_t tda, size_t size1,
  257. size_t size2, gsl_wavelet_workspace * work)
  258. {
  259. return gsl_wavelet2d_transform (w, data, tda, size1, size2, gsl_wavelet_forward, work);
  260. }
  261. int
  262. gsl_wavelet2d_transform_inverse (const gsl_wavelet * w,
  263. double *data, size_t tda, size_t size1,
  264. size_t size2, gsl_wavelet_workspace * work)
  265. {
  266. return gsl_wavelet2d_transform (w, data, tda, size1, size2, gsl_wavelet_backward, work);
  267. }
  268. int
  269. gsl_wavelet2d_nstransform_forward (const gsl_wavelet * w,
  270. double *data, size_t tda, size_t size1,
  271. size_t size2, gsl_wavelet_workspace * work)
  272. {
  273. return gsl_wavelet2d_nstransform (w, data, tda, size1, size2, gsl_wavelet_forward, work);
  274. }
  275. int
  276. gsl_wavelet2d_nstransform_inverse (const gsl_wavelet * w,
  277. double *data, size_t tda, size_t size1,
  278. size_t size2, gsl_wavelet_workspace * work)
  279. {
  280. return gsl_wavelet2d_nstransform (w, data, tda, size1, size2, gsl_wavelet_backward, work);
  281. }
  282. int
  283. gsl_wavelet2d_transform_matrix (const gsl_wavelet * w,
  284. gsl_matrix * a,
  285. gsl_wavelet_direction dir,
  286. gsl_wavelet_workspace * work)
  287. {
  288. return gsl_wavelet2d_transform (w, a->data,
  289. a->tda, a->size1, a->size2,
  290. dir, work);
  291. }
  292. int
  293. gsl_wavelet2d_transform_matrix_forward (const gsl_wavelet * w,
  294. gsl_matrix * a,
  295. gsl_wavelet_workspace * work)
  296. {
  297. return gsl_wavelet2d_transform (w, a->data,
  298. a->tda, a->size1, a->size2,
  299. gsl_wavelet_forward, work);
  300. }
  301. int
  302. gsl_wavelet2d_transform_matrix_inverse (const gsl_wavelet * w,
  303. gsl_matrix * a,
  304. gsl_wavelet_workspace * work)
  305. {
  306. return gsl_wavelet2d_transform (w, a->data,
  307. a->tda, a->size1, a->size2,
  308. gsl_wavelet_backward, work);
  309. }
  310. int
  311. gsl_wavelet2d_nstransform_matrix (const gsl_wavelet * w,
  312. gsl_matrix * a,
  313. gsl_wavelet_direction dir,
  314. gsl_wavelet_workspace * work)
  315. {
  316. return gsl_wavelet2d_nstransform (w, a->data,
  317. a->tda, a->size1, a->size2,
  318. dir, work);
  319. }
  320. int
  321. gsl_wavelet2d_nstransform_matrix_forward (const gsl_wavelet * w,
  322. gsl_matrix * a,
  323. gsl_wavelet_workspace * work)
  324. {
  325. return gsl_wavelet2d_nstransform (w, a->data,
  326. a->tda, a->size1, a->size2,
  327. gsl_wavelet_forward, work);
  328. }
  329. int
  330. gsl_wavelet2d_nstransform_matrix_inverse (const gsl_wavelet * w,
  331. gsl_matrix * a,
  332. gsl_wavelet_workspace * work)
  333. {
  334. return gsl_wavelet2d_nstransform (w, a->data,
  335. a->tda, a->size1, a->size2,
  336. gsl_wavelet_backward, work);
  337. }