summaryrefslogtreecommitdiff
path: root/gsl-1.9/specfunc/mathieu_charv.c
blob: bf284135b16f900110060b719f6884ea5f092e50 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
/* specfunc/mathieu_charv.c
 * 
 * Copyright (C) 2002 Lowell Johnson
 * 
 * This program is free software; you can redistribute it and/or modify
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation; either version 2 of the License, or (at
 * your option) any later version.
 * 
 * This program is distributed in the hope that it will be useful, but
 * WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
 * General Public License for more details.
 * 
 * You should have received a copy of the GNU General Public License
 * along with this program; if not, write to the Free Software
 * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
 */

/* Author:  L. Johnson */

#include <config.h>
#include <stdlib.h>
#include <stdio.h>
#include <math.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_eigen.h>
#include <gsl/gsl_errno.h>
#include <gsl/gsl_sf_mathieu.h>


/* prototypes */
static double solve_cubic(double c2, double c1, double c0);


static double ceer(int order, double qq, double aa, int nterms)
{

  double term, term1;
  int ii, n1;
  

  if (order == 0)
      term = 0.0;
  else
  {      
      term = 2.0*qq*qq/aa;

      if (order != 2)
      {
          n1 = order/2 - 1;

          for (ii=0; ii<n1; ii++)
              term = qq*qq/(aa - 4.0*(ii+1)*(ii+1) - term);
      }
  }
  
  term += order*order;

  term1 = 0.0;

  for (ii=0; ii<nterms; ii++)
      term1 = qq*qq/
        (aa - (order + 2.0*(nterms - ii))*(order + 2.0*(nterms - ii)) - term1);

  if (order == 0)
      term1 *= 2.0;
  
  return (term + term1 - aa);
}


static double ceor(int order, double qq, double aa, int nterms)
{
  double term, term1;
  int ii, n1;

  term = qq;
  n1 = (int)((float)order/2.0 - 0.5);

  for (ii=0; ii<n1; ii++)
      term = qq*qq/(aa - (2.0*ii + 1.0)*(2.0*ii + 1.0) - term);
  term += order*order;

  term1 = 0.0;
  for (ii=0; ii<nterms; ii++)
      term1 = qq*qq/
        (aa - (order + 2.0*(nterms - ii))*(order + 2.0*(nterms - ii)) - term1);

  return (term + term1 - aa);
}


static double seer(int order, double qq, double aa, int nterms)
{
  double term, term1;
  int ii, n1;

  term = 0.0;
  n1 = order/2 - 1;

  for (ii=0; ii<n1; ii++)
      term = qq*qq/(aa - 4.0*(ii + 1)*(ii + 1) - term);
  term += order*order;

  term1 = 0.0;
  for (ii=0; ii<nterms; ii++)
      term1 = qq*qq/
        (aa - (order + 2.0*(nterms - ii))*(order + 2.0*(nterms - ii)) - term1);

  return (term + term1 - aa);
}


static double seor(int order, double qq, double aa, int nterms)
{
  double term, term1;
  int ii, n1;


  term = -1.0*qq;
  n1 = (int)((float)order/2.0 - 0.5);
  for (ii=0; ii<n1; ii++)
      term = qq*qq/(aa - (2.0*ii + 1.0)*(2.0*ii + 1.0) - term);
  term += order*order;

  term1 = 0.0;
  for (ii=0; ii<nterms; ii++)
      term1 = qq*qq/
        (aa - (order + 2.0*(nterms - ii))*(order + 2.0*(nterms - ii)) - term1);

  return (term + term1 - aa);
}


/*----------------------------------------------------------------------------
 * Asymptotic and approximation routines for the characteristic value.
 *
 * Adapted from F.A. Alhargan's paper,
 * "Algorithms for the Computation of All Mathieu Functions of Integer
 * Orders," ACM Transactions on Mathematical Software, Vol. 26, No. 3,
 * September 2000, pp. 390-407.
 *--------------------------------------------------------------------------*/
