zbesh.cpp 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455
  1. /* zbesh.f -- translated by f2c (version 20100827).
  2. This file no longer depends on f2c.
  3. */
  4. #include "slatec-internal.hpp"
  5. /* Table of constant values */
  6. static integer const c__4 = 4;
  7. static integer const c__15 = 15;
  8. static integer const c__16 = 16;
  9. static integer const c__5 = 5;
  10. static integer const c__14 = 14;
  11. static integer const c__9 = 9;
  12. static integer const c__1 = 1;
  13. static integer const c__2 = 2;
  14. int zbesh_(double *zr, double *zi, double const *fnu,
  15. integer const *kode, integer const *m, integer const *n, double *cyr, double *
  16. cyi, integer *nz, integer *ierr)
  17. {
  18. /* Initialized data */
  19. static double const hpi = 1.57079632679489662;
  20. /* System generated locals */
  21. integer i__1, i__2;
  22. double d__1, d__2;
  23. /* Local variables */
  24. integer i__, k, k1, k2;
  25. double aa, bb, fn;
  26. integer mm;
  27. double az;
  28. integer ir, nn;
  29. double rl;
  30. integer mr, nw;
  31. double dig, arg, aln, fmm, r1m5, ufl, sgn;
  32. integer nuf, inu;
  33. double tol, sti, zni, zti, str, znr, alim, elim;
  34. double atol, rhpi;
  35. integer inuh;
  36. double fnul, rtol, ascle, csgni;
  37. double csgnr;
  38. /* ***BEGIN PROLOGUE ZBESH */
  39. /* ***PURPOSE Compute a sequence of the Hankel functions H(m,a,z) */
  40. /* for superscript m=1 or 2, real nonnegative orders a=b, */
  41. /* b+1,... where b>0, and nonzero complex argument z. A */
  42. /* scaling option is available to help avoid overflow. */
  43. /* ***LIBRARY SLATEC */
  44. /* ***CATEGORY C10A4 */
  45. /* ***TYPE COMPLEX (CBESH-C, ZBESH-C) */
  46. /* ***KEYWORDS BESSEL FUNCTIONS OF COMPLEX ARGUMENT, */
  47. /* BESSEL FUNCTIONS OF THE THIRD KIND, H BESSEL FUNCTIONS, */
  48. /* HANKEL FUNCTIONS */
  49. /* ***AUTHOR Amos, D. E., (SNL) */
  50. /* ***DESCRIPTION */
  51. /* ***A DOUBLE PRECISION ROUTINE*** */
  52. /* On KODE=1, ZBESH computes an N member sequence of complex */
  53. /* Hankel (Bessel) functions CY(L)=H(M,FNU+L-1,Z) for super- */
  54. /* script M=1 or 2, real nonnegative orders FNU+L-1, L=1,..., */
  55. /* N, and complex nonzero Z in the cut plane -pi<arg(Z)<=pi */
  56. /* where Z=ZR+i*ZI. On KODE=2, CBESH returns the scaled */
  57. /* functions */
  58. /* CY(L) = H(M,FNU+L-1,Z)*exp(-(3-2*M)*Z*i), i**2=-1 */
  59. /* which removes the exponential behavior in both the upper */
  60. /* and lower half planes. Definitions and notation are found */
  61. /* in the NBS Handbook of Mathematical Functions (Ref. 1). */
  62. /* Input */
  63. /* ZR - DOUBLE PRECISION REAL part of nonzero argument Z */
  64. /* ZI - DOUBLE PRECISION imag part of nonzero argument Z */
  65. /* FNU - DOUBLE PRECISION initial order, FNU>=0 */
  66. /* KODE - A parameter to indicate the scaling option */
  67. /* KODE=1 returns */
  68. /* CY(L)=H(M,FNU+L-1,Z), L=1,...,N */
  69. /* =2 returns */
  70. /* CY(L)=H(M,FNU+L-1,Z)*exp(-(3-2M)*Z*i), */
  71. /* L=1,...,N */
  72. /* M - Superscript of Hankel function, M=1 or 2 */
  73. /* N - Number of terms in the sequence, N>=1 */
  74. /* Output */
  75. /* CYR - DOUBLE PRECISION REAL part of result vector */
  76. /* CYI - DOUBLE PRECISION imag part of result vector */
  77. /* NZ - Number of underflows set to zero */
  78. /* NZ=0 Normal return */
  79. /* NZ>0 CY(L)=0 for NZ values of L (if M=1 and */
  80. /* Im(Z)>0 or if M=2 and Im(Z)<0, then */
  81. /* CY(L)=0 for L=1,...,NZ; in the com- */
  82. /* plementary half planes, the underflows */
  83. /* may not be in an uninterrupted sequence) */
  84. /* IERR - Error flag */
  85. /* IERR=0 Normal return - COMPUTATION COMPLETED */
  86. /* IERR=1 Input error - NO COMPUTATION */
  87. /* IERR=2 Overflow - NO COMPUTATION */
  88. /* (abs(Z) too small and/or FNU+N-1 */
  89. /* too large) */
  90. /* IERR=3 Precision warning - COMPUTATION COMPLETED */
  91. /* (Result has half precision or less */
  92. /* because abs(Z) or FNU+N-1 is large) */
  93. /* IERR=4 Precision error - NO COMPUTATION */
  94. /* (Result has no precision because */
  95. /* abs(Z) or FNU+N-1 is too large) */
  96. /* IERR=5 Algorithmic error - NO COMPUTATION */
  97. /* (Termination condition not met) */
  98. /* *Long Description: */
  99. /* The computation is carried out by the formula */
  100. /* H(m,a,z) = (1/t)*exp(-a*t)*K(a,z*exp(-t)) */
  101. /* t = (3-2*m)*i*pi/2 */
  102. /* where the K Bessel function is computed as described in the */
  103. /* prologue to CBESK. */
  104. /* Exponential decay of H(m,a,z) occurs in the upper half z */
  105. /* plane for m=1 and the lower half z plane for m=2. Exponential */
  106. /* growth occurs in the complementary half planes. Scaling */
  107. /* by exp(-(3-2*m)*z*i) removes the exponential behavior in the */
  108. /* whole z plane as z goes to infinity. */
  109. /* For negative orders, the formula */
  110. /* H(m,-a,z) = H(m,a,z)*exp((3-2*m)*a*pi*i) */
  111. /* can be used. */
  112. /* In most complex variable computation, one must evaluate ele- */
  113. /* mentary functions. When the magnitude of Z or FNU+N-1 is */
  114. /* large, losses of significance by argument reduction occur. */
  115. /* Consequently, if either one exceeds U1=SQRT(0.5/UR), then */
  116. /* losses exceeding half precision are likely and an error flag */
  117. /* IERR=3 is triggered where UR=MAX(D1MACH(4),1.0D-18) is double */
  118. /* precision unit roundoff limited to 18 digits precision. Also, */
  119. /* if either is larger than U2=0.5/UR, then all significance is */
  120. /* lost and IERR=4. In order to use the INT function, arguments */
  121. /* must be further restricted not to exceed the largest machine */
  122. /* integer, U3=I1MACH(9). Thus, the magnitude of Z and FNU+N-1 */
  123. /* is restricted by MIN(U2,U3). In IEEE arithmetic, U1,U2, and */
  124. /* U3 approximate 2.0E+3, 4.2E+6, 2.1E+9 in single precision */
  125. /* and 4.7E+7, 2.3E+15 and 2.1E+9 in double precision. This */
  126. /* makes U2 limiting in single precision and U3 limiting in */
  127. /* double precision. This means that one can expect to retain, */
  128. /* in the worst cases on IEEE machines, no digits in single pre- */
  129. /* cision and only 6 digits in double precision. Similar con- */
  130. /* siderations hold for other machines. */
  131. /* The approximate relative error in the magnitude of a complex */
  132. /* Bessel function can be expressed as P*10**S where P=MAX(UNIT */
  133. /* ROUNDOFF,1.0E-18) is the nominal precision and 10**S repre- */
  134. /* sents the increase in error due to argument reduction in the */
  135. /* elementary functions. Here, S=MAX(1,ABS(LOG10(ABS(Z))), */
  136. /* ABS(LOG10(FNU))) approximately (i.e., S=MAX(1,ABS(EXPONENT OF */
  137. /* ABS(Z),ABS(EXPONENT OF FNU)) ). However, the phase angle may */
  138. /* have only absolute accuracy. This is most likely to occur */
  139. /* when one component (in magnitude) is larger than the other by */
  140. /* several orders of magnitude. If one component is 10**K larger */
  141. /* than the other, then one can expect only MAX(ABS(LOG10(P))-K, */
  142. /* 0) significant digits; or, stated another way, when K exceeds */
  143. /* the exponent of P, no significant digits remain in the smaller */
  144. /* component. However, the phase angle retains absolute accuracy */
  145. /* because, in complex arithmetic with precision P, the smaller */
  146. /* component will not (as a rule) decrease below P times the */
  147. /* magnitude of the larger component. In these extreme cases, */
  148. /* the principal phase angle is on the order of +P, -P, PI/2-P, */
  149. /* or -PI/2+P. */
  150. /* ***REFERENCES 1. M. Abramowitz and I. A. Stegun, Handbook of Mathe- */
  151. /* matical Functions, National Bureau of Standards */
  152. /* Applied Mathematics Series 55, U. S. Department */
  153. /* of Commerce, Tenth Printing (1972) or later. */
  154. /* 2. D. E. Amos, Computation of Bessel Functions of */
  155. /* Complex Argument, Report SAND83-0086, Sandia National */
  156. /* Laboratories, Albuquerque, NM, May 1983. */
  157. /* 3. D. E. Amos, Computation of Bessel Functions of */
  158. /* Complex Argument and Large Order, Report SAND83-0643, */
  159. /* Sandia National Laboratories, Albuquerque, NM, May */
  160. /* 1983. */
  161. /* 4. D. E. Amos, A Subroutine Package for Bessel Functions */
  162. /* of a Complex Argument and Nonnegative Order, Report */
  163. /* SAND85-1018, Sandia National Laboratory, Albuquerque, */
  164. /* NM, May 1985. */
  165. /* 5. D. E. Amos, A portable package for Bessel functions */
  166. /* of a complex argument and nonnegative order, ACM */
  167. /* Transactions on Mathematical Software, 12 (September */
  168. /* 1986), pp. 265-273. */
  169. /* ***ROUTINES CALLED D1MACH, I1MACH, ZABS, ZACON, ZBKNU, ZBUNK, ZUOIK */
  170. /* ***REVISION HISTORY (YYMMDD) */
  171. /* 830501 DATE WRITTEN */
  172. /* 890801 REVISION DATE from Version 3.2 */
  173. /* 910415 Prologue converted to Version 4.0 format. (BAB) */
  174. /* 920128 Category corrected. (WRB) */
  175. /* 920811 Prologue revised. (DWL) */
  176. /* ***END PROLOGUE ZBESH */
  177. /* COMPLEX CY,Z,ZN,ZT,CSGN */
  178. /* Parameter adjustments */
  179. --cyi;
  180. --cyr;
  181. /* Function Body */
  182. /* ***FIRST EXECUTABLE STATEMENT ZBESH */
  183. *ierr = 0;
  184. *nz = 0;
  185. if (*zr == 0. && *zi == 0.) {
  186. *ierr = 1;
  187. }
  188. if (*fnu < 0.) {
  189. *ierr = 1;
  190. }
  191. if (*m < 1 || *m > 2) {
  192. *ierr = 1;
  193. }
  194. if (*kode < 1 || *kode > 2) {
  195. *ierr = 1;
  196. }
  197. if (*n < 1) {
  198. *ierr = 1;
  199. }
  200. if (*ierr != 0) {
  201. return 0;
  202. }
  203. nn = *n;
  204. /* ----------------------------------------------------------------------- */
  205. /* SET PARAMETERS RELATED TO MACHINE CONSTANTS. */
  206. /* TOL IS THE APPROXIMATE UNIT ROUNDOFF LIMITED TO 1.0E-18. */
  207. /* ELIM IS THE APPROXIMATE EXPONENTIAL OVER- AND UNDERFLOW LIMIT. */
  208. /* EXP(-ELIM).LT.EXP(-ALIM)=EXP(-ELIM)/TOL AND */
  209. /* EXP(ELIM).GT.EXP(ALIM)=EXP(ELIM)*TOL ARE INTERVALS NEAR */
  210. /* UNDERFLOW AND OVERFLOW LIMITS WHERE SCALED ARITHMETIC IS DONE. */
  211. /* RL IS THE LOWER BOUNDARY OF THE ASYMPTOTIC EXPANSION FOR LARGE Z. */
  212. /* DIG = NUMBER OF BASE 10 DIGITS IN TOL = 10**(-DIG). */
  213. /* FNUL IS THE LOWER BOUNDARY OF THE ASYMPTOTIC SERIES FOR LARGE FNU */
  214. /* ----------------------------------------------------------------------- */
  215. /* Computing MAX */
  216. d__1 = d1mach_(4);
  217. tol = max(d__1,1e-18);
  218. k1 = i1mach_(15);
  219. k2 = i1mach_(16);
  220. r1m5 = d1mach_(5);
  221. /* Computing MIN */
  222. i__1 = abs(k1), i__2 = abs(k2);
  223. k = min(i__1,i__2);
  224. elim = (k * r1m5 - 3.) * 2.303;
  225. k1 = i1mach_(14) - 1;
  226. aa = r1m5 * k1;
  227. dig = min(aa,18.);
  228. aa *= 2.303;
  229. /* Computing MAX */
  230. d__1 = -aa;
  231. alim = elim + max(d__1,-41.45);
  232. fnul = (dig - 3.) * 6. + 10.;
  233. rl = dig * 1.2 + 3.;
  234. fn = *fnu + (nn - 1);
  235. mm = 3 - *m - *m;
  236. fmm = (double) mm;
  237. znr = fmm * *zi;
  238. zni = -fmm * *zr;
  239. /* ----------------------------------------------------------------------- */
  240. /* TEST FOR PROPER RANGE */
  241. /* ----------------------------------------------------------------------- */
  242. az = zabs_(zr, zi);
  243. aa = .5 / tol;
  244. bb = i1mach_(9) * .5;
  245. aa = min(aa,bb);
  246. if (az > aa) {
  247. goto L260;
  248. }
  249. if (fn > aa) {
  250. goto L260;
  251. }
  252. aa = sqrt(aa);
  253. if (az > aa) {
  254. *ierr = 3;
  255. }
  256. if (fn > aa) {
  257. *ierr = 3;
  258. }
  259. /* ----------------------------------------------------------------------- */
  260. /* OVERFLOW TEST ON THE LAST MEMBER OF THE SEQUENCE */
  261. /* ----------------------------------------------------------------------- */
  262. ufl = d1mach_(1) * 1e3;
  263. if (az < ufl) {
  264. goto L230;
  265. }
  266. if (*fnu > fnul) {
  267. goto L90;
  268. }
  269. if (fn <= 1.) {
  270. goto L70;
  271. }
  272. if (fn > 2.) {
  273. goto L60;
  274. }
  275. if (az > tol) {
  276. goto L70;
  277. }
  278. arg = az * .5;
  279. aln = -fn * log(arg);
  280. if (aln > elim) {
  281. goto L230;
  282. }
  283. goto L70;
  284. L60:
  285. zuoik_(&znr, &zni, fnu, kode, &c__2, &nn, &cyr[1], &cyi[1], &nuf, &tol, &
  286. elim, &alim);
  287. if (nuf < 0) {
  288. goto L230;
  289. }
  290. *nz += nuf;
  291. nn -= nuf;
  292. /* ----------------------------------------------------------------------- */
  293. /* HERE NN=N OR NN=0 SINCE NUF=0,NN, OR -1 ON RETURN FROM CUOIK */
  294. /* IF NUF=NN, THEN CY(I)=CZERO FOR ALL I */
  295. /* ----------------------------------------------------------------------- */
  296. if (nn == 0) {
  297. goto L140;
  298. }
  299. L70:
  300. if (znr < 0. || znr == 0. && zni < 0. && *m == 2) {
  301. goto L80;
  302. }
  303. /* ----------------------------------------------------------------------- */
  304. /* RIGHT HALF PLANE COMPUTATION, XN.GE.0. .AND. (XN.NE.0. .OR. */
  305. /* YN.GE.0. .OR. M=1) */
  306. /* ----------------------------------------------------------------------- */
  307. zbknu_(&znr, &zni, fnu, kode, &nn, &cyr[1], &cyi[1], nz, &tol, &elim, &
  308. alim);
  309. goto L110;
  310. /* ----------------------------------------------------------------------- */
  311. /* LEFT HALF PLANE COMPUTATION */
  312. /* ----------------------------------------------------------------------- */
  313. L80:
  314. mr = -mm;
  315. zacon_(&znr, &zni, fnu, kode, &mr, &nn, &cyr[1], &cyi[1], &nw, &rl, &fnul,
  316. &tol, &elim, &alim);
  317. if (nw < 0) {
  318. goto L240;
  319. }
  320. *nz = nw;
  321. goto L110;
  322. L90:
  323. /* ----------------------------------------------------------------------- */
  324. /* UNIFORM ASYMPTOTIC EXPANSIONS FOR FNU.GT.FNUL */
  325. /* ----------------------------------------------------------------------- */
  326. mr = 0;
  327. if (znr >= 0. && (znr != 0. || zni >= 0. || *m != 2)) {
  328. goto L100;
  329. }
  330. mr = -mm;
  331. if (znr != 0. || zni >= 0.) {
  332. goto L100;
  333. }
  334. znr = -znr;
  335. zni = -zni;
  336. L100:
  337. zbunk_(&znr, &zni, fnu, kode, &mr, &nn, &cyr[1], &cyi[1], &nw, &tol, &
  338. elim, &alim);
  339. if (nw < 0) {
  340. goto L240;
  341. }
  342. *nz += nw;
  343. L110:
  344. /* ----------------------------------------------------------------------- */
  345. /* H(M,FNU,Z) = -FMM*(I/HPI)*(ZT**FNU)*K(FNU,-Z*ZT) */
  346. /* ZT=EXP(-FMM*HPI*I) = CMPLX(0.0,-FMM), FMM=3-2*M, M=1,2 */
  347. /* ----------------------------------------------------------------------- */
  348. d__1 = -fmm;
  349. sgn = f2c::d_sign(&hpi, &d__1);
  350. /* ----------------------------------------------------------------------- */
  351. /* CALCULATE EXP(FNU*HPI*I) TO MINIMIZE LOSSES OF SIGNIFICANCE */
  352. /* WHEN FNU IS LARGE */
  353. /* ----------------------------------------------------------------------- */
  354. inu = (integer) (*fnu);
  355. inuh = inu / 2;
  356. ir = inu - (inuh << 1);
  357. arg = (*fnu - (inu - ir)) * sgn;
  358. rhpi = 1. / sgn;
  359. /* ZNI = RHPI*COS(ARG) */
  360. /* ZNR = -RHPI*SIN(ARG) */
  361. csgni = rhpi * cos(arg);
  362. csgnr = -rhpi * sin(arg);
  363. if (inuh % 2 == 0) {
  364. goto L120;
  365. }
  366. /* ZNR = -ZNR */
  367. /* ZNI = -ZNI */
  368. csgnr = -csgnr;
  369. csgni = -csgni;
  370. L120:
  371. zti = -fmm;
  372. rtol = 1. / tol;
  373. ascle = ufl * rtol;
  374. i__1 = nn;
  375. for (i__ = 1; i__ <= i__1; ++i__) {
  376. /* STR = CYR(I)*ZNR - CYI(I)*ZNI */
  377. /* CYI(I) = CYR(I)*ZNI + CYI(I)*ZNR */
  378. /* CYR(I) = STR */
  379. /* STR = -ZNI*ZTI */
  380. /* ZNI = ZNR*ZTI */
  381. /* ZNR = STR */
  382. aa = cyr[i__];
  383. bb = cyi[i__];
  384. atol = 1.;
  385. /* Computing MAX */
  386. d__1 = abs(aa), d__2 = abs(bb);
  387. if (max(d__1,d__2) > ascle) {
  388. goto L135;
  389. }
  390. aa *= rtol;
  391. bb *= rtol;
  392. atol = tol;
  393. L135:
  394. str = aa * csgnr - bb * csgni;
  395. sti = aa * csgni + bb * csgnr;
  396. cyr[i__] = str * atol;
  397. cyi[i__] = sti * atol;
  398. str = -csgni * zti;
  399. csgni = csgnr * zti;
  400. csgnr = str;
  401. /* L130: */
  402. }
  403. return 0;
  404. L140:
  405. if (znr < 0.) {
  406. goto L230;
  407. }
  408. return 0;
  409. L230:
  410. *nz = 0;
  411. *ierr = 2;
  412. return 0;
  413. L240:
  414. if (nw == -1) {
  415. goto L230;
  416. }
  417. *nz = 0;
  418. *ierr = 5;
  419. return 0;
  420. L260:
  421. *nz = 0;
  422. *ierr = 4;
  423. return 0;
  424. } /* zbesh_ */