gsl_specfunc__mathieu_charv.c 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850
  1. /* specfunc/mathieu_charv.c
  2. *
  3. * Copyright (C) 2002 Lowell Johnson
  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., 675 Mass Ave, Cambridge, MA 02139, USA.
  18. */
  19. /* Author: L. Johnson */
  20. #include "gsl__config.h"
  21. #include <stdlib.h>
  22. #include <stdio.h>
  23. #include <math.h>
  24. #include "gsl_math.h"
  25. #include "gsl_eigen.h"
  26. #include "gsl_errno.h"
  27. #include "gsl_sf_mathieu.h"
  28. /* prototypes */
  29. static double solve_cubic(double c2, double c1, double c0);
  30. static double ceer(int order, double qq, double aa, int nterms)
  31. {
  32. double term, term1;
  33. int ii, n1;
  34. if (order == 0)
  35. term = 0.0;
  36. else
  37. {
  38. term = 2.0*qq*qq/aa;
  39. if (order != 2)
  40. {
  41. n1 = order/2 - 1;
  42. for (ii=0; ii<n1; ii++)
  43. term = qq*qq/(aa - 4.0*(ii+1)*(ii+1) - term);
  44. }
  45. }
  46. term += order*order;
  47. term1 = 0.0;
  48. for (ii=0; ii<nterms; ii++)
  49. term1 = qq*qq/
  50. (aa - (order + 2.0*(nterms - ii))*(order + 2.0*(nterms - ii)) - term1);
  51. if (order == 0)
  52. term1 *= 2.0;
  53. return (term + term1 - aa);
  54. }
  55. static double ceor(int order, double qq, double aa, int nterms)
  56. {
  57. double term, term1;
  58. int ii, n1;
  59. term = qq;
  60. n1 = (int)((float)order/2.0 - 0.5);
  61. for (ii=0; ii<n1; ii++)
  62. term = qq*qq/(aa - (2.0*ii + 1.0)*(2.0*ii + 1.0) - term);
  63. term += order*order;
  64. term1 = 0.0;
  65. for (ii=0; ii<nterms; ii++)
  66. term1 = qq*qq/
  67. (aa - (order + 2.0*(nterms - ii))*(order + 2.0*(nterms - ii)) - term1);
  68. return (term + term1 - aa);
  69. }
  70. static double seer(int order, double qq, double aa, int nterms)
  71. {
  72. double term, term1;
  73. int ii, n1;
  74. term = 0.0;
  75. n1 = order/2 - 1;
  76. for (ii=0; ii<n1; ii++)
  77. term = qq*qq/(aa - 4.0*(ii + 1)*(ii + 1) - term);
  78. term += order*order;
  79. term1 = 0.0;
  80. for (ii=0; ii<nterms; ii++)
  81. term1 = qq*qq/
  82. (aa - (order + 2.0*(nterms - ii))*(order + 2.0*(nterms - ii)) - term1);
  83. return (term + term1 - aa);
  84. }
  85. static double seor(int order, double qq, double aa, int nterms)
  86. {
  87. double term, term1;
  88. int ii, n1;
  89. term = -1.0*qq;
  90. n1 = (int)((float)order/2.0 - 0.5);
  91. for (ii=0; ii<n1; ii++)
  92. term = qq*qq/(aa - (2.0*ii + 1.0)*(2.0*ii + 1.0) - term);
  93. term += order*order;
  94. term1 = 0.0;
  95. for (ii=0; ii<nterms; ii++)
  96. term1 = qq*qq/
  97. (aa - (order + 2.0*(nterms - ii))*(order + 2.0*(nterms - ii)) - term1);
  98. return (term + term1 - aa);
  99. }
  100. /*----------------------------------------------------------------------------
  101. * Asymptotic and approximation routines for the characteristic value.
  102. *
  103. * Adapted from F.A. Alhargan's paper,
  104. * "Algorithms for the Computation of All Mathieu Functions of Integer
  105. * Orders," ACM Transactions on Mathematical Software, Vol. 26, No. 3,
  106. * September 2000, pp. 390-407.
  107. *--------------------------------------------------------------------------*/
  108. static double asymptotic(int order, double qq)
  109. {
  110. double asymp;
  111. double nn, n2, n4, n6;
  112. double hh, ah, ah2, ah3, ah4, ah5;
  113. /* Set up temporary variables to simplify the readability. */
  114. nn = 2*order + 1;
  115. n2 = nn*nn;
  116. n4 = n2*n2;
  117. n6 = n4*n2;
  118. hh = 2*sqrt(qq);
  119. ah = 16*hh;
  120. ah2 = ah*ah;
  121. ah3 = ah2*ah;
  122. ah4 = ah3*ah;
  123. ah5 = ah4*ah;
  124. /* Equation 38, p. 397 of Alhargan's paper. */
  125. asymp = -2*qq + nn*hh - 0.125*(n2 + 1);
  126. asymp -= 0.25*nn*( n2 + 3)/ah;
  127. asymp -= 0.25* ( 5*n4 + 34*n2 + 9)/ah2;
  128. asymp -= 0.25*nn*( 33*n4 + 410*n2 + 405)/ah3;
  129. asymp -= ( 63*n6 + 1260*n4 + 2943*n2 + 486)/ah4;
  130. asymp -= nn*(527*n6 + 15617*n4 + 69001*n2 + 41607)/ah5;
  131. return asymp;
  132. }
  133. /* Solve the cubic x^3 + c2*x^2 + c1*x + c0 = 0 */
  134. static double solve_cubic(double c2, double c1, double c0)
  135. {
  136. double qq, rr, ww, ss, tt;
  137. qq = (3*c1 - c2*c2)/9;
  138. rr = (9*c2*c1 - 27*c0 - 2*c2*c2*c2)/54;
  139. ww = qq*qq*qq + rr*rr;
  140. if (ww >= 0)
  141. {
  142. double t1 = rr + sqrt(ww);
  143. ss = fabs(t1)/t1*pow(fabs(t1), 1/3.);
  144. t1 = rr - sqrt(ww);
  145. tt = fabs(t1)/t1*pow(fabs(t1), 1/3.);
  146. }
  147. else
  148. {
  149. double theta = acos(rr/sqrt(-qq*qq*qq));
  150. ss = 2*sqrt(-qq)*cos((theta + 4*M_PI)/3.);
  151. tt = 0.0;
  152. }
  153. return (ss + tt - c2/3);
  154. }
  155. /* Compute an initial approximation for the characteristic value. */
  156. static double approx_c(int order, double qq)
  157. {
  158. double approx;
  159. double c0, c1, c2;
  160. if (order < 0)
  161. {
  162. GSL_ERROR_VAL("Undefined order for Mathieu function", GSL_EINVAL, 0.0);
  163. }
  164. switch (order)
  165. {
  166. case 0:
  167. if (qq <= 4)
  168. return (2 - sqrt(4 + 2*qq*qq)); /* Eqn. 31 */
  169. else
  170. return asymptotic(order, qq);
  171. break;
  172. case 1:
  173. if (qq <= 4)
  174. return (5 + 0.5*(qq - sqrt(5*qq*qq - 16*qq + 64))); /* Eqn. 32 */
  175. else
  176. return asymptotic(order, qq);
  177. break;
  178. case 2:
  179. if (qq <= 3)
  180. {
  181. c2 = -8.0; /* Eqn. 33 */
  182. c1 = -48 - 3*qq*qq;
  183. c0 = 20*qq*qq;
  184. }
  185. else
  186. return asymptotic(order, qq);
  187. break;
  188. case 3:
  189. if (qq <= 6.25)
  190. {
  191. c2 = -qq - 8; /* Eqn. 34 */
  192. c1 = 16*qq - 128 - 2*qq*qq;
  193. c0 = qq*qq*(qq + 8);
  194. }
  195. else
  196. return asymptotic(order, qq);
  197. break;
  198. default:
  199. if (order < 70)
  200. {
  201. if (1.7*order > 2*sqrt(qq))
  202. {
  203. /* Eqn. 30 */
  204. double n2 = (double)(order*order);
  205. double n22 = (double)((n2 - 1)*(n2 - 1));
  206. double q2 = qq*qq;
  207. double q4 = q2*q2;
  208. approx = n2 + 0.5*q2/(n2 - 1);
  209. approx += (5*n2 + 7)*q4/(32*n22*(n2 - 1)*(n2 - 4));
  210. approx += (9*n2*n2 + 58*n2 + 29)*q4*q2/
  211. (64*n22*n22*(n2 - 1)*(n2 - 4)*(n2 - 9));
  212. if (1.4*order < 2*sqrt(qq))
  213. {
  214. approx += asymptotic(order, qq);
  215. approx *= 0.5;
  216. }
  217. }
  218. else
  219. approx = asymptotic(order, qq);
  220. return approx;
  221. }
  222. else
  223. return order*order;
  224. }
  225. /* Solve the cubic x^3 + c2*x^2 + c1*x + c0 = 0 */
  226. approx = solve_cubic(c2, c1, c0);
  227. if ( approx < 0 && sqrt(qq) > 0.1*order )
  228. return asymptotic(order-1, qq);
  229. else
  230. return (order*order + fabs(approx));
  231. }
  232. static double approx_s(int order, double qq)
  233. {
  234. double approx;
  235. double c0, c1, c2;
  236. if (order < 1)
  237. {
  238. GSL_ERROR_VAL("Undefined order for Mathieu function", GSL_EINVAL, 0.0);
  239. }
  240. switch (order)
  241. {
  242. case 1:
  243. if (qq <= 4)
  244. return (5 - 0.5*(qq + sqrt(5*qq*qq + 16*qq + 64))); /* Eqn. 35 */
  245. else
  246. return asymptotic(order-1, qq);
  247. break;
  248. case 2:
  249. if (qq <= 5)
  250. return (10 - sqrt(36 + qq*qq)); /* Eqn. 36 */
  251. else
  252. return asymptotic(order-1, qq);
  253. break;
  254. case 3:
  255. if (qq <= 6.25)
  256. {
  257. c2 = qq - 8; /* Eqn. 37 */
  258. c1 = -128 - 16*qq - 2*qq*qq;
  259. c0 = qq*qq*(8 - qq);
  260. }
  261. else
  262. return asymptotic(order-1, qq);
  263. break;
  264. default:
  265. if (order < 70)
  266. {
  267. if (1.7*order > 2*sqrt(qq))
  268. {
  269. /* Eqn. 30 */
  270. double n2 = (double)(order*order);
  271. double n22 = (double)((n2 - 1)*(n2 - 1));
  272. double q2 = qq*qq;
  273. double q4 = q2*q2;
  274. approx = n2 + 0.5*q2/(n2 - 1);
  275. approx += (5*n2 + 7)*q4/(32*n22*(n2 - 1)*(n2 - 4));
  276. approx += (9*n2*n2 + 58*n2 + 29)*q4*q2/
  277. (64*n22*n22*(n2 - 1)*(n2 - 4)*(n2 - 9));
  278. if (1.4*order < 2*sqrt(qq))
  279. {
  280. approx += asymptotic(order-1, qq);
  281. approx *= 0.5;
  282. }
  283. }
  284. else
  285. approx = asymptotic(order-1, qq);
  286. return approx;
  287. }
  288. else
  289. return order*order;
  290. }
  291. /* Solve the cubic x^3 + c2*x^2 + c1*x + c0 = 0 */
  292. approx = solve_cubic(c2, c1, c0);
  293. if ( approx < 0 && sqrt(qq) > 0.1*order )
  294. return asymptotic(order-1, qq);
  295. else
  296. return (order*order + fabs(approx));
  297. }
  298. int gsl_sf_mathieu_a(int order, double qq, gsl_sf_result *result)
  299. {
  300. int even_odd, nterms = 50, ii, counter = 0, maxcount = 200;
  301. double a1, a2, fa, fa1, dela, aa_orig, da = 0.025, aa;
  302. even_odd = 0;
  303. if (order % 2 != 0)
  304. even_odd = 1;
  305. /* If the argument is 0, then the coefficient is simply the square of
  306. the order. */
  307. if (qq == 0)
  308. {
  309. result->val = order*order;
  310. result->err = 0.0;
  311. return GSL_SUCCESS;
  312. }
  313. /* Use symmetry characteristics of the functions to handle cases with
  314. negative order and/or argument q. See Abramowitz & Stegun, 20.8.3. */
  315. if (order < 0)
  316. order *= -1;
  317. if (qq < 0.0)
  318. {
  319. if (even_odd == 0)
  320. return gsl_sf_mathieu_a(order, -qq, result);
  321. else
  322. return gsl_sf_mathieu_b(order, -qq, result);
  323. }
  324. /* Compute an initial approximation for the characteristic value. */
  325. aa = approx_c(order, qq);
  326. /* Save the original approximation for later comparison. */
  327. aa_orig = aa;
  328. /* Loop as long as the final value is not near the approximate value
  329. (with a max limit to avoid potential infinite loop). */
  330. while (counter < maxcount)
  331. {
  332. a1 = aa + 0.001;
  333. ii = 0;
  334. if (even_odd == 0)
  335. fa1 = ceer(order, qq, a1, nterms);
  336. else
  337. fa1 = ceor(order, qq, a1, nterms);
  338. for (;;)
  339. {
  340. if (even_odd == 0)
  341. fa = ceer(order, qq, aa, nterms);
  342. else
  343. fa = ceor(order, qq, aa, nterms);
  344. a2 = a1;
  345. a1 = aa;
  346. if (fa == fa1)
  347. {
  348. result->err = GSL_DBL_EPSILON;
  349. break;
  350. }
  351. aa -= (aa - a2)/(fa - fa1)*fa;
  352. dela = fabs(aa - a2);
  353. if (dela < GSL_DBL_EPSILON)
  354. {
  355. result->err = GSL_DBL_EPSILON;
  356. break;
  357. }
  358. if (ii > 20)
  359. {
  360. result->err = dela;
  361. break;
  362. }
  363. fa1 = fa;
  364. ii++;
  365. }
  366. /* If the solution found is not near the original approximation,
  367. tweak the approximate value, and try again. */
  368. if (fabs(aa - aa_orig) > (3 + 0.01*order*fabs(aa_orig)))
  369. {
  370. counter++;
  371. if (counter == maxcount)
  372. {
  373. result->err = fabs(aa - aa_orig);
  374. break;
  375. }
  376. if (aa > aa_orig)
  377. aa = aa_orig - da*counter;
  378. else
  379. aa = aa_orig + da*counter;
  380. continue;
  381. }
  382. else
  383. break;
  384. }
  385. result->val = aa;
  386. /* If we went through the maximum number of retries and still didn't
  387. find the solution, let us know. */
  388. if (counter == maxcount)
  389. {
  390. GSL_ERROR("Wrong characteristic Mathieu value", GSL_EFAILED);
  391. }
  392. return GSL_SUCCESS;
  393. }
  394. int gsl_sf_mathieu_b(int order, double qq, gsl_sf_result *result)
  395. {
  396. int even_odd, nterms = 50, ii, counter = 0, maxcount = 200;
  397. double a1, a2, fa, fa1, dela, aa_orig, da = 0.025, aa;
  398. even_odd = 0;
  399. if (order % 2 != 0)
  400. even_odd = 1;
  401. /* The order cannot be 0. */
  402. if (order == 0)
  403. {
  404. GSL_ERROR("Characteristic value undefined for order 0", GSL_EFAILED);
  405. }
  406. /* If the argument is 0, then the coefficient is simply the square of
  407. the order. */
  408. if (qq == 0)
  409. {
  410. result->val = order*order;
  411. result->err = 0.0;
  412. return GSL_SUCCESS;
  413. }
  414. /* Use symmetry characteristics of the functions to handle cases with
  415. negative order and/or argument q. See Abramowitz & Stegun, 20.8.3. */
  416. if (order < 0)
  417. order *= -1;
  418. if (qq < 0.0)
  419. {
  420. if (even_odd == 0)
  421. return gsl_sf_mathieu_b(order, -qq, result);
  422. else
  423. return gsl_sf_mathieu_a(order, -qq, result);
  424. }
  425. /* Compute an initial approximation for the characteristic value. */
  426. aa = approx_s(order, qq);
  427. /* Save the original approximation for later comparison. */
  428. aa_orig = aa;
  429. /* Loop as long as the final value is not near the approximate value
  430. (with a max limit to avoid potential infinite loop). */
  431. while (counter < maxcount)
  432. {
  433. a1 = aa + 0.001;
  434. ii = 0;
  435. if (even_odd == 0)
  436. fa1 = seer(order, qq, a1, nterms);
  437. else
  438. fa1 = seor(order, qq, a1, nterms);
  439. for (;;)
  440. {
  441. if (even_odd == 0)
  442. fa = seer(order, qq, aa, nterms);
  443. else
  444. fa = seor(order, qq, aa, nterms);
  445. a2 = a1;
  446. a1 = aa;
  447. if (fa == fa1)
  448. {
  449. result->err = GSL_DBL_EPSILON;
  450. break;
  451. }
  452. aa -= (aa - a2)/(fa - fa1)*fa;
  453. dela = fabs(aa - a2);
  454. if (dela < 1e-18)
  455. {
  456. result->err = GSL_DBL_EPSILON;
  457. break;
  458. }
  459. if (ii > 20)
  460. {
  461. result->err = dela;
  462. break;
  463. }
  464. fa1 = fa;
  465. ii++;
  466. }
  467. /* If the solution found is not near the original approximation,
  468. tweak the approximate value, and try again. */
  469. if (fabs(aa - aa_orig) > (3 + 0.01*order*fabs(aa_orig)))
  470. {
  471. counter++;
  472. if (counter == maxcount)
  473. {
  474. result->err = fabs(aa - aa_orig);
  475. break;
  476. }
  477. if (aa > aa_orig)
  478. aa = aa_orig - da*counter;
  479. else
  480. aa = aa_orig + da*counter;
  481. continue;
  482. }
  483. else
  484. break;
  485. }
  486. result->val = aa;
  487. /* If we went through the maximum number of retries and still didn't
  488. find the solution, let us know. */
  489. if (counter == maxcount)
  490. {
  491. GSL_ERROR("Wrong characteristic Mathieu value", GSL_EFAILED);
  492. }
  493. return GSL_SUCCESS;
  494. }
  495. /* Eigenvalue solutions for characteristic values below. */
  496. /* figi.c converted from EISPACK Fortran FIGI.F.
  497. *
  498. * given a nonsymmetric tridiagonal matrix such that the products
  499. * of corresponding pairs of off-diagonal elements are all
  500. * non-negative, this subroutine reduces it to a symmetric
  501. * tridiagonal matrix with the same eigenvalues. if, further,
  502. * a zero product only occurs when both factors are zero,
  503. * the reduced matrix is similar to the original matrix.
  504. *
  505. * on input
  506. *
  507. * n is the order of the matrix.
  508. *
  509. * t contains the input matrix. its subdiagonal is
  510. * stored in the last n-1 positions of the first column,
  511. * its diagonal in the n positions of the second column,
  512. * and its superdiagonal in the first n-1 positions of
  513. * the third column. t(1,1) and t(n,3) are arbitrary.
  514. *
  515. * on output
  516. *
  517. * t is unaltered.
  518. *
  519. * d contains the diagonal elements of the symmetric matrix.
  520. *
  521. * e contains the subdiagonal elements of the symmetric
  522. * matrix in its last n-1 positions. e(1) is not set.
  523. *
  524. * e2 contains the squares of the corresponding elements of e.
  525. * e2 may coincide with e if the squares are not needed.
  526. *
  527. * ierr is set to
  528. * zero for normal return,
  529. * n+i if t(i,1)*t(i-1,3) is negative,
  530. * -(3*n+i) if t(i,1)*t(i-1,3) is zero with one factor
  531. * non-zero. in this case, the eigenvectors of
  532. * the symmetric matrix are not simply related
  533. * to those of t and should not be sought.
  534. *
  535. * questions and comments should be directed to burton s. garbow,
  536. * mathematics and computer science div, argonne national laboratory
  537. *
  538. * this version dated august 1983.
  539. */
  540. static int figi(int nn, double *tt, double *dd, double *ee,
  541. double *e2)
  542. {
  543. int ii;
  544. for (ii=0; ii<nn; ii++)
  545. {
  546. if (ii != 0)
  547. {
  548. e2[ii] = tt[3*ii]*tt[3*(ii-1)+2];
  549. if (e2[ii] < 0.0)
  550. {
  551. /* set error -- product of some pair of off-diagonal
  552. elements is negative */
  553. return (nn + ii);
  554. }
  555. if (e2[ii] == 0.0 && (tt[3*ii] != 0.0 || tt[3*(ii-1)+2] != 0.0))
  556. {
  557. /* set error -- product of some pair of off-diagonal
  558. elements is zero with one member non-zero */
  559. return (-1*(3*nn + ii));
  560. }
  561. ee[ii] = sqrt(e2[ii]);
  562. }
  563. dd[ii] = tt[3*ii+1];
  564. }
  565. return 0;
  566. }
  567. int gsl_sf_mathieu_a_array(int order_min, int order_max, double qq, gsl_sf_mathieu_workspace *work, double result_array[])
  568. {
  569. unsigned int even_order = work->even_order, odd_order = work->odd_order,
  570. extra_values = work->extra_values, ii, jj;
  571. int status;
  572. double *tt = work->tt, *dd = work->dd, *ee = work->ee, *e2 = work->e2,
  573. *zz = work->zz, *aa = work->aa;
  574. gsl_matrix_view mat, evec;
  575. gsl_vector_view eval;
  576. gsl_eigen_symmv_workspace *wmat = work->wmat;
  577. if (order_max > work->size || order_max <= order_min || order_min < 0)
  578. {
  579. GSL_ERROR ("invalid range [order_min,order_max]", GSL_EINVAL);
  580. }
  581. /* Convert the nonsymmetric tridiagonal matrix to a symmetric tridiagonal
  582. form. */
  583. tt[0] = 0.0;
  584. tt[1] = 0.0;
  585. tt[2] = qq;
  586. for (ii=1; ii<even_order-1; ii++)
  587. {
  588. tt[3*ii] = qq;
  589. tt[3*ii+1] = 4*ii*ii;
  590. tt[3*ii+2] = qq;
  591. }
  592. tt[3*even_order-3] = qq;
  593. tt[3*even_order-2] = 4*(even_order - 1)*(even_order - 1);
  594. tt[3*even_order-1] = 0.0;
  595. tt[3] *= 2;
  596. status = figi((signed int)even_order, tt, dd, ee, e2);
  597. if (status)
  598. {
  599. GSL_ERROR("Internal error in tridiagonal Mathieu matrix", GSL_EFAILED);
  600. }
  601. /* Fill the period \pi matrix. */
  602. for (ii=0; ii<even_order*even_order; ii++)
  603. zz[ii] = 0.0;
  604. zz[0] = dd[0];
  605. zz[1] = ee[1];
  606. for (ii=1; ii<even_order-1; ii++)
  607. {
  608. zz[ii*even_order+ii-1] = ee[ii];
  609. zz[ii*even_order+ii] = dd[ii];
  610. zz[ii*even_order+ii+1] = ee[ii+1];
  611. }
  612. zz[even_order*(even_order-1)+even_order-2] = ee[even_order-1];
  613. zz[even_order*even_order-1] = dd[even_order-1];
  614. /* Compute (and sort) the eigenvalues of the matrix. */
  615. mat = gsl_matrix_view_array(zz, even_order, even_order);
  616. eval = gsl_vector_subvector(work->eval, 0, even_order);
  617. evec = gsl_matrix_submatrix(work->evec, 0, 0, even_order, even_order);
  618. gsl_eigen_symmv(&mat.matrix, &eval.vector, &evec.matrix, wmat);
  619. gsl_eigen_symmv_sort(&eval.vector, &evec.matrix, GSL_EIGEN_SORT_VAL_ASC);
  620. for (ii=0; ii<even_order-extra_values; ii++)
  621. aa[2*ii] = gsl_vector_get(&eval.vector, ii);
  622. /* Fill the period 2\pi matrix. */
  623. for (ii=0; ii<odd_order*odd_order; ii++)
  624. zz[ii] = 0.0;
  625. for (ii=0; ii<odd_order; ii++)
  626. for (jj=0; jj<odd_order; jj++)
  627. {
  628. if (ii == jj)
  629. zz[ii*odd_order+jj] = (2*ii + 1)*(2*ii + 1);
  630. else if (ii == jj + 1 || ii + 1 == jj)
  631. zz[ii*odd_order+jj] = qq;
  632. }
  633. zz[0] += qq;
  634. /* Compute (and sort) the eigenvalues of the matrix. */
  635. mat = gsl_matrix_view_array(zz, odd_order, odd_order);
  636. eval = gsl_vector_subvector(work->eval, 0, odd_order);
  637. evec = gsl_matrix_submatrix(work->evec, 0, 0, odd_order, odd_order);
  638. gsl_eigen_symmv(&mat.matrix, &eval.vector, &evec.matrix, wmat);
  639. gsl_eigen_symmv_sort(&eval.vector, &evec.matrix, GSL_EIGEN_SORT_VAL_ASC);
  640. for (ii=0; ii<odd_order-extra_values; ii++)
  641. aa[2*ii+1] = gsl_vector_get(&eval.vector, ii);
  642. for (ii = order_min ; ii <= order_max ; ii++)
  643. {
  644. result_array[ii - order_min] = aa[ii];
  645. }
  646. return GSL_SUCCESS;
  647. }
  648. int gsl_sf_mathieu_b_array(int order_min, int order_max, double qq, gsl_sf_mathieu_workspace *work, double result_array[])
  649. {
  650. unsigned int even_order = work->even_order-1, odd_order = work->odd_order,
  651. extra_values = work->extra_values, ii, jj;
  652. double *zz = work->zz, *bb = work->bb;
  653. gsl_matrix_view mat, evec;
  654. gsl_vector_view eval;
  655. gsl_eigen_symmv_workspace *wmat = work->wmat;
  656. if (order_max > work->size || order_max <= order_min || order_min < 0)
  657. {
  658. GSL_ERROR ("invalid range [order_min,order_max]", GSL_EINVAL);
  659. }
  660. /* Fill the period \pi matrix. */
  661. for (ii=0; ii<even_order*even_order; ii++)
  662. zz[ii] = 0.0;
  663. for (ii=0; ii<even_order; ii++)
  664. for (jj=0; jj<even_order; jj++)
  665. {
  666. if (ii == jj)
  667. zz[ii*even_order+jj] = 4*(ii + 1)*(ii + 1);
  668. else if (ii == jj + 1 || ii + 1 == jj)
  669. zz[ii*even_order+jj] = qq;
  670. }
  671. /* Compute (and sort) the eigenvalues of the matrix. */
  672. mat = gsl_matrix_view_array(zz, even_order, even_order);
  673. eval = gsl_vector_subvector(work->eval, 0, even_order);
  674. evec = gsl_matrix_submatrix(work->evec, 0, 0, even_order, even_order);
  675. gsl_eigen_symmv(&mat.matrix, &eval.vector, &evec.matrix, wmat);
  676. gsl_eigen_symmv_sort(&eval.vector, &evec.matrix, GSL_EIGEN_SORT_VAL_ASC);
  677. bb[0] = 0.0;
  678. for (ii=0; ii<even_order-extra_values; ii++)
  679. bb[2*(ii+1)] = gsl_vector_get(&eval.vector, ii);
  680. /* Fill the period 2\pi matrix. */
  681. for (ii=0; ii<odd_order*odd_order; ii++)
  682. zz[ii] = 0.0;
  683. for (ii=0; ii<odd_order; ii++)
  684. for (jj=0; jj<odd_order; jj++)
  685. {
  686. if (ii == jj)
  687. zz[ii*odd_order+jj] = (2*ii + 1)*(2*ii + 1);
  688. else if (ii == jj + 1 || ii + 1 == jj)
  689. zz[ii*odd_order+jj] = qq;
  690. }
  691. zz[0] -= qq;
  692. /* Compute (and sort) the eigenvalues of the matrix. */
  693. mat = gsl_matrix_view_array(zz, odd_order, odd_order);
  694. eval = gsl_vector_subvector(work->eval, 0, odd_order);
  695. evec = gsl_matrix_submatrix(work->evec, 0, 0, odd_order, odd_order);
  696. gsl_eigen_symmv(&mat.matrix, &eval.vector, &evec.matrix, wmat);
  697. gsl_eigen_symmv_sort(&eval.vector, &evec.matrix, GSL_EIGEN_SORT_VAL_ASC);
  698. for (ii=0; ii<odd_order-extra_values; ii++)
  699. bb[2*ii+1] = gsl_vector_get(&eval.vector, ii);
  700. for (ii = order_min ; ii <= order_max ; ii++)
  701. {
  702. result_array[ii - order_min] = bb[ii];
  703. }
  704. return GSL_SUCCESS;
  705. }