static double asymptotic(int order, double qq)
{
  double asymp;
  double nn, n2, n4, n6;
  double hh, ah, ah2, ah3, ah4, ah5;


  /* Set up temporary variables to simplify the readability. */
  nn = 2*order + 1;
  n2 = nn*nn;
  n4 = n2*n2;
  n6 = n4*n2;
  
  hh = 2*sqrt(qq);
  ah = 16*hh;
  ah2 = ah*ah;
  ah3 = ah2*ah;
  ah4 = ah3*ah;
  ah5 = ah4*ah;

  /* Equation 38, p. 397 of Alhargan's paper. */
  asymp = -2*qq + nn*hh - 0.125*(n2 + 1);
  asymp -= 0.25*nn*(                          n2 +     3)/ah;
  asymp -= 0.25*   (             5*n4 +    34*n2 +     9)/ah2;
  asymp -= 0.25*nn*(            33*n4 +   410*n2 +   405)/ah3;
  asymp -=         ( 63*n6 +  1260*n4 +  2943*n2 +   486)/ah4;
  asymp -=      nn*(527*n6 + 15617*n4 + 69001*n2 + 41607)/ah5;

  return asymp;
}


/* Solve the cubic x^3 + c2*x^2 + c1*x + c0 = 0 */
static double solve_cubic(double c2, double c1, double c0)
{
  double qq, rr, ww, ss, tt;

  
  qq = (3*c1 - c2*c2)/9;
  rr = (9*c2*c1 - 27*c0 - 2*c2*c2*c2)/54;
  ww = qq*qq*qq + rr*rr;
  
  if (ww >= 0)
  {
      double t1 = rr + sqrt(ww);
      ss = fabs(t1)/t1*pow(fabs(t1), 1/3.);
      t1 = rr - sqrt(ww);
      tt = fabs(t1)/t1*pow(fabs(t1), 1/3.);
  }
  else
  {
      double theta = acos(rr/sqrt(-qq*qq*qq));
      ss = 2*sqrt(-qq)*cos((theta + 4*M_PI)/3.);
      tt = 0.0;
  }
  
  return (ss + tt - c2/3);
}


/* Compute an initial approximation for the characteristic value. */
static double approx_c(int order, double qq)
{
  double approx;
  double c0, c1, c2;


  if (order < 0)
  {
    GSL_ERROR_VAL("Undefined order for Mathieu function", GSL_EINVAL, 0.0);
  }
  
  switch (order)
  {
      case 0:
          if (qq <= 4)
              return (2 - sqrt(4 + 2*qq*qq)); /* Eqn. 31 */
          else
              return asymptotic(order, qq);
          break;

      case 1:
          if (qq <= 4)
              return (5 + 0.5*(qq - sqrt(5*qq*qq - 16*qq + 64))); /* Eqn. 32 */
          else
              return asymptotic(order, qq);
          break;

      case 2:
          if (qq <= 3)
          {
              c2 = -8.0;  /* Eqn. 33 */
              c1 = -48 - 3*qq*qq;
              c0 = 20*qq*qq;
          }
          else
              return asymptotic(order, qq);
          break;

      case 3:
          if (qq <= 6.25)
          {
              c2 = -qq - 8;  /* Eqn. 34 */
              c1 = 16*qq - 128 - 2*qq*qq;
              c0 = qq*qq*(qq + 8);
          }
          else
              return asymptotic(order, qq);
          break;

      default:
          if (order < 70)
          {
              if (1.7*order > 2*sqrt(qq))
              {
                  /* Eqn. 30 */
                  double n2 = (double)(order*order);
                  double n22 = (double)((n2 - 1)*(n2 - 1));
                  double q2 = qq*qq;
                  double q4 = q2*q2;
                  approx = n2 + 0.5*q2/(n2 - 1);
                  approx += (5*n2 + 7)*q4/(32*n22*(n2 - 1)*(n2 - 4));
                  approx += (9*n2*n2 + 58*n2 + 29)*q4*q2/
                      (64*n22*n22*(n2 - 1)*(n2 - 4)*(n2 - 9));
                  if (1.4*order < 2*sqrt(qq))
                  {
                      approx += asymptotic(order, qq);
                      approx *= 0.5;
                  }
              }
              else
                  approx = asymptotic(order, qq);

              return approx;
          }
          else
              return order*order;
  }

  /* Solve the cubic x^3 + c2*x^2 + c1*x + c0 = 0 */
  approx = solve_cubic(c2, c1, c0);
      
  if ( approx < 0 && sqrt(qq) > 0.1*order )
      return asymptotic(order-1, qq);
  else
      return (order*order + fabs(approx));
}

  
static double approx_s(int order, double qq)
{
  double approx;
  double c0, c1, c2;

  
  if (order < 1)
  {
    GSL_ERROR_VAL("Undefined order for Mathieu function", GSL_EINVAL, 0.0);
  }
  
  switch (order)
  {
      case 1:
          if (qq <= 4)
              return (5 - 0.5*(qq + sqrt(5*qq*qq + 16*qq + 64))); /* Eqn. 35 */
          else
              return asymptotic(order-1, qq);
          break;

      case 2:
          if (qq <= 5)
              return (10 - sqrt(36 + qq*qq)); /* Eqn. 36 */
          else
              return asymptotic(order-1, qq);
          break;

      case 3:
          if (qq <= 6.25)
          {
              c2 = qq - 8; /* Eqn. 37 */
              c1 = -128 - 16*qq - 2*qq*qq;
              c0 = qq*qq*(8 - qq);
          }
          else
              return asymptotic(order-1, qq);
          break;

      default:
          if (order < 70)
          {
              if (1.7*order > 2*sqrt(qq))
              {
                  /* Eqn. 30 */
                  double n2 = (double)(order*order);
                  double n22 = (double)((n2 - 1)*(n2 - 1));
                  double q2 = qq*qq;
                  double q4 = q2*q2;
                  approx = n2 + 0.5*q2/(n2 - 1);
                  approx += (5*n2 + 7)*q4/(32*n22*(n2 - 1)*(n2 - 4));
                  approx += (9*n2*n2 + 58*n2 + 29)*q4*q2/
                      (64*n22*n22*(n2 - 1)*(n2 - 4)*(n2 - 9));
                  if (1.4*order < 2*sqrt(qq))
                  {
                      approx += asymptotic(order-1, qq);
                      approx *= 0.5;
                  }
              }
              else
                  approx = asymptotic(order-1, qq);

              return approx;
          }
          else
              return order*order;
  }

  /* Solve the cubic x^3 + c2*x^2 + c1*x + c0 = 0 */
  approx = solve_cubic(c2, c1, c0);
      
  if ( approx < 0 && sqrt(qq) > 0.1*order )
      return asymptotic(order-1, qq);
  else
      return (order*order + fabs(approx));
}


