bessel.cc 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325
  1. //**********************************************************************************//
  2. // Copyright (C) 2009-2015 Ovidio Pena <ovidio@bytesfall.com> //
  3. // Copyright (C) 2013-2015 Konstantin Ladutenko <kostyfisik@gmail.com> //
  4. // //
  5. // This file is part of scattnlay //
  6. // //
  7. // This program is free software: you can redistribute it and/or modify //
  8. // it under the terms of the GNU General Public License as published by //
  9. // the Free Software Foundation, either version 3 of the License, or //
  10. // (at your option) any later version. //
  11. // //
  12. // This program is distributed in the hope that it will be useful, //
  13. // but WITHOUT ANY WARRANTY; without even the implied warranty of //
  14. // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the //
  15. // GNU General Public License for more details. //
  16. // //
  17. // The only additional remark is that we expect that all publications //
  18. // describing work using this software, or all commercial products //
  19. // using it, cite the following reference: //
  20. // [1] O. Pena and U. Pal, "Scattering of electromagnetic radiation by //
  21. // a multilayered sphere," Computer Physics Communications, //
  22. // vol. 180, Nov. 2009, pp. 2348-2354. //
  23. // //
  24. // You should have received a copy of the GNU General Public License //
  25. // along with this program. If not, see <http://www.gnu.org/licenses/>. //
  26. //**********************************************************************************//
  27. #include "bessel.h"
  28. #include <algorithm>
  29. #include <complex>
  30. #include <cmath>
  31. #include <limits>
  32. #include <stdexcept>
  33. #include <vector>
  34. namespace nmie {
  35. namespace bessel {
  36. // !*****************************************************************************80
  37. //
  38. // C++ port of fortran code
  39. //
  40. // !! CSPHJY: spherical Bessel functions jn(z) and yn(z) for complex argument.
  41. // !
  42. // ! Discussion:
  43. // !
  44. // ! This procedure computes spherical Bessel functions jn(z) and yn(z)
  45. // ! and their derivatives for a complex argument.
  46. // !
  47. // ! Licensing:
  48. // !
  49. // ! This routine is copyrighted by Shanjie Zhang and Jianming Jin. However,
  50. // ! they give permission to incorporate this routine into a user program
  51. // ! provided that the copyright is acknowledged.
  52. // !
  53. // ! Modified:
  54. // !
  55. // ! 01 August 2012
  56. // !
  57. // ! Author:
  58. // !
  59. // ! Shanjie Zhang, Jianming Jin
  60. // !
  61. // ! Reference:
  62. // !
  63. // ! Shanjie Zhang, Jianming Jin,
  64. // ! Computation of Special Functions,
  65. // ! Wiley, 1996,
  66. // ! ISBN: 0-471-11963-6,
  67. // ! LC: QA351.C45.
  68. // !
  69. // ! Parameters:
  70. // !
  71. // ! Input, integer ( kind = 4 ) N, the order of jn(z) and yn(z).
  72. // !
  73. // ! Input, complex ( kind = 8 ) Z, the argument.
  74. // !
  75. // ! Output, integer ( kind = 4 ) NM, the highest order computed.
  76. // !
  77. // ! Output, complex ( kind = 8 ) CSJ(0:N0, CDJ(0:N), CSY(0:N), CDY(0:N),
  78. // ! the values of jn(z), jn'(z), yn(z), yn'(z).
  79. // !
  80. void csphjy (int n, std::complex<double>z, int& nm,
  81. std::vector< std::complex<double> >& csj,
  82. std::vector< std::complex<double> >& cdj,
  83. std::vector< std::complex<double> >& csy,
  84. std::vector< std::complex<double> >& cdy ) {
  85. double a0;
  86. csj.resize(n+1);
  87. cdj.resize(n+1);
  88. csy.resize(n+1);
  89. cdy.resize(n+1);
  90. std::complex<double> cf, cf0, cf1, cs, csa, csb;
  91. int k, m;
  92. a0 = std::abs(z);
  93. nm = n;
  94. if (a0 < 1.0e-60) {
  95. for (int k = 0; k < n+1; ++k) {
  96. csj[k] = 0.0;
  97. cdj[k] = 0.0;
  98. csy[k] = -1.0e+300;
  99. cdy[k] = 1.0e+300;
  100. }
  101. csj[0] = std::complex<double>( 1.0, 0.0);
  102. cdj[1] = std::complex<double>( 0.333333333333333, 0.0);
  103. return;
  104. }
  105. csj[0] = std::sin ( z ) / z;
  106. csj[1] = ( csj[0] - std::cos ( z ) ) / z;
  107. if ( 2 <= n ) {
  108. csa = csj[0];
  109. csb = csj[1];
  110. m = msta1 ( a0, 200 );
  111. if ( m < n ) nm = m;
  112. else m = msta2 ( a0, n, 15 );
  113. cf0 = 0.0;
  114. cf1 = 1.0e-100;
  115. for (int k = m; k>=0; --k) {
  116. cf = ( 2.0 * k + 3.0 ) * cf1 / z - cf0;
  117. if ( k <= nm ) csj[k] = cf;
  118. cf0 = cf1;
  119. cf1 = cf;
  120. }
  121. if ( std::abs ( csa ) <= std::abs ( csb ) ) cs = csb / cf0;
  122. else cs = csa / cf;
  123. for (int k = 0; k <= nm; ++k) {
  124. csj[k] = cs * csj[k];
  125. }
  126. }
  127. cdj[0] = ( std::cos ( z ) - std::sin ( z ) / z ) / z;
  128. for (int k = 1; k <=nm; ++k) {
  129. cdj[k] = csj[k-1] - ( k + 1.0 ) * csj[k] / z;
  130. }
  131. csy[0] = - std::cos ( z ) / z;
  132. csy[1] = ( csy[0] - std::sin ( z ) ) / z;
  133. cdy[0] = ( std::sin ( z ) + std::cos ( z ) / z ) / z;
  134. cdy[1] = ( 2.0 * cdy[0] - std::cos ( z ) ) / z;
  135. for (int k = 2; k<=nm; ++k) {
  136. if ( std::abs ( csj[k-2] ) < std::abs ( csj[k-1] ) ) {
  137. csy[k] = ( csj[k] * csy[k-1] - 1.0 / ( z * z ) ) / csj[k-1];
  138. } else {
  139. csy[k] = ( csj[k] * csy[k-2] - ( 2.0 * k - 1.0 ) / std::pow(z,3) )
  140. / csj[k-2];
  141. }
  142. }
  143. for (int k = 2; k<=nm; ++k) {
  144. cdy[k] = csy[k-1] - ( k + 1.0 ) * csy[k] / z;
  145. }
  146. return;
  147. }
  148. // function msta2 ( x, n, mp )
  149. // !*****************************************************************************80
  150. // !
  151. // !! MSTA2 determines a backward recurrence starting point for Jn(x).
  152. // !
  153. // ! Discussion:
  154. // !
  155. // ! This procedure determines the starting point for a backward
  156. // ! recurrence such that all Jn(x) has MP significant digits.
  157. // !
  158. // ! Licensing:
  159. // !
  160. // ! This routine is copyrighted by Shanjie Zhang and Jianming Jin. However,
  161. // ! they give permission to incorporate this routine into a user program
  162. // ! provided that the copyright is acknowledged.
  163. // !
  164. // ! Modified:
  165. // !
  166. // ! 08 July 2012
  167. // !
  168. // ! Author:
  169. // !
  170. // ! Shanjie Zhang, Jianming Jin
  171. // !
  172. // ! Reference:
  173. // !
  174. // ! Shanjie Zhang, Jianming Jin,
  175. // ! Computation of Special Functions,
  176. // ! Wiley, 1996,
  177. // ! ISBN: 0-471-11963-6,
  178. // ! LC: QA351.C45.
  179. // !
  180. // ! Parameters:
  181. // !
  182. // ! Input, real ( kind = 8 ) X, the argument of Jn(x).
  183. // !
  184. // ! Input, integer ( kind = 4 ) N, the order of Jn(x).
  185. // !
  186. // ! Input, integer ( kind = 4 ) MP, the number of significant digits.
  187. // !
  188. // ! Output, integer ( kind = 4 ) MSTA2, the starting point.
  189. // !
  190. int msta2 ( double x, int n, int mp ) {
  191. double a0, ejn, f, f0, f1, hmp;
  192. int it, n0, n1, nn;
  193. double obj;
  194. a0 = std::abs ( x );
  195. hmp = 0.5 * mp;
  196. ejn = envj ( n, a0 );
  197. if ( ejn <= hmp ) {
  198. obj = mp;
  199. n0 = static_cast<int> ( 1.1 * a0 );
  200. } else {
  201. obj = hmp + ejn;
  202. n0 = n;
  203. }
  204. f0 = envj ( n0, a0 ) - obj;
  205. n1 = n0 + 5;
  206. f1 = envj ( n1, a0 ) - obj;
  207. for (int it = 1; it < 21; ++it) {
  208. nn = n1 - ( n1 - n0 ) / ( 1.0 - f0 / f1 );
  209. f = envj ( nn, a0 ) - obj;
  210. if ( std::abs ( nn - n1 ) < 1 ) break;
  211. n0 = n1;
  212. f0 = f1;
  213. n1 = nn;
  214. f1 = f;
  215. }
  216. return nn + 10;
  217. }
  218. // !*****************************************************************************80
  219. // !
  220. // !! MSTA1 determines a backward recurrence starting point for Jn(x).
  221. // !
  222. // ! Discussion:
  223. // !
  224. // ! This procedure determines the starting point for backward
  225. // ! recurrence such that the magnitude of
  226. // ! Jn(x) at that point is about 10^(-MP).
  227. // !
  228. // ! Licensing:
  229. // !
  230. // ! This routine is copyrighted by Shanjie Zhang and Jianming Jin. However,
  231. // ! they give permission to incorporate this routine into a user program
  232. // ! provided that the copyright is acknowledged.
  233. // !
  234. // ! Modified:
  235. // !
  236. // ! 08 July 2012
  237. // !
  238. // ! Author:
  239. // !
  240. // ! Shanjie Zhang, Jianming Jin
  241. // !
  242. // ! Reference:
  243. // !
  244. // ! Shanjie Zhang, Jianming Jin,
  245. // ! Computation of Special Functions,
  246. // ! Wiley, 1996,
  247. // ! ISBN: 0-471-11963-6,
  248. // ! LC: QA351.C45.
  249. // !
  250. // ! Parameters:
  251. // !
  252. // ! Input, real ( kind = 8 ) X, the argument.
  253. // !
  254. // ! Input, integer ( kind = 4 ) MP, the negative logarithm of the
  255. // ! desired magnitude.
  256. // !
  257. // ! Output, integer ( kind = 4 ) MSTA1, the starting point.
  258. // !
  259. int msta1 ( double x, int mp ) {
  260. double a0, f, f0, f1;
  261. int it, n0, n1, nn;
  262. a0 = std::abs ( x );
  263. n0 = static_cast<int> (1.1 * a0 ) + 1;
  264. f0 = envj ( n0, a0 ) - mp;
  265. n1 = n0 + 5;
  266. f1 = envj ( n1, a0 ) - mp;
  267. for (int it = 1; it <= 20; ++it) {
  268. nn = n1 - ( n1 - n0 ) / ( 1.0 - f0 / f1 );
  269. f = envj ( nn, a0 ) - mp;
  270. if ( abs ( nn - n1 ) < 1 ) break;
  271. n0 = n1;
  272. f0 = f1;
  273. n1 = nn;
  274. f1 = f;
  275. }
  276. return nn;
  277. }
  278. // function envj ( n, x )
  279. // !*****************************************************************************80
  280. // !
  281. // !! ENVJ is a utility function used by MSTA1 and MSTA2.
  282. // !
  283. // ! Licensing:
  284. // !
  285. // ! This routine is copyrighted by Shanjie Zhang and Jianming Jin. However,
  286. // ! they give permission to incorporate this routine into a user program
  287. // ! provided that the copyright is acknowledged.
  288. // !
  289. // ! Modified:
  290. // !
  291. // ! 14 March 2012
  292. // !
  293. // ! Author:
  294. // !
  295. // ! Shanjie Zhang, Jianming Jin
  296. // !
  297. // ! Reference:
  298. // !
  299. // ! Shanjie Zhang, Jianming Jin,
  300. // ! Computation of Special Functions,
  301. // ! Wiley, 1996,
  302. // ! ISBN: 0-471-11963-6,
  303. // ! LC: QA351.C45.
  304. // !
  305. // ! Parameters:
  306. // !
  307. // ! Input, integer ( kind = 4 ) N, ?
  308. // !
  309. // ! Input, real ( kind = 8 ) X, ?
  310. // !
  311. // ! Output, real ( kind = 8 ) ENVJ, ?
  312. // !
  313. double envj (int n, double x ) {
  314. return 0.5 * std::log10(6.28 * n) - n * std::log10(1.36 * x / static_cast<double>(n) );
  315. }
  316. } // end of namespace bessel
  317. } // end of namespace nmie