bessel_data.cpp 9.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355
  1. // Copyright (c) 2007 John Maddock
  2. // Use, modification and distribution are subject to the
  3. // Boost Software License, Version 1.0. (See accompanying file
  4. // LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
  5. //
  6. // Computes test data for the various bessel functions using
  7. // archived - deliberately naive - version of the code.
  8. // We'll rely on the high precision of mp_t to get us out of
  9. // trouble and not worry about how long the calculations take.
  10. // This provides a reasonably independent set of test data to
  11. // compare against newly added asymptotic expansions etc.
  12. //
  13. #include <fstream>
  14. #include <boost/math/tools/test_data.hpp>
  15. #include <boost/math/special_functions/bessel.hpp>
  16. #include "mp_t.hpp"
  17. using namespace boost::math::tools;
  18. using namespace boost::math;
  19. using namespace boost::math::detail;
  20. using namespace std;
  21. // Compute J(v, x) and Y(v, x) simultaneously by Steed's method, see
  22. // Barnett et al, Computer Physics Communications, vol 8, 377 (1974)
  23. template <typename T>
  24. int bessel_jy_bare(T v, T x, T* J, T* Y, int kind = need_j|need_y)
  25. {
  26. // Jv1 = J_(v+1), Yv1 = Y_(v+1), fv = J_(v+1) / J_v
  27. // Ju1 = J_(u+1), Yu1 = Y_(u+1), fu = J_(u+1) / J_u
  28. T u, Jv, Ju, Yv, Yv1, Yu, Yu1, fv, fu;
  29. T W, p, q, gamma, current, prev, next;
  30. bool reflect = false;
  31. int n, k, s;
  32. using namespace std;
  33. using namespace boost::math::tools;
  34. using namespace boost::math::constants;
  35. if (v < 0)
  36. {
  37. reflect = true;
  38. v = -v; // v is non-negative from here
  39. kind = need_j|need_y; // need both for reflection formula
  40. }
  41. n = real_cast<int>(v + 0.5L);
  42. u = v - n; // -1/2 <= u < 1/2
  43. if (x < 0)
  44. {
  45. *J = *Y = policies::raise_domain_error<T>("",
  46. "Real argument x=%1% must be non-negative, complex number result not supported", x, policies::policy<>());
  47. return 1;
  48. }
  49. if (x == 0)
  50. {
  51. *J = *Y = policies::raise_overflow_error<T>(
  52. "", 0, policies::policy<>());
  53. return 1;
  54. }
  55. // x is positive until reflection
  56. W = T(2) / (x * pi<T>()); // Wronskian
  57. if (x <= 2) // x in (0, 2]
  58. {
  59. if(temme_jy(u, x, &Yu, &Yu1, policies::policy<>())) // Temme series
  60. {
  61. // domain error:
  62. *J = *Y = Yu;
  63. return 1;
  64. }
  65. prev = Yu;
  66. current = Yu1;
  67. for (k = 1; k <= n; k++) // forward recurrence for Y
  68. {
  69. next = 2 * (u + k) * current / x - prev;
  70. prev = current;
  71. current = next;
  72. }
  73. Yv = prev;
  74. Yv1 = current;
  75. CF1_jy(v, x, &fv, &s, policies::policy<>()); // continued fraction CF1
  76. Jv = W / (Yv * fv - Yv1); // Wronskian relation
  77. }
  78. else // x in (2, \infty)
  79. {
  80. // Get Y(u, x):
  81. CF1_jy(v, x, &fv, &s, policies::policy<>());
  82. // tiny initial value to prevent overflow
  83. T init = sqrt(tools::min_value<T>());
  84. prev = fv * s * init;
  85. current = s * init;
  86. for (k = n; k > 0; k--) // backward recurrence for J
  87. {
  88. next = 2 * (u + k) * current / x - prev;
  89. prev = current;
  90. current = next;
  91. }
  92. T ratio = (s * init) / current; // scaling ratio
  93. // can also call CF1() to get fu, not much difference in precision
  94. fu = prev / current;
  95. CF2_jy(u, x, &p, &q, policies::policy<>()); // continued fraction CF2
  96. T t = u / x - fu; // t = J'/J
  97. gamma = (p - t) / q;
  98. Ju = sign(current) * sqrt(W / (q + gamma * (p - t)));
  99. Jv = Ju * ratio; // normalization
  100. Yu = gamma * Ju;
  101. Yu1 = Yu * (u/x - p - q/gamma);
  102. // compute Y:
  103. prev = Yu;
  104. current = Yu1;
  105. for (k = 1; k <= n; k++) // forward recurrence for Y
  106. {
  107. next = 2 * (u + k) * current / x - prev;
  108. prev = current;
  109. current = next;
  110. }
  111. Yv = prev;
  112. }
  113. if (reflect)
  114. {
  115. T z = (u + n % 2) * pi<T>();
  116. *J = cos(z) * Jv - sin(z) * Yv; // reflection formula
  117. *Y = sin(z) * Jv + cos(z) * Yv;
  118. }
  119. else
  120. {
  121. *J = Jv;
  122. *Y = Yv;
  123. }
  124. return 0;
  125. }
  126. int progress = 0;
  127. template <class T>
  128. T cyl_bessel_j_bare(T v, T x)
  129. {
  130. T j, y;
  131. bessel_jy_bare(v, x, &j, &y);
  132. std::cout << progress++ << ": J(" << v << ", " << x << ") = " << j << std::endl;
  133. if(fabs(j) > 1e30)
  134. throw std::domain_error("");
  135. return j;
  136. }
  137. template <class T>
  138. T cyl_bessel_i_bare(T v, T x)
  139. {
  140. using namespace std;
  141. if(x < 0)
  142. {
  143. // better have integer v:
  144. if(floor(v) == v)
  145. {
  146. T r = cyl_bessel_i_bare(v, -x);
  147. if(tools::real_cast<int>(v) & 1)
  148. r = -r;
  149. return r;
  150. }
  151. else
  152. return policies::raise_domain_error<T>(
  153. "",
  154. "Got x = %1%, but we need x >= 0", x, policies::policy<>());
  155. }
  156. if(x == 0)
  157. {
  158. return (v == 0) ? 1 : 0;
  159. }
  160. T I, K;
  161. boost::math::detail::bessel_ik(v, x, &I, &K, 0xffff, policies::policy<>());
  162. std::cout << progress++ << ": I(" << v << ", " << x << ") = " << I << std::endl;
  163. if(fabs(I) > 1e30)
  164. throw std::domain_error("");
  165. return I;
  166. }
  167. template <class T>
  168. T cyl_bessel_k_bare(T v, T x)
  169. {
  170. using namespace std;
  171. if(x < 0)
  172. {
  173. return policies::raise_domain_error<T>(
  174. "",
  175. "Got x = %1%, but we need x > 0", x, policies::policy<>());
  176. }
  177. if(x == 0)
  178. {
  179. return (v == 0) ? policies::raise_overflow_error<T>("", 0, policies::policy<>())
  180. : policies::raise_domain_error<T>(
  181. "",
  182. "Got x = %1%, but we need x > 0", x, policies::policy<>());
  183. }
  184. T I, K;
  185. bessel_ik(v, x, &I, &K, 0xFFFF, policies::policy<>());
  186. std::cout << progress++ << ": K(" << v << ", " << x << ") = " << K << std::endl;
  187. if(fabs(K) > 1e30)
  188. throw std::domain_error("");
  189. return K;
  190. }
  191. template <class T>
  192. T cyl_neumann_bare(T v, T x)
  193. {
  194. T j, y;
  195. bessel_jy(v, x, &j, &y, 0xFFFF, policies::policy<>());
  196. std::cout << progress++ << ": Y(" << v << ", " << x << ") = " << y << std::endl;
  197. if(fabs(y) > 1e30)
  198. throw std::domain_error("");
  199. return y;
  200. }
  201. template <class T>
  202. T sph_bessel_j_bare(T v, T x)
  203. {
  204. std::cout << progress++ << ": j(" << v << ", " << x << ") = ";
  205. if((v < 0) || (floor(v) != v))
  206. throw std::domain_error("");
  207. T r = sqrt(constants::pi<T>() / (2 * x)) * cyl_bessel_j_bare(v+0.5, x);
  208. std::cout << r << std::endl;
  209. return r;
  210. }
  211. template <class T>
  212. T sph_bessel_y_bare(T v, T x)
  213. {
  214. std::cout << progress++ << ": y(" << v << ", " << x << ") = ";
  215. if((v < 0) || (floor(v) != v))
  216. throw std::domain_error("");
  217. T r = sqrt(constants::pi<T>() / (2 * x)) * cyl_neumann_bare(v+0.5, x);
  218. std::cout << r << std::endl;
  219. return r;
  220. }
  221. enum
  222. {
  223. func_J = 0,
  224. func_Y,
  225. func_I,
  226. func_K,
  227. func_j,
  228. func_y
  229. };
  230. int main(int argc, char* argv[])
  231. {
  232. std::cout << std::setprecision(17) << std::scientific;
  233. std::cout << sph_bessel_j_bare(0., 0.1185395751953125e4) << std::endl;
  234. std::cout << sph_bessel_j_bare(22., 0.6540834903717041015625) << std::endl;
  235. std::cout << std::setprecision(40) << std::scientific;
  236. parameter_info<mp_t> arg1, arg2;
  237. test_data<mp_t> data;
  238. int functype = 0;
  239. std::string letter = "J";
  240. if(argc == 2)
  241. {
  242. if(std::strcmp(argv[1], "--Y") == 0)
  243. {
  244. functype = func_Y;
  245. letter = "Y";
  246. }
  247. else if(std::strcmp(argv[1], "--I") == 0)
  248. {
  249. functype = func_I;
  250. letter = "I";
  251. }
  252. else if(std::strcmp(argv[1], "--K") == 0)
  253. {
  254. functype = func_K;
  255. letter = "K";
  256. }
  257. else if(std::strcmp(argv[1], "--j") == 0)
  258. {
  259. functype = func_j;
  260. letter = "j";
  261. }
  262. else if(std::strcmp(argv[1], "--y") == 0)
  263. {
  264. functype = func_y;
  265. letter = "y";
  266. }
  267. else
  268. BOOST_ASSERT(0);
  269. }
  270. bool cont;
  271. std::string line;
  272. std::cout << "Welcome.\n"
  273. "This program will generate spot tests for the Bessel " << letter << " function\n\n";
  274. do{
  275. get_user_parameter_info(arg1, "v");
  276. get_user_parameter_info(arg2, "x");
  277. mp_t (*fp)(mp_t, mp_t);
  278. if(functype == func_J)
  279. fp = cyl_bessel_j_bare;
  280. else if(functype == func_I)
  281. fp = cyl_bessel_i_bare;
  282. else if(functype == func_K)
  283. fp = cyl_bessel_k_bare;
  284. else if(functype == func_Y)
  285. fp = cyl_neumann_bare;
  286. else if(functype == func_j)
  287. fp = sph_bessel_j_bare;
  288. else if(functype == func_y)
  289. fp = sph_bessel_y_bare;
  290. else
  291. BOOST_ASSERT(0);
  292. data.insert(fp, arg1, arg2);
  293. std::cout << "Any more data [y/n]?";
  294. std::getline(std::cin, line);
  295. boost::algorithm::trim(line);
  296. cont = (line == "y");
  297. }while(cont);
  298. std::cout << "Enter name of test data file [default=bessel_j_data.ipp]";
  299. std::getline(std::cin, line);
  300. boost::algorithm::trim(line);
  301. if(line == "")
  302. line = "bessel_j_data.ipp";
  303. std::ofstream ofs(line.c_str());
  304. line.erase(line.find('.'));
  305. ofs << std::scientific << std::setprecision(40);
  306. write_code(ofs, data, line.c_str());
  307. return 0;
  308. }