int gsl_sf_mathieu_a(int order, double qq, gsl_sf_result *result)
{
  int even_odd, nterms = 50, ii, counter = 0, maxcount = 200;
  double a1, a2, fa, fa1, dela, aa_orig, da = 0.025, aa;


  even_odd = 0;
  if (order % 2 != 0)
      even_odd = 1;

  /* If the argument is 0, then the coefficient is simply the square of
     the order. */
  if (qq == 0)
  {
      result->val = order*order;
      result->err = 0.0;
      return GSL_SUCCESS;
  }

  /* Use symmetry characteristics of the functions to handle cases with
     negative order and/or argument q.  See Abramowitz & Stegun, 20.8.3. */
  if (order < 0)
      order *= -1;
  if (qq < 0.0)
  {
      if (even_odd == 0)
          return gsl_sf_mathieu_a(order, -qq, result);
      else
          return gsl_sf_mathieu_b(order, -qq, result);
  }
  
  /* Compute an initial approximation for the characteristic value. */
  aa = approx_c(order, qq);

  /* Save the original approximation for later comparison. */
  aa_orig = aa;
  
  /* Loop as long as the final value is not near the approximate value
     (with a max limit to avoid potential infinite loop). */
  while (counter < maxcount)
  {
      a1 = aa + 0.001;
      ii = 0;
      if (even_odd == 0)
          fa1 = ceer(order, qq, a1, nterms);
      else
          fa1 = ceor(order, qq, a1, nterms);

      for (;;)
      {
          if (even_odd == 0)
              fa = ceer(order, qq, aa, nterms);
          else
              fa = ceor(order, qq, aa, nterms);
      
          a2 = a1;
          a1 = aa;

          if (fa == fa1)
          {
              result->err = GSL_DBL_EPSILON;
              break;
          }
          aa -= (aa - a2)/(fa - fa1)*fa;
          dela = fabs(aa - a2);
          if (dela < GSL_DBL_EPSILON)
          {
              result->err = GSL_DBL_EPSILON;
              break;
          }
          if (ii > 20)
          {
              result->err = dela;
              break;
          }
          fa1 = fa;
          ii++;
      }

      /* If the solution found is not near the original approximation,
         tweak the approximate value, and try again. */
      if (fabs(aa - aa_orig) > (3 + 0.01*order*fabs(aa_orig)))
      {
          counter++;
          if (counter == maxcount)
          {
              result->err = fabs(aa - aa_orig);
              break;
          }
          if (aa > aa_orig)
              aa = aa_orig - da*counter;
          else
              aa = aa_orig + da*counter;

          continue;
      }
      else
          break;
  }

  result->val = aa;
      
  /* If we went through the maximum number of retries and still didn't
     find the solution, let us know. */
  if (counter == maxcount)
  {
      GSL_ERROR("Wrong characteristic Mathieu value", GSL_EFAILED);
  }
  
  return GSL_SUCCESS;
}


