bessel_functions.cpp 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2016 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #include "main.h"
  10. #include "../Eigen/SpecialFunctions"
  11. template<typename X, typename Y>
  12. void verify_component_wise(const X& x, const Y& y)
  13. {
  14. for(Index i=0; i<x.size(); ++i)
  15. {
  16. if((numext::isfinite)(y(i))) {
  17. VERIFY_IS_APPROX( x(i), y(i) );
  18. }
  19. else if((numext::isnan)(y(i)))
  20. VERIFY((numext::isnan)(x(i)));
  21. else
  22. VERIFY_IS_EQUAL( x(i), y(i) );
  23. }
  24. }
  25. template<typename ArrayType> void array_bessel_functions()
  26. {
  27. // Test Bessel function i0. Reference results obtained with SciPy.
  28. {
  29. ArrayType x(21);
  30. ArrayType expected(21);
  31. ArrayType res(21);
  32. x << -20.0, -18.0, -16.0, -14.0, -12.0, -10.0, -8.0, -6.0, -4.0, -2.0, 0.0,
  33. 2.0, 4.0, 6.0, 8.0, 10.0, 12.0, 14.0, 16.0, 18.0, 20.0;
  34. expected << 4.35582826e+07, 6.21841242e+06, 8.93446228e+05, 1.29418563e+05,
  35. 1.89489253e+04, 2.81571663e+03, 4.27564116e+02, 6.72344070e+01,
  36. 1.13019220e+01, 2.27958530e+00, 1.00000000e+00, 2.27958530e+00,
  37. 1.13019220e+01, 6.72344070e+01, 4.27564116e+02, 2.81571663e+03,
  38. 1.89489253e+04, 1.29418563e+05, 8.93446228e+05, 6.21841242e+06,
  39. 4.35582826e+07;
  40. CALL_SUBTEST(res = bessel_i0(x);
  41. verify_component_wise(res, expected););
  42. }
  43. // Test Bessel function i0e. Reference results obtained with SciPy.
  44. {
  45. ArrayType x(21);
  46. ArrayType expected(21);
  47. ArrayType res(21);
  48. x << -20.0, -18.0, -16.0, -14.0, -12.0, -10.0, -8.0, -6.0, -4.0, -2.0, 0.0,
  49. 2.0, 4.0, 6.0, 8.0, 10.0, 12.0, 14.0, 16.0, 18.0, 20.0;
  50. expected << 0.0897803118848, 0.0947062952128, 0.100544127361,
  51. 0.107615251671, 0.116426221213, 0.127833337163, 0.143431781857,
  52. 0.16665743264, 0.207001921224, 0.308508322554, 1.0, 0.308508322554,
  53. 0.207001921224, 0.16665743264, 0.143431781857, 0.127833337163,
  54. 0.116426221213, 0.107615251671, 0.100544127361, 0.0947062952128,
  55. 0.0897803118848;
  56. CALL_SUBTEST(res = bessel_i0e(x);
  57. verify_component_wise(res, expected););
  58. }
  59. // Test Bessel function i1. Reference results obtained with SciPy.
  60. {
  61. ArrayType x(21);
  62. ArrayType expected(21);
  63. ArrayType res(21);
  64. x << -20.0, -18.0, -16.0, -14.0, -12.0, -10.0, -8.0, -6.0, -4.0, -2.0, 0.0,
  65. 2.0, 4.0, 6.0, 8.0, 10.0, 12.0, 14.0, 16.0, 18.0, 20.0;
  66. expected << -4.24549734e+07, -6.04313324e+06, -8.65059436e+05, -1.24707259e+05,
  67. -1.81413488e+04, -2.67098830e+03, -3.99873137e+02, -6.13419368e+01,
  68. -9.75946515e+00, -1.59063685e+00, 0.00000000e+00, 1.59063685e+00,
  69. 9.75946515e+00, 6.13419368e+01, 3.99873137e+02, 2.67098830e+03,
  70. 1.81413488e+04, 1.24707259e+05, 8.65059436e+05, 6.04313324e+06,
  71. 4.24549734e+07;
  72. CALL_SUBTEST(res = bessel_i1(x);
  73. verify_component_wise(res, expected););
  74. }
  75. // Test Bessel function i1e. Reference results obtained with SciPy.
  76. {
  77. ArrayType x(21);
  78. ArrayType expected(21);
  79. ArrayType res(21);
  80. x << -20.0, -18.0, -16.0, -14.0, -12.0, -10.0, -8.0, -6.0, -4.0, -2.0, 0.0,
  81. 2.0, 4.0, 6.0, 8.0, 10.0, 12.0, 14.0, 16.0, 18.0, 20.0;
  82. expected << -0.0875062221833, -0.092036796872, -0.0973496147565,
  83. -0.103697667463, -0.11146429929, -0.121262681384, -0.134142493293,
  84. -0.152051459309, -0.178750839502, -0.215269289249, 0.0, 0.215269289249,
  85. 0.178750839502, 0.152051459309, 0.134142493293, 0.121262681384,
  86. 0.11146429929, 0.103697667463, 0.0973496147565, 0.092036796872,
  87. 0.0875062221833;
  88. CALL_SUBTEST(res = bessel_i1e(x);
  89. verify_component_wise(res, expected););
  90. }
  91. // Test Bessel function j0. Reference results obtained with SciPy.
  92. {
  93. ArrayType x(77);
  94. ArrayType expected(77);
  95. ArrayType res(77);
  96. x << -38., -37., -36., -35., -34., -33., -32., -31., -30.,
  97. -29., -28., -27., -26., -25., -24., -23., -22., -21., -20., -19.,
  98. -18., -17., -16., -15., -14., -13., -12., -11., -10., -9., -8.,
  99. -7., -6., -5., -4., -3., -2., -1., 0., 1., 2., 3.,
  100. 4., 5., 6., 7., 8., 9., 10., 11., 12., 13., 14.,
  101. 15., 16., 17., 18., 19., 20., 21., 22., 23., 24., 25.,
  102. 26., 27., 28., 29., 30., 31., 32., 33., 34., 35., 36.,
  103. 37., 38.;
  104. expected << 0.11433274, 0.01086237, -0.10556738,
  105. -0.12684568, -0.03042119, 0.09727067, 0.13807901, 0.05120815,
  106. -0.08636798, -0.14784876, -0.07315701, 0.07274192, 0.15599932,
  107. 0.09626678, -0.05623027, -0.16241278, -0.12065148, 0.03657907,
  108. 0.16702466, 0.14662944, -0.01335581, -0.16985425, -0.17489907,
  109. -0.01422447, 0.17107348, 0.2069261 , 0.04768931, -0.1711903 ,
  110. -0.24593576, -0.09033361, 0.17165081, 0.30007927, 0.15064526,
  111. -0.17759677, -0.39714981, -0.26005195, 0.22389078, 0.76519769,
  112. 1. , 0.76519769, 0.22389078, -0.26005195, -0.39714981,
  113. -0.17759677, 0.15064526, 0.30007927, 0.17165081, -0.09033361,
  114. -0.24593576, -0.1711903 , 0.04768931, 0.2069261 , 0.17107348,
  115. -0.01422447, -0.17489907, -0.16985425, -0.01335581, 0.14662944,
  116. 0.16702466, 0.03657907, -0.12065148, -0.16241278, -0.05623027,
  117. 0.09626678, 0.15599932, 0.07274192, -0.07315701, -0.14784876,
  118. -0.08636798, 0.05120815, 0.13807901, 0.09727067, -0.03042119,
  119. -0.12684568, -0.10556738, 0.01086237, 0.11433274;
  120. CALL_SUBTEST(res = bessel_j0(x);
  121. verify_component_wise(res, expected););
  122. }
  123. // Test Bessel function j1. Reference results obtained with SciPy.
  124. {
  125. ArrayType x(81);
  126. ArrayType expected(81);
  127. ArrayType res(81);
  128. x << -40., -39., -38., -37., -36., -35., -34., -33., -32., -31., -30.,
  129. -29., -28., -27., -26., -25., -24., -23., -22., -21., -20., -19.,
  130. -18., -17., -16., -15., -14., -13., -12., -11., -10., -9., -8.,
  131. -7., -6., -5., -4., -3., -2., -1., 0., 1., 2., 3.,
  132. 4., 5., 6., 7., 8., 9., 10., 11., 12., 13., 14.,
  133. 15., 16., 17., 18., 19., 20., 21., 22., 23., 24., 25.,
  134. 26., 27., 28., 29., 30., 31., 32., 33., 34., 35., 36.,
  135. 37., 38., 39., 40.;
  136. expected << -0.12603832, -0.0640561 , 0.05916189, 0.13058004, 0.08232981,
  137. -0.04399094, -0.13297118, -0.10061965, 0.02658903, 0.13302432,
  138. 0.11875106, -0.0069342 , -0.13055149, -0.13658472, -0.01504573,
  139. 0.12535025, 0.15403807, 0.03951932, -0.11717779, -0.17112027,
  140. -0.06683312, 0.10570143, 0.18799489, 0.09766849, -0.09039718,
  141. -0.20510404, -0.13337515, 0.07031805, 0.2234471 , 0.1767853 ,
  142. -0.04347275, -0.24531179, -0.23463635, 0.00468282, 0.27668386,
  143. 0.32757914, 0.06604333, -0.33905896, -0.57672481, -0.44005059,
  144. 0. , 0.44005059, 0.57672481, 0.33905896, -0.06604333,
  145. -0.32757914, -0.27668386, -0.00468282, 0.23463635, 0.24531179,
  146. 0.04347275, -0.1767853 , -0.2234471 , -0.07031805, 0.13337515,
  147. 0.20510404, 0.09039718, -0.09766849, -0.18799489, -0.10570143,
  148. 0.06683312, 0.17112027, 0.11717779, -0.03951932, -0.15403807,
  149. -0.12535025, 0.01504573, 0.13658472, 0.13055149, 0.0069342 ,
  150. -0.11875106, -0.13302432, -0.02658903, 0.10061965, 0.13297118,
  151. 0.04399094, -0.08232981, -0.13058004, -0.05916189, 0.0640561 ,
  152. 0.12603832;
  153. CALL_SUBTEST(res = bessel_j1(x);
  154. verify_component_wise(res, expected););
  155. }
  156. // Test Bessel function k0e. Reference results obtained with SciPy.
  157. {
  158. ArrayType x(42);
  159. ArrayType expected(42);
  160. ArrayType res(42);
  161. x << 0.25, 0.5, 1., 2., 3., 4., 5., 6., 7., 8., 9., 10., 11., 12.,
  162. 13., 14., 15., 16., 17., 18., 19., 20., 21., 22., 23., 24., 25.,
  163. 26., 27., 28., 29., 30., 31., 32., 33., 34., 35., 36., 37., 38.,
  164. 39., 40.;
  165. expected << 1.97933385, 1.52410939, 1.14446308, 0.84156822,
  166. 0.6977616 , 0.60929767, 0.54780756, 0.50186313, 0.4658451 ,
  167. 0.43662302, 0.41229555, 0.39163193, 0.3737955 , 0.35819488,
  168. 0.34439865, 0.33208364, 0.32100235, 0.31096159, 0.30180802,
  169. 0.29341821, 0.28569149, 0.27854488, 0.2719092 , 0.26572635,
  170. 0.25994703, 0.25452917, 0.2494366 , 0.24463801, 0.24010616,
  171. 0.23581722, 0.23175022, 0.22788667, 0.22421014, 0.22070602,
  172. 0.21736123, 0.21416406, 0.21110397, 0.20817141, 0.20535778,
  173. 0.20265524, 0.20005668, 0.19755558;
  174. CALL_SUBTEST(res = bessel_k0e(x);
  175. verify_component_wise(res, expected););
  176. }
  177. // Test Bessel function k0. Reference results obtained with SciPy.
  178. {
  179. ArrayType x(42);
  180. ArrayType expected(42);
  181. ArrayType res(42);
  182. x << 0.25, 0.5, 1., 2., 3., 4., 5., 6., 7., 8., 9., 10., 11., 12.,
  183. 13., 14., 15., 16., 17., 18., 19., 20., 21., 22., 23., 24., 25.,
  184. 26., 27., 28., 29., 30., 31., 32., 33., 34., 35., 36., 37., 38.,
  185. 39., 40.;
  186. expected << 1.54150675, 0.92441907, 4.21024438e-01, 1.13893873e-01,
  187. 3.47395044e-02, 1.11596761e-02, 3.69109833e-03, 1.24399433e-03,
  188. 4.24795742e-04, 1.46470705e-04, 5.08813130e-05, 1.77800623e-05,
  189. 6.24302055e-06, 2.20082540e-06, 7.78454386e-07, 2.76137082e-07,
  190. 9.81953648e-08, 3.49941166e-08, 1.24946640e-08, 4.46875334e-09,
  191. 1.60067129e-09, 5.74123782e-10, 2.06176797e-10, 7.41235161e-11,
  192. 2.66754511e-11, 9.60881878e-12, 3.46416156e-12, 1.24987740e-12,
  193. 4.51286453e-13, 1.63053459e-13, 5.89495073e-14, 2.13247750e-14,
  194. 7.71838266e-15, 2.79505752e-15, 1.01266123e-15, 3.67057597e-16,
  195. 1.33103515e-16, 4.82858338e-17, 1.75232770e-17, 6.36161716e-18,
  196. 2.31029936e-18, 8.39286110e-19;
  197. CALL_SUBTEST(res = bessel_k0(x);
  198. verify_component_wise(res, expected););
  199. }
  200. // Test Bessel function k0e. Reference results obtained with SciPy.
  201. {
  202. ArrayType x(42);
  203. ArrayType expected(42);
  204. ArrayType res(42);
  205. x << 0.25, 0.5, 1., 2., 3., 4., 5., 6., 7., 8., 9., 10., 11., 12.,
  206. 13., 14., 15., 16., 17., 18., 19., 20., 21., 22., 23., 24., 25.,
  207. 26., 27., 28., 29., 30., 31., 32., 33., 34., 35., 36., 37., 38.,
  208. 39., 40.;
  209. expected << 1.97933385, 1.52410939, 1.14446308, 0.84156822,
  210. 0.6977616 , 0.60929767, 0.54780756, 0.50186313,
  211. 0.4658451 , 0.43662302, 0.41229555, 0.39163193,
  212. 0.3737955 , 0.35819488, 0.34439865, 0.33208364,
  213. 0.32100235, 0.31096159, 0.30180802, 0.29341821,
  214. 0.28569149, 0.27854488, 0.2719092 , 0.26572635,
  215. 0.25994703, 0.25452917, 0.2494366 , 0.24463801,
  216. 0.24010616, 0.23581722, 0.23175022, 0.22788667,
  217. 0.22421014, 0.22070602, 0.21736123, 0.21416406,
  218. 0.21110397, 0.20817141, 0.20535778, 0.20265524,
  219. 0.20005668, 0.19755558;
  220. CALL_SUBTEST(res = bessel_k0e(x);
  221. verify_component_wise(res, expected););
  222. }
  223. // Test Bessel function k1. Reference results obtained with SciPy.
  224. {
  225. ArrayType x(42);
  226. ArrayType expected(42);
  227. ArrayType res(42);
  228. x << 0.25, 0.5, 1., 2., 3., 4., 5., 6., 7., 8., 9., 10., 11., 12.,
  229. 13., 14., 15., 16., 17., 18., 19., 20., 21., 22., 23., 24., 25.,
  230. 26., 27., 28., 29., 30., 31., 32., 33., 34., 35., 36., 37., 38.,
  231. 39., 40.;
  232. expected << 3.74702597, 1.65644112, 6.01907230e-01, 1.39865882e-01,
  233. 4.01564311e-02, 1.24834989e-02, 4.04461345e-03, 1.34391972e-03,
  234. 4.54182487e-04, 1.55369212e-04, 5.36370164e-05, 1.86487735e-05,
  235. 6.52086067e-06, 2.29075746e-06, 8.07858841e-07, 2.85834365e-07,
  236. 1.01417294e-07, 3.60715712e-08, 1.28570417e-08, 4.59124963e-09,
  237. 1.64226697e-09, 5.88305797e-10, 2.11029922e-10, 7.57898116e-11,
  238. 2.72493059e-11, 9.80699893e-12, 3.53277807e-12, 1.27369078e-12,
  239. 4.59568940e-13, 1.65940011e-13, 5.99574032e-14, 2.16773200e-14,
  240. 7.84189960e-15, 2.83839927e-15, 1.02789171e-15, 3.72416929e-16,
  241. 1.34991783e-16, 4.89519373e-17, 1.77585196e-17, 6.44478588e-18,
  242. 2.33973340e-18, 8.49713195e-19;
  243. CALL_SUBTEST(res = bessel_k1(x);
  244. verify_component_wise(res, expected););
  245. }
  246. // Test Bessel function k1e. Reference results obtained with SciPy.
  247. {
  248. ArrayType x(42);
  249. ArrayType expected(42);
  250. ArrayType res(42);
  251. x << 0.25, 0.5, 1., 2., 3., 4., 5., 6., 7., 8., 9., 10., 11., 12.,
  252. 13., 14., 15., 16., 17., 18., 19., 20., 21., 22., 23., 24., 25.,
  253. 26., 27., 28., 29., 30., 31., 32., 33., 34., 35., 36., 37., 38.,
  254. 39., 40.;
  255. expected << 4.81127659, 2.73100971, 1.63615349, 1.03347685,
  256. 0.80656348, 0.68157595, 0.60027386, 0.54217591,
  257. 0.49807158, 0.46314909, 0.43462525, 0.41076657,
  258. 0.39043094, 0.37283175, 0.35740757, 0.34374563,
  259. 0.33153489, 0.32053597, 0.31056123, 0.30146131,
  260. 0.29311559, 0.2854255 , 0.27830958, 0.27169987,
  261. 0.26553913, 0.25977879, 0.25437733, 0.249299 ,
  262. 0.24451285, 0.23999191, 0.2357126 , 0.23165413,
  263. 0.22779816, 0.22412841, 0.22063036, 0.21729103,
  264. 0.21409878, 0.21104314, 0.20811462, 0.20530466,
  265. 0.20260547, 0.20000997;
  266. CALL_SUBTEST(res = bessel_k1e(x);
  267. verify_component_wise(res, expected););
  268. }
  269. // Test Bessel function y0. Reference results obtained with SciPy.
  270. {
  271. ArrayType x(42);
  272. ArrayType expected(42);
  273. ArrayType res(42);
  274. x << 0.25, 0.5, 1., 2., 3., 4., 5., 6., 7., 8., 9., 10., 11., 12.,
  275. 13., 14., 15., 16., 17., 18., 19., 20., 21., 22., 23., 24., 25.,
  276. 26., 27., 28., 29., 30., 31., 32., 33., 34., 35., 36., 37., 38.,
  277. 39., 40.;
  278. expected << -0.93157302, -0.44451873, 0.08825696, 0.51037567, 0.37685001,
  279. -0.01694074, -0.30851763, -0.28819468, -0.02594974, 0.22352149,
  280. 0.2499367 , 0.05567117, -0.16884732, -0.22523731, -0.07820786,
  281. 0.12719257, 0.2054643 , 0.095811 , -0.0926372 , -0.18755216,
  282. -0.10951969, 0.0626406 , 0.17020176, 0.1198876 , -0.03598179,
  283. -0.15283403, -0.12724943, 0.01204463, 0.13521498, 0.13183647,
  284. 0.00948116, -0.11729573, -0.13383266, -0.02874248, 0.09913483,
  285. 0.13340405, 0.04579799, -0.08085609, -0.13071488, -0.06066076,
  286. 0.06262353, 0.12593642;
  287. CALL_SUBTEST(res = bessel_y0(x);
  288. verify_component_wise(res, expected););
  289. }
  290. // Test Bessel function y1. Reference results obtained with SciPy.
  291. {
  292. ArrayType x(42);
  293. ArrayType expected(42);
  294. ArrayType res(42);
  295. x << 0.25, 0.5, 1., 2., 3., 4., 5., 6., 7., 8., 9., 10., 11., 12.,
  296. 13., 14., 15., 16., 17., 18., 19., 20., 21., 22., 23., 24., 25.,
  297. 26., 27., 28., 29., 30., 31., 32., 33., 34., 35., 36., 37., 38.,
  298. 39., 40.;
  299. expected << -2.70410523, -1.47147239, -0.78121282, -0.10703243,
  300. 0.32467442, 0.39792571, 0.14786314, -0.17501034, -0.30266724,
  301. -0.15806046, 0.10431458, 0.24901542, 0.16370554, -0.05709922,
  302. -0.21008141, -0.16664484, 0.02107363, 0.17797517, 0.16720504,
  303. 0.00815513, -0.14956011, -0.16551161, -0.03253926, 0.12340586,
  304. 0.1616692 , 0.05305978, -0.09882996, -0.15579655, -0.07025124,
  305. 0.07552213, 0.14803412, 0.08442557, -0.05337283, -0.13854483,
  306. -0.09578012, 0.03238588, 0.12751273, 0.10445477, -0.01262946,
  307. -0.11514066, -0.11056411, -0.00579351;
  308. CALL_SUBTEST(res = bessel_y1(x);
  309. verify_component_wise(res, expected););
  310. }
  311. }
  312. EIGEN_DECLARE_TEST(bessel_functions)
  313. {
  314. CALL_SUBTEST_1(array_bessel_functions<ArrayXf>());
  315. CALL_SUBTEST_2(array_bessel_functions<ArrayXd>());
  316. }