int gsl_sf_mathieu_b(int order, double qq, gsl_sf_result *result)
{
  int even_odd, nterms = 50, ii, counter = 0, maxcount = 200;
  double a1, a2, fa, fa1, dela, aa_orig, da = 0.025, aa;


  even_odd = 0;
  if (order % 2 != 0)
      even_odd = 1;

  /* The order cannot be 0. */
  if (order == 0)
  {
      GSL_ERROR("Characteristic value undefined for order 0", GSL_EFAILED);
  }

  /* If the argument is 0, then the coefficient is simply the square of
     the order. */
  if (qq == 0)
  {
      result->val = order*order;
      result->err = 0.0;
      return GSL_SUCCESS;
  }

  /* Use symmetry characteristics of the functions to handle cases with
     negative order and/or argument q.  See Abramowitz & Stegun, 20.8.3. */
  if (order < 0)
      order *= -1;
  if (qq < 0.0)
  {
      if (even_odd == 0)
          return gsl_sf_mathieu_b(order, -qq, result);
      else
          return gsl_sf_mathieu_a(order, -qq, result);
  }
  
  /* Compute an initial approximation for the characteristic value. */
  aa = approx_s(order, qq);
  
  /* Save the original approximation for later comparison. */
  aa_orig = aa;
  
  /* Loop as long as the final value is not near the approximate value
     (with a max limit to avoid potential infinite loop). */
  while (counter < maxcount)
  {
      a1 = aa + 0.001;
      ii = 0;
      if (even_odd == 0)
          fa1 = seer(order, qq, a1, nterms);
      else
          fa1 = seor(order, qq, a1, nterms);

      for (;;)
      {
          if (even_odd == 0)
              fa = seer(order, qq, aa, nterms);
          else
              fa = seor(order, qq, aa, nterms);
      
          a2 = a1;
          a1 = aa;

          if (fa == fa1)
          {
              result->err = GSL_DBL_EPSILON;
              break;
          }
          aa -= (aa - a2)/(fa - fa1)*fa;
          dela = fabs(aa - a2);
          if (dela < 1e-18)
          {
              result->err = GSL_DBL_EPSILON;
              break;
          }
          if (ii > 20)
          {
              result->err = dela;
              break;
          }
          fa1 = fa;
          ii++;
      }
      
      /* If the solution found is not near the original approximation,
         tweak the approximate value, and try again. */
      if (fabs(aa - aa_orig) > (3 + 0.01*order*fabs(aa_orig)))
      {
          counter++;
          if (counter == maxcount)
          {
              result->err = fabs(aa - aa_orig);
              break;
          }
          if (aa > aa_orig)
              aa = aa_orig - da*counter;
          else
              aa = aa_orig + da*counter;
          
          continue;
      }
      else
          break;
  }
  
  result->val = aa;
      
  /* If we went through the maximum number of retries and still didn't
     find the solution, let us know. */
  if (counter == maxcount)
  {
      GSL_ERROR("Wrong characteristic Mathieu value", GSL_EFAILED);
  }
  
  return GSL_SUCCESS;
}


/* Eigenvalue solutions for characteristic values below. */


/*  figi.c converted from EISPACK Fortran FIGI.F.
 *
 *   given a nonsymmetric tridiagonal matrix such that the products
 *    of corresponding pairs of off-diagonal elements are all
 *    non-negative, this subroutine reduces it to a symmetric
 *    tridiagonal matrix with the same eigenvalues.  if, further,
 *    a zero product only occurs when both factors are zero,
 *    the reduced matrix is similar to the original matrix.
 *
 *    on input
 *
 *       n is the order of the matrix.
 *
 *       t contains the input matrix.  its subdiagonal is
 *         stored in the last n-1 positions of the first column,
 *         its diagonal in the n positions of the second column,
 *         and its superdiagonal in the first n-1 positions of
 *         the third column.  t(1,1) and t(n,3) are arbitrary.
 *
 *    on output
 *
 *       t is unaltered.
 *
 *       d contains the diagonal elements of the symmetric matrix.
 *
 *       e contains the subdiagonal elements of the symmetric
 *         matrix in its last n-1 positions.  e(1) is not set.
 *
 *       e2 contains the squares of the corresponding elements of e.
 *         e2 may coincide with e if the squares are not needed.
 *
 *       ierr is set to
 *         zero       for normal return,
 *         n+i        if t(i,1)*t(i-1,3) is negative,
 *         -(3*n+i)   if t(i,1)*t(i-1,3) is zero with one factor
 *                    non-zero.  in this case, the eigenvectors of
 *                    the symmetric matrix are not simply related
 *                    to those of  t  and should not be sought.
 *
 *    questions and comments should be directed to burton s. garbow,
 *    mathematics and computer science div, argonne national laboratory
 *
 *    this version dated august 1983.
 */
static int figi(int nn, double *tt, double *dd, double *ee,
                double *e2)
{
  int ii;

  for (ii=0; ii<nn; ii++)
  {
      if (ii != 0)
      {
          e2[ii] = tt[3*ii]*tt[3*(ii-1)+2];

          if (e2[ii] < 0.0)
          {
              /* set error -- product of some pair of off-diagonal
                 elements is negative */
              return (nn + ii);
          }

          if (e2[ii] == 0.0 && (tt[3*ii] != 0.0 || tt[3*(ii-1)+2] != 0.0))
          {
              /* set error -- product of some pair of off-diagonal
                 elements is zero with one member non-zero */
              return (-1*(3*nn + ii));
          }

          ee[ii] = sqrt(e2[ii]);
      }

      dd[ii] = tt[3*ii+1];
  }

  return 0;
}


int gsl_sf_mathieu_a_array(int order_min, int order_max, double qq, gsl_sf_mathieu_workspace *work, double result_array[])
{
  unsigned int even_order = work->even_order, odd_order = work->odd_order,
      extra_values = work->extra_values, ii, jj;
  int status;
  double *tt = work->tt, *dd = work->dd, *ee = work->ee, *e2 = work->e2,
         *zz = work->zz, *aa = work->aa;
  gsl_matrix_view mat, evec;
  gsl_vector_view eval;
  gsl_eigen_symmv_workspace *wmat = work->wmat;
  
  if (order_max > work->size || order_max <= order_min || order_min < 0)
    {
      GSL_ERROR ("invalid range [order_min,order_max]", GSL_EINVAL);
    }
  
  /* Convert the nonsymmetric tridiagonal matrix to a symmetric tridiagonal
     form. */

  tt[0] = 0.0;
  tt[1] = 0.0;
  tt[2] = qq;
  for (ii=1; ii<even_order-1; ii++)
  {
      tt[3*ii] = qq;
      tt[3*ii+1] = 4*ii*ii;
      tt[3*ii+2] = qq;
  }
  tt[3*even_order-3] = qq;
  tt[3*even_order-2] = 4*(even_order - 1)*(even_order - 1);
  tt[3*even_order-1] = 0.0;

  tt[3] *= 2;
  
  status = figi((signed int)even_order, tt, dd, ee, e2);

  if (status) 
    {
      GSL_ERROR("Internal error in tridiagonal Mathieu matrix", GSL_EFAILED);
    }

  /* Fill the period \pi matrix. */
  for (ii=0; ii<even_order*even_order; ii++)
      zz[ii] = 0.0;

  zz[0] = dd[0];
  zz[1] = ee[1];
  for (ii=1; ii<even_order-1; ii++)
  {
      zz[ii*even_order+ii-1] = ee[ii];
      zz[ii*even_order+ii] = dd[ii];
      zz[ii*even_order+ii+1] = ee[ii+1];
  }
  zz[even_order*(even_order-1)+even_order-2] = ee[even_order-1];
  zz[even_order*even_order-1] = dd[even_order-1];
  
  /* Compute (and sort) the eigenvalues of the matrix. */
  mat = gsl_matrix_view_array(zz, even_order, even_order);
  eval = gsl_vector_subvector(work->eval, 0, even_order);
  evec = gsl_matrix_submatrix(work->evec, 0, 0, even_order, even_order);
  gsl_eigen_symmv(&mat.matrix, &eval.vector, &evec.matrix, wmat);
  gsl_eigen_symmv_sort(&eval.vector, &evec.matrix, GSL_EIGEN_SORT_VAL_ASC);
  
  for (ii=0; ii<even_order-extra_values; ii++)
      aa[2*ii] = gsl_vector_get(&eval.vector, ii);
  
  /* Fill the period 2\pi matrix. */
  for (ii=0; ii<odd_order*odd_order; ii++)
      zz[ii] = 0.0;
  for (ii=0; ii<odd_order; ii++)
      for (jj=0; jj<odd_order; jj++)
      {
          if (ii == jj)
              zz[ii*odd_order+jj] = (2*ii + 1)*(2*ii + 1);
          else if (ii == jj + 1 || ii + 1 == jj)
              zz[ii*odd_order+jj] = qq;
      }
  zz[0] += qq;

  /* Compute (and sort) the eigenvalues of the matrix. */
  mat = gsl_matrix_view_array(zz, odd_order, odd_order);
  eval = gsl_vector_subvector(work->eval, 0, odd_order);
  evec = gsl_matrix_submatrix(work->evec, 0, 0, odd_order, odd_order);
  gsl_eigen_symmv(&mat.matrix, &eval.vector, &evec.matrix, wmat);
  gsl_eigen_symmv_sort(&eval.vector, &evec.matrix, GSL_EIGEN_SORT_VAL_ASC);

  for (ii=0; ii<odd_order-extra_values; ii++)
      aa[2*ii+1] = gsl_vector_get(&eval.vector, ii);

  for (ii = order_min ; ii <= order_max ; ii++)
    {
      result_array[ii - order_min] = aa[ii];
    }
  
  return GSL_SUCCESS;
}


int gsl_sf_mathieu_b_array(int order_min, int order_max, double qq, gsl_sf_mathieu_workspace *work, double result_array[])
{
  unsigned int even_order = work->even_order-1, odd_order = work->odd_order,
      extra_values = work->extra_values, ii, jj;
  double *zz = work->zz, *bb = work->bb;
  gsl_matrix_view mat, evec;
  gsl_vector_view eval;
  gsl_eigen_symmv_workspace *wmat = work->wmat;

  if (order_max > work->size || order_max <= order_min || order_min < 0)
    {
      GSL_ERROR ("invalid range [order_min,order_max]", GSL_EINVAL);
    }

  /* Fill the period \pi matrix. */
  for (ii=0; ii<even_order*even_order; ii++)
      zz[ii] = 0.0;
  for (ii=0; ii<even_order; ii++)
      for (jj=0; jj<even_order; jj++)
      {
          if (ii == jj)
              zz[ii*even_order+jj] = 4*(ii + 1)*(ii + 1);
          else if (ii == jj + 1 || ii + 1 == jj)
              zz[ii*even_order+jj] = qq;
      }

  /* Compute (and sort) the eigenvalues of the matrix. */
  mat = gsl_matrix_view_array(zz, even_order, even_order);
  eval = gsl_vector_subvector(work->eval, 0, even_order);
  evec = gsl_matrix_submatrix(work->evec, 0, 0, even_order, even_order);
  gsl_eigen_symmv(&mat.matrix, &eval.vector, &evec.matrix, wmat);
  gsl_eigen_symmv_sort(&eval.vector, &evec.matrix, GSL_EIGEN_SORT_VAL_ASC);

  bb[0] = 0.0;
  for (ii=0; ii<even_order-extra_values; ii++)
      bb[2*(ii+1)] = gsl_vector_get(&eval.vector, ii);
  
  /* Fill the period 2\pi matrix. */
  for (ii=0; ii<odd_order*odd_order; ii++)
      zz[ii] = 0.0;
  for (ii=0; ii<odd_order; ii++)
      for (jj=0; jj<odd_order; jj++)
      {
          if (ii == jj)
              zz[ii*odd_order+jj] = (2*ii + 1)*(2*ii + 1);
          else if (ii == jj + 1 || ii + 1 == jj)
              zz[ii*odd_order+jj] = qq;
      }

  zz[0] -= qq;

  /* Compute (and sort) the eigenvalues of the matrix. */
  mat = gsl_matrix_view_array(zz, odd_order, odd_order);
  eval = gsl_vector_subvector(work->eval, 0, odd_order);
  evec = gsl_matrix_submatrix(work->evec, 0, 0, odd_order, odd_order);
  gsl_eigen_symmv(&mat.matrix, &eval.vector, &evec.matrix, wmat);
  gsl_eigen_symmv_sort(&eval.vector, &evec.matrix, GSL_EIGEN_SORT_VAL_ASC);
  
  for (ii=0; ii<odd_order-extra_values; ii++)
      bb[2*ii+1] = gsl_vector_get(&eval.vector, ii);  

  for (ii = order_min ; ii <= order_max ; ii++)
    {
      result_array[ii - order_min] = bb[ii];
    }

  return GSL_SUCCESS;
}