Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Polylib.cpp
Go to the documentation of this file.
1 
2 #include <stdlib.h>
3 
4 #include <stdio.h>
5 
6 #include <math.h>
7 
8 #include "Polylib.h"
9 
10 #include <float.h>
11 
12 #include <complex>
13 
15 
16 
17 /// Maximum number of iterations in polynomial defalation routine Jacobz
18 
19 #define STOP 30
20 
21 /// Precision tolerance for two points to be similar
22 
23 #define EPS 100*DBL_EPSILON
24 
25 /// return the sign(b)*a
26 
27 #define sign(a,b) ((b)<0 ? -fabs(a) : fabs(a))
28 
29 
30 
31 namespace Polylib {
32 
33 
34 
35  /// Define whether to use polynomial deflation (1) or tridiagonal solver (0).
36 
37 #define POLYNOMIAL_DEFLATION 0
38 
39 
40 
41 #ifdef POLYNOMIAL_DEFLATION
42 
43  /// zero determination using Newton iteration with polynomial deflation
44 
45 #define jacobz(n,z,alpha,beta) Jacobz(n,z,alpha,beta)
46 
47 #else
48 
49  /// zero determination using eigenvalues of tridiagaonl matrix
50 
51 #define jacobz(n,z,alpha,beta) JacZeros(n,z,alpha,beta)
52 
53 #endif
54 
55 
56 
57 
58 
59 
60 
61  /* local functions */
62 
63  static void Jacobz (const int n, double *z, const double alpha,
64 
65  const double beta);
66 
67  // static void JacZeros (const int n, double *a, const double alpha,
68 
69  // const double beta);
70 
71  //static void TriQL (const int n, double *d, double *e);
72 
73  static void TriQL(const int, double *,double *, double **);
74 
75  double gammaF (const double);
76 
77  static void RecCoeff(const int, double *, double *,const double,
78 
79  const double);
80 
81  void JKMatrix(int, double *, double *);
82 
83  void chri1(int,double*,double*,double*,double*,double);
84 
85 
86 
87  /**
88 
89  \brief Gauss-Jacobi zeros and weights.
90 
91 
92 
93  \li Generate \a np Gauss Jacobi zeros, \a z, and weights,\a w,
94 
95  associated with the Jacobi polynomial \f$ P^{\alpha,\beta}_{np}(z)
96 
97  \f$,
98 
99 
100 
101  \li Exact for polynomials of order \a 2np-1 or less
102 
103  */
104 
105 
106 
107  void zwgj (double *z, double *w, const int np, const double alpha,
108 
109  const double beta)
110 
111  {
112 
113  register int i;
114 
115  double fac, one = 1.0, two = 2.0, apb = alpha + beta;
116 
117 
118 
119  jacobz (np,z,alpha,beta);
120 
121  jacobd (np,z,w,np,alpha,beta);
122 
123 
124 
125  fac = pow(two,apb + one)*gammaF(alpha + np + one)*gammaF(beta + np + one);
126 
127  fac /= gammaF(np + one)*gammaF(apb + np + one);
128 
129 
130 
131  for(i = 0; i < np; ++i) w[i] = fac/(w[i]*w[i]*(one-z[i]*z[i]));
132 
133 
134 
135  return;
136 
137  }
138 
139 
140 
141 
142 
143  /**
144 
145  \brief Gauss-Radau-Jacobi zeros and weights with end point at \a z=-1.
146 
147 
148 
149  \li Generate \a np Gauss-Radau-Jacobi zeros, \a z, and weights,\a w,
150 
151  associated with the polynomial \f$(1+z) P^{\alpha,\beta+1}_{np-1}(z)
152 
153  \f$.
154 
155 
156 
157  \li Exact for polynomials of order \a 2np-2 or less
158 
159  */
160 
161 
162 
163  void zwgrjm(double *z, double *w, const int np, const double alpha,
164 
165  const double beta)
166 
167  {
168 
169 
170 
171  if(np == 1){
172 
173  z[0] = 0.0;
174 
175  w[0] = 2.0;
176 
177  }
178 
179  else{
180 
181  register int i;
182 
183  double fac, one = 1.0, two = 2.0, apb = alpha + beta;
184 
185 
186 
187  z[0] = -one;
188 
189  jacobz (np-1,z+1,alpha,beta+1);
190 
191  jacobfd (np,z,w,NULL,np-1,alpha,beta);
192 
193 
194 
195  fac = pow(two,apb)*gammaF(alpha + np)*gammaF(beta + np);
196 
197  fac /= gammaF(np)*(beta + np)*gammaF(apb + np + 1);
198 
199 
200 
201  for(i = 0; i < np; ++i) w[i] = fac*(1-z[i])/(w[i]*w[i]);
202 
203  w[0] *= (beta + one);
204 
205  }
206 
207 
208 
209  return;
210 
211  }
212 
213 
214 
215 
216 
217  /**
218 
219  \brief Gauss-Radau-Jacobi zeros and weights with end point at \a z=1
220 
221 
222 
223 
224 
225  \li Generate \a np Gauss-Radau-Jacobi zeros, \a z, and weights,\a w,
226 
227  associated with the polynomial \f$(1-z) P^{\alpha+1,\beta}_{np-1}(z)
228 
229  \f$.
230 
231 
232 
233  \li Exact for polynomials of order \a 2np-2 or less
234 
235  */
236 
237 
238 
239  void zwgrjp(double *z, double *w, const int np, const double alpha,
240 
241  const double beta)
242 
243  {
244 
245 
246 
247  if(np == 1){
248 
249  z[0] = 0.0;
250 
251  w[0] = 2.0;
252 
253  }
254 
255  else{
256 
257  register int i;
258 
259  double fac, one = 1.0, two = 2.0, apb = alpha + beta;
260 
261 
262 
263  jacobz (np-1,z,alpha+1,beta);
264 
265  z[np-1] = one;
266 
267  jacobfd (np,z,w,NULL,np-1,alpha,beta);
268 
269 
270 
271  fac = pow(two,apb)*gammaF(alpha + np)*gammaF(beta + np);
272 
273  fac /= gammaF(np)*(alpha + np)*gammaF(apb + np + 1);
274 
275 
276 
277  for(i = 0; i < np; ++i) w[i] = fac*(1+z[i])/(w[i]*w[i]);
278 
279  w[np-1] *= (alpha + one);
280 
281  }
282 
283 
284 
285  return;
286 
287  }
288 
289 
290 
291 
292 
293  /**
294 
295  \brief Gauss-Lobatto-Jacobi zeros and weights with end point at \a z=-1,\a 1
296 
297 
298 
299 
300 
301  \li Generate \a np Gauss-Lobatto-Jacobi points, \a z, and weights, \a w,
302 
303  associated with polynomial \f$ (1-z)(1+z) P^{\alpha+1,\beta+1}_{np-2}(z) \f$
304 
305  \li Exact for polynomials of order \a 2np-3 or less
306 
307  */
308 
309 
310 
311  void zwglj(double *z, double *w, const int np, const double alpha,
312 
313  const double beta)
314 
315  {
316 
317 
318 
319  if( np == 1 ){
320 
321  z[0] = 0.0;
322 
323  w[0] = 2.0;
324 
325  }
326 
327  else{
328 
329  register int i;
330 
331  double fac, one = 1.0, apb = alpha + beta, two = 2.0;
332 
333 
334 
335  z[0] = -one;
336 
337  z[np-1] = one;
338 
339  jacobz (np-2,z + 1,alpha + one,beta + one);
340 
341  jacobfd (np,z,w,NULL,np-1,alpha,beta);
342 
343 
344 
345  fac = pow(two,apb + 1)*gammaF(alpha + np)*gammaF(beta + np);
346 
347  fac /= (np-1)*gammaF(np)*gammaF(alpha + beta + np + one);
348 
349 
350 
351  for(i = 0; i < np; ++i) w[i] = fac/(w[i]*w[i]);
352 
353  w[0] *= (beta + one);
354 
355  w[np-1] *= (alpha + one);
356 
357  }
358 
359 
360 
361  return;
362 
363  }
364 
365 
366 
367  /**
368 
369  \brief Gauss-Kronrod-Jacobi zeros and weights.
370 
371 
372 
373  \li Generate \a npt=2*np+1 Gauss-Kronrod Jacobi zeros, \a z, and weights,\a w,
374 
375  associated with the Jacobi polynomial \f$ P^{\alpha,\beta}_{np}(z)
376 
377  \f$,
378 
379 
380 
381  \li Exact for polynomials of order \a 3np+1 or less
382 
383  */
384 
385  void zwgk(double *z, double *w, const int npt , const double alpha,
386 
387  const double beta)
388 
389  {
390 
391 
392 
393  int np = (npt-1)/2;
394 
395 
396 
397  int i,j;
398 
399 
400 
401  // number of kronrod points associated with the np gauss rule
402 
403  int kpoints = 2*np + 1;
404 
405 
406 
407  // Define the number of required recurrence coefficents
408 
409  int ncoeffs = (int)floor(3.0*(np+1)/2);
410 
411 
412 
413  // Define arrays for the recurrence coefficients
414 
415  // We will use these arrays for the Kronrod results too, hence the
416 
417  // reason for the size of the arrays
418 
419  double *a = new double[kpoints];
420 
421  double *b = new double[kpoints];
422 
423 
424 
425  // Initialize a and b to zero
426 
427  for(i = 0; i < kpoints; i++)
428 
429  {
430 
431  a[i] = 0.0;
432 
433  b[i] = 0.0;
434 
435  }
436 
437 
438 
439  // Call the routine to calculate the recurrence coefficients
440 
441  RecCoeff(ncoeffs,a,b,alpha,beta);
442 
443 
444 
445  // Call the routine to caluclate the jacobi-Kronrod matrix
446 
447  JKMatrix(np,a,b);
448 
449 
450 
451  // Set up the identity matrix
452 
453  double** zmatrix = new double*[kpoints];
454 
455  for(i = 0; i < kpoints; i++)
456 
457  {
458 
459  zmatrix[i] = new double[kpoints];
460 
461  for(j = 0; j < kpoints; j++)
462 
463  {
464 
465  zmatrix[i][j] = 0.0;
466 
467  }
468 
469  }
470 
471  for(i = 0; i < kpoints; i++)
472 
473  {
474 
475  zmatrix[i][i] = 1.0;
476 
477  }
478 
479 
480 
481  // Calculte the points and weights
482 
483  TriQL(kpoints, a, b, zmatrix);
484 
485 
486 
487  for(i = 0; i < kpoints; i++)
488 
489  {
490 
491  z[i] = a[i];
492 
493  w[i] = b[i];
494 
495  }
496 
497 
498 
499 
500 
501  }
502 
503 
504 
505  /**
506 
507  \brief Gauss-Radau-Kronrod-Jacobi zeros and weights.
508 
509 
510 
511  \li Generate \a npt=2*np Radau-Kronrod Jacobi zeros, \a z, and weights,\a w,
512 
513  associated with the Jacobi polynomial \f$ P^{\alpha,\beta}_{np}(z)
514 
515  \f$,
516 
517  */
518 
519  void zwrk(double* z, double* w, const int npt ,const double alpha,
520 
521  const double beta)
522 
523  {
524 
525 
526 
527  int np = npt/2;
528 
529 
530 
531  if(np < 2)
532 
533  {
534 
535  fprintf(stderr,"too few points in formula\n");
536 
537  return;
538 
539  }
540 
541 
542 
543  double end0 = -1;
544 
545 
546 
547  int i,j;
548 
549 
550 
551  // number of kronrod points associated with the np gauss rule
552 
553  int kpoints = 2*np;
554 
555 
556 
557  // Define the number of required recurrence coefficents
558 
559  int ncoeffs = (int)ceil(3.0*np/2);
560 
561 
562 
563  // Define arrays for the recurrence coefficients
564 
565  double *a = new double[ncoeffs+1];
566 
567  double *b = new double[ncoeffs+1];
568 
569 
570 
571  // Initialize a and b to zero
572 
573  for(i = 0; i < ncoeffs+1; i++)
574 
575  {
576 
577  a[i] = 0.0;
578 
579  b[i] = 0.0;
580 
581  }
582 
583 
584 
585  // Call the routine to calculate the recurrence coefficients
586 
587  RecCoeff(ncoeffs,a,b,alpha,beta);
588 
589 
590 
591  double* a0 = new double[ncoeffs];
592 
593  double* b0 = new double[ncoeffs];
594 
595 
596 
597  chri1(ncoeffs,a,b,a0,b0,end0);
598 
599 
600 
601  double s = b0[0]/fabs(b0[0]);
602 
603  b0[0] = s*b0[0];
604 
605 
606 
607  // Finding the 2*np-1 gauss-kronrod points
608 
609  double* z1 = new double[2*np-1];
610 
611  double* w1 = new double[2*np-1];
612 
613  for(i = 0; i < ncoeffs; i++)
614 
615  {
616 
617  z1[i] = a0[i];
618 
619  w1[i] = b0[i];
620 
621  }
622 
623  JKMatrix(np-1,z1,w1);
624 
625  // Set up the identity matrix
626 
627  double** zmatrix = new double*[2*np-1];
628 
629  for(i = 0; i < 2*np-1; i++)
630 
631  {
632 
633  zmatrix[i] = new double[2*np-1];
634 
635  for(j = 0; j < 2*np-1; j++)
636 
637  {
638 
639  zmatrix[i][j] = 0.0;
640 
641  }
642 
643  }
644 
645  for(i = 0; i < 2*np-1; i++)
646 
647  {
648 
649  zmatrix[i][i] = 1.0;
650 
651  }
652 
653 
654 
655  // Calculate the points and weights
656 
657  TriQL(2*np-1, z1, w1, zmatrix);
658 
659 
660 
661  double sumW1 = 0.0;
662 
663  for(i = 0; i < 2*np-1; i++)
664 
665  {
666 
667  w1[i] = s*w1[i]/(z1[i]-end0);
668 
669  sumW1 += w1[i];
670 
671  }
672 
673 
674 
675  z[0] = end0;
676 
677  w[0] = b[0]- sumW1;
678 
679  for(i = 1; i < kpoints; i++)
680 
681  {
682 
683  z[i] = z1[i-1];
684 
685  w[i] = w1[i-1];
686 
687  }
688 
689 
690 
691 
692 
693  }
694 
695 
696 
697  /**
698 
699  \brief Gauss-Lobatto-Kronrod-Jacobi zeros and weights.
700 
701 
702 
703  \li Generate \a npt=2*np-1 Lobatto-Kronrod Jacobi zeros, \a z, and weights,\a w,
704 
705  associated with the Jacobi polynomial \f$ P^{\alpha,\beta}_{np}(z)
706 
707  \f$,
708 
709  */
710 
711  void zwlk(double* z, double* w, const int npt,
712 
713  const double alpha, const double beta)
714 
715  {
716 
717 
718 
719  int np = (npt+1)/2;
720 
721 
722 
723  if(np < 4)
724 
725  {
726 
727  fprintf (stderr,"too few points in formula\n");
728 
729  return;
730 
731  }
732 
733 
734 
735  double endl = -1;
736 
737  double endr = 1;
738 
739  int i,j;
740 
741 
742 
743  // number of kronrod points associated with the np gauss rule
744 
745  int kpoints = 2*np-1;
746 
747 
748 
749  // Define the number of required recurrence coefficents
750 
751  int ncoeffs = (int)ceil(3.0*np/2)-1;
752 
753 
754 
755  // Define arrays for the recurrence coefficients
756 
757  double *a = new double[ncoeffs+1];
758 
759  double *b = new double[ncoeffs+1];
760 
761 
762 
763  // Initialize a and b to zero
764 
765  for(i = 0; i < ncoeffs+1; i++)
766 
767  {
768 
769  a[i] = 0.0;
770 
771  b[i] = 0.0;
772 
773  }
774 
775 
776 
777  // Call the routine to calculate the recurrence coefficients
778 
779  RecCoeff(ncoeffs,a,b,alpha,beta);
780 
781 
782 
783 
784 
785  double* a0 = new double[ncoeffs];
786 
787  double* b0 = new double[ncoeffs];
788 
789 
790 
791  chri1(ncoeffs,a,b,a0,b0,endl);
792 
793 
794 
795  double* a1 = new double[ncoeffs-1];
796 
797  double* b1 = new double[ncoeffs-1];
798 
799 
800 
801  chri1(ncoeffs-1,a0,b0,a1,b1,endr);
802 
803 
804 
805 
806 
807  double s = b1[0]/fabs(b1[0]);
808 
809  b1[0] = s*b1[0];
810 
811 
812 
813  // Finding the 2*np-1 gauss-kronrod points
814 
815  double* z1 = new double[2*np-3];
816 
817  double* w1 = new double[2*np-3];
818 
819  for(i = 0; i < ncoeffs; i++)
820 
821  {
822 
823  z1[i] = a1[i];
824 
825  w1[i] = b1[i];
826 
827  }
828 
829  JKMatrix(np-2,z1,w1);
830 
831  // Set up the identity matrix
832 
833  double** zmatrix = new double*[2*np-3];
834 
835  for(i = 0; i < 2*np-3; i++)
836 
837  {
838 
839  zmatrix[i] = new double[2*np-3];
840 
841  for(j = 0; j < 2*np-3; j++)
842 
843  {
844 
845  zmatrix[i][j] = 0.0;
846 
847  }
848 
849  }
850 
851  for(i = 0; i < 2*np-3; i++)
852 
853  {
854 
855  zmatrix[i][i] = 1.0;
856 
857  }
858 
859 
860 
861  // Calculate the points and weights
862 
863  TriQL(2*np-3, z1, w1, zmatrix);
864 
865 
866 
867  double sumW1 = 0.0;
868 
869  double sumW1Z1 = 0.0;
870 
871  for(i = 0; i < 2*np-3; i++)
872 
873  {
874 
875  w1[i] = s*w1[i]/(z1[i]-endl)/(z1[i]-endr);
876 
877  sumW1 += w1[i];
878 
879  sumW1Z1 += z1[i]*w1[i];
880 
881  }
882 
883 
884 
885  double c0 = b[0]-sumW1;
886 
887  double c1 = a[0]*b[0]-sumW1Z1;
888 
889 
890 
891  z[0] = endl;
892 
893  z[2*np-2] = endr;
894 
895  w[0] = (c0*endr-c1)/(endr-endl);
896 
897  w[2*np-2] = (c1-c0*endl)/(endr-endl);
898 
899 
900 
901  for(i = 1; i < kpoints-1; i++)
902 
903  {
904 
905  z[i] = z1[i-1];
906 
907  w[i] = w1[i-1];
908 
909  }
910 
911  }
912 
913 
914 
915  /**
916 
917  \brief Compute the Derivative Matrix and its transpose associated
918 
919  with the Gauss-Jacobi zeros.
920 
921 
922 
923  \li Compute the derivative matrix, \a d, and its transpose, \a dt,
924 
925  associated with the n_th order Lagrangian interpolants through the
926 
927  \a np Gauss-Jacobi points \a z such that \n
928 
929  \f$ \frac{du}{dz}(z[i]) = \sum_{j=0}^{np-1} D[i*np+j] u(z[j]) \f$
930 
931 
932 
933  */
934 
935 
936 
937  void Dgj(double *D, const double *z, const int np, const double alpha,
938 
939  const double beta)
940 
941  {
942 
943 
944 
945  double one = 1.0, two = 2.0;
946 
947 
948 
949  if (np <= 0){
950 
951  D[0] = 0.0;
952 
953  }
954 
955  else{
956 
957  register int i,j;
958 
959  double *pd;
960 
961 
962 
963  pd = (double *)malloc(np*sizeof(double));
964 
965  jacobd(np,z,pd,np,alpha,beta);
966 
967 
968 
969  for (i = 0; i < np; i++){
970 
971  for (j = 0; j < np; j++){
972 
973 
974 
975  if (i != j)
976 
977  D[i*np+j] = pd[j]/(pd[i]*(z[j]-z[i]));
978 
979  else
980 
981  D[i*np+j] = (alpha - beta + (alpha + beta + two)*z[j])/
982 
983  (two*(one - z[j]*z[j]));
984 
985  }
986 
987  }
988 
989  free(pd);
990 
991  }
992 
993  return;
994 
995  }
996 
997 
998 
999 
1000 
1001  /**
1002 
1003  \brief Compute the Derivative Matrix and its transpose associated
1004 
1005  with the Gauss-Radau-Jacobi zeros with a zero at \a z=-1.
1006 
1007 
1008 
1009  \li Compute the derivative matrix, \a d, associated with the n_th
1010 
1011  order Lagrangian interpolants through the \a np Gauss-Radau-Jacobi
1012 
1013  points \a z such that \n \f$ \frac{du}{dz}(z[i]) =
1014 
1015  \sum_{j=0}^{np-1} D[i*np+j] u(z[j]) \f$
1016 
1017 
1018 
1019  */
1020 
1021 
1022 
1023  void Dgrjm(double *D, const double *z, const int np, const double alpha,
1024 
1025  const double beta)
1026 
1027  {
1028 
1029 
1030 
1031  if (np <= 0){
1032 
1033  D[0] = 0.0;
1034 
1035  }
1036 
1037  else{
1038 
1039  register int i, j;
1040 
1041  double one = 1.0, two = 2.0;
1042 
1043  double *pd;
1044 
1045 
1046 
1047  pd = (double *)malloc(np*sizeof(double));
1048 
1049 
1050 
1051  pd[0] = pow(-one,np-1)*gammaF(np+beta+one);
1052 
1053  pd[0] /= gammaF(np)*gammaF(beta+two);
1054 
1055  jacobd(np-1,z+1,pd+1,np-1,alpha,beta+1);
1056 
1057  for(i = 1; i < np; ++i) pd[i] *= (1+z[i]);
1058 
1059 
1060 
1061  for (i = 0; i < np; i++)
1062 
1063  for (j = 0; j < np; j++){
1064 
1065  if (i != j)
1066 
1067  D[i*np+j] = pd[j]/(pd[i]*(z[j]-z[i]));
1068 
1069  else {
1070 
1071  if(j == 0)
1072 
1073  D[i*np+j] = -(np + alpha + beta + one)*(np - one)/
1074 
1075  (two*(beta + two));
1076 
1077  else
1078 
1079  D[i*np+j] = (alpha - beta + one + (alpha + beta + one)*z[j])/
1080 
1081  (two*(one - z[j]*z[j]));
1082 
1083  }
1084 
1085  }
1086 
1087  free(pd);
1088 
1089  }
1090 
1091 
1092 
1093  return;
1094 
1095  }
1096 
1097 
1098 
1099 
1100 
1101  /**
1102 
1103  \brief Compute the Derivative Matrix associated with the
1104 
1105  Gauss-Radau-Jacobi zeros with a zero at \a z=1.
1106 
1107 
1108 
1109  \li Compute the derivative matrix, \a d, associated with the n_th
1110 
1111  order Lagrangian interpolants through the \a np Gauss-Radau-Jacobi
1112 
1113  points \a z such that \n \f$ \frac{du}{dz}(z[i]) =
1114 
1115  \sum_{j=0}^{np-1} D[i*np+j] u(z[j]) \f$
1116 
1117  */
1118 
1119 
1120 
1121  void Dgrjp(double *D, const double *z, const int np, const double alpha,
1122 
1123  const double beta)
1124 
1125  {
1126 
1127 
1128 
1129  if (np <= 0){
1130 
1131  D[0] = 0.0;
1132 
1133  }
1134 
1135  else{
1136 
1137  register int i, j;
1138 
1139  double one = 1.0, two = 2.0;
1140 
1141  double *pd;
1142 
1143 
1144 
1145  pd = (double *)malloc(np*sizeof(double));
1146 
1147 
1148 
1149 
1150 
1151  jacobd(np-1,z,pd,np-1,alpha+1,beta);
1152 
1153  for(i = 0; i < np-1; ++i) pd[i] *= (1-z[i]);
1154 
1155  pd[np-1] = -gammaF(np+alpha+one);
1156 
1157  pd[np-1] /= gammaF(np)*gammaF(alpha+two);
1158 
1159 
1160 
1161  for (i = 0; i < np; i++)
1162 
1163  for (j = 0; j < np; j++){
1164 
1165  if (i != j)
1166 
1167  D[i*np+j] = pd[j]/(pd[i]*(z[j]-z[i]));
1168 
1169  else {
1170 
1171  if(j == np-1)
1172 
1173  D[i*np+j] = (np + alpha + beta + one)*(np - one)/
1174 
1175  (two*(alpha + two));
1176 
1177  else
1178 
1179  D[i*np+j] = (alpha - beta - one + (alpha + beta + one)*z[j])/
1180 
1181  (two*(one - z[j]*z[j]));
1182 
1183  }
1184 
1185  }
1186 
1187  free(pd);
1188 
1189  }
1190 
1191 
1192 
1193  return;
1194 
1195  }
1196 
1197 
1198 
1199  /**
1200 
1201  \brief Compute the Derivative Matrix associated with the
1202 
1203  Gauss-Lobatto-Jacobi zeros.
1204 
1205 
1206 
1207  \li Compute the derivative matrix, \a d, associated with the n_th
1208 
1209  order Lagrange interpolants through the \a np
1210 
1211  Gauss-Lobatto-Jacobi points \a z such that \n \f$
1212 
1213  \frac{du}{dz}(z[i]) = \sum_{j=0}^{np-1} D[i*np+j] u(z[j]) \f$
1214 
1215 
1216 
1217  */
1218 
1219 
1220 
1221  void Dglj(double *D, const double *z, const int np, const double alpha,
1222 
1223  const double beta)
1224 
1225  {
1226 
1227 
1228 
1229  if (np <= 0){
1230 
1231  D[0] = 0.0;
1232 
1233  }
1234 
1235  else{
1236 
1237  register int i, j;
1238 
1239  double one = 1.0, two = 2.0;
1240 
1241  double *pd;
1242 
1243 
1244 
1245  pd = (double *)malloc(np*sizeof(double));
1246 
1247 
1248 
1249  pd[0] = two*pow(-one,np)*gammaF(np + beta);
1250 
1251  pd[0] /= gammaF(np - one)*gammaF(beta + two);
1252 
1253  jacobd(np-2,z+1,pd+1,np-2,alpha+1,beta+1);
1254 
1255  for(i = 1; i < np-1; ++i) pd[i] *= (one-z[i]*z[i]);
1256 
1257  pd[np-1] = -two*gammaF(np + alpha);
1258 
1259  pd[np-1] /= gammaF(np - one)*gammaF(alpha + two);
1260 
1261 
1262 
1263  for (i = 0; i < np; i++)
1264 
1265  for (j = 0; j < np; j++){
1266 
1267  if (i != j)
1268 
1269  D[i*np+j] = pd[j]/(pd[i]*(z[j]-z[i]));
1270 
1271  else {
1272 
1273  if (j == 0)
1274 
1275  D[i*np+j] = (alpha - (np-1)*(np + alpha + beta))/(two*(beta+ two));
1276 
1277  else if (j == np-1)
1278 
1279  D[i*np+j] =-(beta - (np-1)*(np + alpha + beta))/(two*(alpha+ two));
1280 
1281  else
1282 
1283  D[i*np+j] = (alpha - beta + (alpha + beta)*z[j])/
1284 
1285  (two*(one - z[j]*z[j]));
1286 
1287  }
1288 
1289  }
1290 
1291  free(pd);
1292 
1293  }
1294 
1295 
1296 
1297  return;
1298 
1299  }
1300 
1301 
1302 
1303 
1304 
1305  /**
1306 
1307  \brief Compute the value of the \a i th Lagrangian interpolant through
1308 
1309  the \a np Gauss-Jacobi points \a zgj at the arbitrary location \a z.
1310 
1311 
1312 
1313  \li \f$ -1 \leq z \leq 1 \f$
1314 
1315 
1316 
1317  \li Uses the defintion of the Lagrangian interpolant:\n
1318 
1319  %
1320 
1321  \f$ \begin{array}{rcl}
1322 
1323  h_j(z) = \left\{ \begin{array}{ll}
1324 
1325  \displaystyle \frac{P_{np}^{\alpha,\beta}(z)}
1326 
1327  {[P_{np}^{\alpha,\beta}(z_j)]^\prime
1328 
1329  (z-z_j)} & \mbox{if $z \ne z_j$}\\
1330 
1331  & \\
1332 
1333  1 & \mbox{if $z=z_j$}
1334 
1335  \end{array}
1336 
1337  \right.
1338 
1339  \end{array} \f$
1340 
1341  */
1342 
1343 
1344 
1345  double hgj (const int i, const double z, const double *zgj,
1346 
1347  const int np, const double alpha, const double beta)
1348 
1349  {
1350 
1351 
1352 
1353  double zi, dz, p, pd, h;
1354 
1355 
1356 
1357  zi = *(zgj+i);
1358 
1359  dz = z - zi;
1360 
1361  if (fabs(dz) < EPS) return 1.0;
1362 
1363 
1364 
1365  jacobd (1, &zi, &pd , np, alpha, beta);
1366 
1367  jacobfd(1, &z , &p, NULL , np, alpha, beta);
1368 
1369  h = p/(pd*dz);
1370 
1371 
1372 
1373  return h;
1374 
1375  }
1376 
1377 
1378 
1379  /**
1380 
1381  \brief Compute the value of the \a i th Lagrangian interpolant through the
1382 
1383  \a np Gauss-Radau-Jacobi points \a zgrj at the arbitrary location
1384 
1385  \a z. This routine assumes \a zgrj includes the point \a -1.
1386 
1387 
1388 
1389  \li \f$ -1 \leq z \leq 1 \f$
1390 
1391 
1392 
1393  \li Uses the defintion of the Lagrangian interpolant:\n
1394 
1395  %
1396 
1397  \f$ \begin{array}{rcl}
1398 
1399  h_j(z) = \left\{ \begin{array}{ll}
1400 
1401  \displaystyle \frac{(1+z) P_{np-1}^{\alpha,\beta+1}(z)}
1402 
1403  {((1+z_j) [P_{np-1}^{\alpha,\beta+1}(z_j)]^\prime +
1404 
1405  P_{np-1}^{\alpha,\beta+1}(z_j) ) (z-z_j)} & \mbox{if $z \ne z_j$}\\
1406 
1407  & \\
1408 
1409  1 & \mbox{if $z=z_j$}
1410 
1411  \end{array}
1412 
1413  \right.
1414 
1415  \end{array} \f$
1416 
1417  */
1418 
1419 
1420 
1421  double hgrjm (const int i, const double z, const double *zgrj, const int np,
1422 
1423  const double alpha, const double beta)
1424 
1425  {
1426 
1427 
1428 
1429  double zi, dz, p, pd, h;
1430 
1431 
1432 
1433  zi = *(zgrj+i);
1434 
1435  dz = z - zi;
1436 
1437  if (fabs(dz) < EPS) return 1.0;
1438 
1439 
1440 
1441  jacobfd (1, &zi, &p , NULL, np-1, alpha, beta + 1);
1442 
1443  // need to use this routine in caes zi = -1 or 1
1444 
1445  jacobd (1, &zi, &pd, np-1, alpha, beta + 1);
1446 
1447  h = (1.0 + zi)*pd + p;
1448 
1449  jacobfd (1, &z, &p, NULL, np-1, alpha, beta + 1);
1450 
1451  h = (1.0 + z )*p/(h*dz);
1452 
1453 
1454 
1455  return h;
1456 
1457  }
1458 
1459 
1460 
1461 
1462 
1463  /**
1464 
1465  \brief Compute the value of the \a i th Lagrangian interpolant through the
1466 
1467  \a np Gauss-Radau-Jacobi points \a zgrj at the arbitrary location
1468 
1469  \a z. This routine assumes \a zgrj includes the point \a +1.
1470 
1471 
1472 
1473  \li \f$ -1 \leq z \leq 1 \f$
1474 
1475 
1476 
1477  \li Uses the defintion of the Lagrangian interpolant:\n
1478 
1479  %
1480 
1481  \f$ \begin{array}{rcl}
1482 
1483  h_j(z) = \left\{ \begin{array}{ll}
1484 
1485  \displaystyle \frac{(1-z) P_{np-1}^{\alpha+1,\beta}(z)}
1486 
1487  {((1-z_j) [P_{np-1}^{\alpha+1,\beta}(z_j)]^\prime -
1488 
1489  P_{np-1}^{\alpha+1,\beta}(z_j) ) (z-z_j)} & \mbox{if $z \ne z_j$}\\
1490 
1491  & \\
1492 
1493  1 & \mbox{if $z=z_j$}
1494 
1495  \end{array}
1496 
1497  \right.
1498 
1499  \end{array} \f$
1500 
1501  */
1502 
1503 
1504 
1505  double hgrjp (const int i, const double z, const double *zgrj, const int np,
1506 
1507  const double alpha, const double beta)
1508 
1509  {
1510 
1511 
1512 
1513  double zi, dz, p, pd, h;
1514 
1515 
1516 
1517  zi = *(zgrj+i);
1518 
1519  dz = z - zi;
1520 
1521  if (fabs(dz) < EPS) return 1.0;
1522 
1523 
1524 
1525  jacobfd (1, &zi, &p , NULL, np-1, alpha+1, beta );
1526 
1527  // need to use this routine in caes z = -1 or 1
1528 
1529  jacobd (1, &zi, &pd, np-1, alpha+1, beta );
1530 
1531  h = (1.0 - zi)*pd - p;
1532 
1533  jacobfd (1, &z, &p, NULL, np-1, alpha+1, beta);
1534 
1535  h = (1.0 - z )*p/(h*dz);
1536 
1537 
1538 
1539  return h;
1540 
1541  }
1542 
1543 
1544 
1545 
1546 
1547  /**
1548 
1549  \brief Compute the value of the \a i th Lagrangian interpolant through the
1550 
1551  \a np Gauss-Lobatto-Jacobi points \a zgrj at the arbitrary location
1552 
1553  \a z.
1554 
1555 
1556 
1557  \li \f$ -1 \leq z \leq 1 \f$
1558 
1559 
1560 
1561  \li Uses the defintion of the Lagrangian interpolant:\n
1562 
1563  %
1564 
1565  \f$ \begin{array}{rcl}
1566 
1567  h_j(z) = \left\{ \begin{array}{ll}
1568 
1569  \displaystyle \frac{(1-z^2) P_{np-2}^{\alpha+1,\beta+1}(z)}
1570 
1571  {((1-z^2_j) [P_{np-2}^{\alpha+1,\beta+1}(z_j)]^\prime -
1572 
1573  2 z_j P_{np-2}^{\alpha+1,\beta+1}(z_j) ) (z-z_j)}&\mbox{if $z \ne z_j$}\\
1574 
1575  & \\
1576 
1577  1 & \mbox{if $z=z_j$}
1578 
1579  \end{array}
1580 
1581  \right.
1582 
1583  \end{array} \f$
1584 
1585  */
1586 
1587 
1588 
1589  double hglj (const int i, const double z, const double *zglj, const int np,
1590 
1591  const double alpha, const double beta)
1592 
1593  {
1594 
1595  double one = 1., two = 2.;
1596 
1597  double zi, dz, p, pd, h;
1598 
1599 
1600 
1601  zi = *(zglj+i);
1602 
1603  dz = z - zi;
1604 
1605  if (fabs(dz) < EPS) return 1.0;
1606 
1607 
1608 
1609  jacobfd(1, &zi, &p , NULL, np-2, alpha + one, beta + one);
1610 
1611  // need to use this routine in caes z = -1 or 1
1612 
1613  jacobd (1, &zi, &pd, np-2, alpha + one, beta + one);
1614 
1615  h = (one - zi*zi)*pd - two*zi*p;
1616 
1617  jacobfd(1, &z, &p, NULL, np-2, alpha + one, beta + one);
1618 
1619  h = (one - z*z)*p/(h*dz);
1620 
1621 
1622 
1623  return h;
1624 
1625  }
1626 
1627 
1628 
1629 
1630 
1631  /**
1632 
1633  \brief Interpolation Operator from Gauss-Jacobi points to an
1634 
1635  arbitrary distribution at points \a zm
1636 
1637 
1638 
1639  \li Computes the one-dimensional interpolation matrix, \a im, to
1640 
1641  interpolate a function from at Gauss-Jacobi distribution of \a nz
1642 
1643  zeros \a zgrj to an arbitrary distribution of \a mz points \a zm, i.e.\n
1644 
1645  \f$
1646 
1647  u(zm[i]) = \sum_{j=0}^{nz-1} im[i*nz+j] \ u(zgj[j])
1648 
1649  \f$
1650 
1651 
1652 
1653  */
1654 
1655 
1656 
1657  void Imgj(double *im, const double *zgj, const double *zm, const int nz,
1658 
1659  const int mz,const double alpha, const double beta){
1660 
1661  double zp;
1662 
1663  register int i, j;
1664 
1665 
1666 
1667  for (i = 0; i < nz; ++i) {
1668 
1669  for (j = 0; j < mz; ++j)
1670 
1671  {
1672 
1673  zp = zm[j];
1674 
1675  im [i*mz+j] = hgj(i, zp, zgj, nz, alpha, beta);
1676 
1677  }
1678 
1679  }
1680 
1681 
1682 
1683  return;
1684 
1685  }
1686 
1687 
1688 
1689  /**
1690 
1691  \brief Interpolation Operator from Gauss-Radau-Jacobi points
1692 
1693  (including \a z=-1) to an arbitrary distrubtion at points \a zm
1694 
1695 
1696 
1697  \li Computes the one-dimensional interpolation matrix, \a im, to
1698 
1699  interpolate a function from at Gauss-Radau-Jacobi distribution of
1700 
1701  \a nz zeros \a zgrj (where \a zgrj[0]=-1) to an arbitrary
1702 
1703  distribution of \a mz points \a zm, i.e.
1704 
1705  \n
1706 
1707  \f$ u(zm[i]) = \sum_{j=0}^{nz-1} im[i*nz+j] \ u(zgj[j]) \f$
1708 
1709 
1710 
1711  */
1712 
1713 
1714 
1715  void Imgrjm(double *im, const double *zgrj, const double *zm, const int nz,
1716 
1717  const int mz, const double alpha, const double beta)
1718 
1719  {
1720 
1721  double zp;
1722 
1723  register int i, j;
1724 
1725 
1726 
1727  for (i = 0; i < nz; i++) {
1728 
1729  for (j = 0; j < mz; j++)
1730 
1731  {
1732 
1733  zp = zm[j];
1734 
1735  im [i*mz+j] = hgrjm(i, zp, zgrj, nz, alpha, beta);
1736 
1737  }
1738 
1739  }
1740 
1741 
1742 
1743  return;
1744 
1745  }
1746 
1747 
1748 
1749  /**
1750 
1751  \brief Interpolation Operator from Gauss-Radau-Jacobi points
1752 
1753  (including \a z=1) to an arbitrary distrubtion at points \a zm
1754 
1755 
1756 
1757  \li Computes the one-dimensional interpolation matrix, \a im, to
1758 
1759  interpolate a function from at Gauss-Radau-Jacobi distribution of
1760 
1761  \a nz zeros \a zgrj (where \a zgrj[nz-1]=1) to an arbitrary
1762 
1763  distribution of \a mz points \a zm, i.e.
1764 
1765  \n
1766 
1767  \f$ u(zm[i]) = \sum_{j=0}^{nz-1} im[i*nz+j] \ u(zgj[j]) \f$
1768 
1769 
1770 
1771  */
1772 
1773 
1774 
1775  void Imgrjp(double *im, const double *zgrj, const double *zm, const int nz,
1776 
1777  const int mz,const double alpha, const double beta)
1778 
1779  {
1780 
1781  double zp;
1782 
1783  register int i, j;
1784 
1785 
1786 
1787  for (i = 0; i < nz; i++) {
1788 
1789  for (j = 0; j < mz; j++)
1790 
1791  {
1792 
1793  zp = zm[j];
1794 
1795  im [i*mz+j] = hgrjp(i, zp, zgrj, nz, alpha, beta);
1796 
1797  }
1798 
1799  }
1800 
1801 
1802 
1803  return;
1804 
1805  }
1806 
1807 
1808 
1809 
1810 
1811  /**
1812 
1813  \brief Interpolation Operator from Gauss-Lobatto-Jacobi points
1814 
1815  to an arbitrary distrubtion at points \a zm
1816 
1817 
1818 
1819  \li Computes the one-dimensional interpolation matrix, \a im, to
1820 
1821  interpolate a function from at Gauss-Lobatto-Jacobi distribution of
1822 
1823  \a nz zeros \a zgrj (where \a zgrj[0]=-1) to an arbitrary
1824 
1825  distribution of \a mz points \a zm, i.e.
1826 
1827  \n
1828 
1829  \f$ u(zm[i]) = \sum_{j=0}^{nz-1} im[i*nz+j] \ u(zgj[j]) \f$
1830 
1831 
1832 
1833  */
1834 
1835 
1836 
1837  void Imglj(double *im, const double *zglj, const double *zm, const int nz,
1838 
1839  const int mz, const double alpha, const double beta)
1840 
1841  {
1842 
1843  double zp;
1844 
1845  register int i, j;
1846 
1847 
1848 
1849  for (i = 0; i < nz; i++) {
1850 
1851  for (j = 0; j < mz; j++)
1852 
1853  {
1854 
1855  zp = zm[j];
1856 
1857  im[i*mz+j] = hglj(i, zp, zglj, nz, alpha, beta);
1858 
1859  }
1860 
1861  }
1862 
1863 
1864 
1865  return;
1866 
1867  }
1868 
1869 
1870 
1871  /**
1872 
1873  \brief Routine to calculate Jacobi polynomials, \f$
1874 
1875  P^{\alpha,\beta}_n(z) \f$, and their first derivative, \f$
1876 
1877  \frac{d}{dz} P^{\alpha,\beta}_n(z) \f$.
1878 
1879 
1880 
1881  \li This function returns the vectors \a poly_in and \a poly_d
1882 
1883  containing the value of the \f$ n^th \f$ order Jacobi polynomial
1884 
1885  \f$ P^{\alpha,\beta}_n(z) \alpha > -1, \beta > -1 \f$ and its
1886 
1887  derivative at the \a np points in \a z[i]
1888 
1889 
1890 
1891  - If \a poly_in = NULL then only calculate derivatice
1892 
1893 
1894 
1895  - If \a polyd = NULL then only calculate polynomial
1896 
1897 
1898 
1899  - To calculate the polynomial this routine uses the recursion
1900 
1901  relationship (see appendix A ref [4]) :
1902 
1903  \f$ \begin{array}{rcl}
1904 
1905  P^{\alpha,\beta}_0(z) &=& 1 \\
1906 
1907  P^{\alpha,\beta}_1(z) &=& \frac{1}{2} [ \alpha-\beta+(\alpha+\beta+2)z] \\
1908 
1909  a^1_n P^{\alpha,\beta}_{n+1}(z) &=& (a^2_n + a^3_n z)
1910 
1911  P^{\alpha,\beta}_n(z) - a^4_n P^{\alpha,\beta}_{n-1}(z) \\
1912 
1913  a^1_n &=& 2(n+1)(n+\alpha + \beta + 1)(2n + \alpha + \beta) \\
1914 
1915  a^2_n &=& (2n + \alpha + \beta + 1)(\alpha^2 - \beta^2) \\
1916 
1917  a^3_n &=& (2n + \alpha + \beta)(2n + \alpha + \beta + 1)
1918 
1919  (2n + \alpha + \beta + 2) \\
1920 
1921  a^4_n &=& 2(n+\alpha)(n+\beta)(2n + \alpha + \beta + 2)
1922 
1923  \end{array} \f$
1924 
1925 
1926 
1927  - To calculate the derivative of the polynomial this routine uses
1928 
1929  the relationship (see appendix A ref [4]) :
1930 
1931  \f$ \begin{array}{rcl}
1932 
1933  b^1_n(z)\frac{d}{dz} P^{\alpha,\beta}_n(z)&=&b^2_n(z)P^{\alpha,\beta}_n(z)
1934 
1935  + b^3_n(z) P^{\alpha,\beta}_{n-1}(z) \hspace{2.2cm} \\
1936 
1937  b^1_n(z) &=& (2n+\alpha + \beta)(1-z^2) \\
1938 
1939  b^2_n(z) &=& n[\alpha - \beta - (2n+\alpha + \beta)z]\\
1940 
1941  b^3_n(z) &=& 2(n+\alpha)(n+\beta)
1942 
1943  \end{array} \f$
1944 
1945 
1946 
1947  - Note the derivative from this routine is only valid for -1 < \a z < 1.
1948 
1949  */
1950 
1951  void jacobfd(const int np, const double *z, double *poly_in, double *polyd,
1952 
1953  const int n, const double alpha, const double beta){
1954 
1955  register int i;
1956 
1957  double zero = 0.0, one = 1.0, two = 2.0;
1958 
1959 
1960 
1961  if(!np)
1962 
1963  return;
1964 
1965 
1966 
1967  if(n == 0){
1968 
1969  if(poly_in)
1970 
1971  for(i = 0; i < np; ++i)
1972 
1973  poly_in[i] = one;
1974 
1975  if(polyd)
1976 
1977  for(i = 0; i < np; ++i)
1978 
1979  polyd[i] = zero;
1980 
1981  }
1982 
1983  else if (n == 1){
1984 
1985  if(poly_in)
1986 
1987  for(i = 0; i < np; ++i)
1988 
1989  poly_in[i] = 0.5*(alpha - beta + (alpha + beta + two)*z[i]);
1990 
1991  if(polyd)
1992 
1993  for(i = 0; i < np; ++i)
1994 
1995  polyd[i] = 0.5*(alpha + beta + two);
1996 
1997  }
1998 
1999  else{
2000 
2001  register int k;
2002 
2003  double a1,a2,a3,a4;
2004 
2005  double two = 2.0, apb = alpha + beta;
2006 
2007  double *poly, *polyn1,*polyn2;
2008 
2009 
2010 
2011  if(poly_in){ // switch for case of no poynomial function return
2012 
2013  polyn1 = (double *)malloc(2*np*sizeof(double));
2014 
2015  polyn2 = polyn1+np;
2016 
2017  poly = poly_in;
2018 
2019  }
2020 
2021  else{
2022 
2023  polyn1 = (double *)malloc(3*np*sizeof(double));
2024 
2025  polyn2 = polyn1+np;
2026 
2027  poly = polyn2+np;
2028 
2029  }
2030 
2031 
2032 
2033  for(i = 0; i < np; ++i){
2034 
2035  polyn2[i] = one;
2036 
2037  polyn1[i] = 0.5*(alpha - beta + (alpha + beta + two)*z[i]);
2038 
2039  }
2040 
2041 
2042 
2043  for(k = 2; k <= n; ++k){
2044 
2045  a1 = two*k*(k + apb)*(two*k + apb - two);
2046 
2047  a2 = (two*k + apb - one)*(alpha*alpha - beta*beta);
2048 
2049  a3 = (two*k + apb - two)*(two*k + apb - one)*(two*k + apb);
2050 
2051  a4 = two*(k + alpha - one)*(k + beta - one)*(two*k + apb);
2052 
2053 
2054 
2055  a2 /= a1;
2056 
2057  a3 /= a1;
2058 
2059  a4 /= a1;
2060 
2061 
2062 
2063  for(i = 0; i < np; ++i){
2064 
2065  poly [i] = (a2 + a3*z[i])*polyn1[i] - a4*polyn2[i];
2066 
2067  polyn2[i] = polyn1[i];
2068 
2069  polyn1[i] = poly [i];
2070 
2071  }
2072 
2073  }
2074 
2075 
2076 
2077  if(polyd){
2078 
2079  a1 = n*(alpha - beta);
2080 
2081  a2 = n*(two*n + alpha + beta);
2082 
2083  a3 = two*(n + alpha)*(n + beta);
2084 
2085  a4 = (two*n + alpha + beta);
2086 
2087  a1 /= a4; a2 /= a4; a3 /= a4;
2088 
2089 
2090 
2091  // note polyn2 points to polyn1 at end of poly iterations
2092 
2093  for(i = 0; i < np; ++i){
2094 
2095  polyd[i] = (a1- a2*z[i])*poly[i] + a3*polyn2[i];
2096 
2097  polyd[i] /= (one - z[i]*z[i]);
2098 
2099  }
2100 
2101  }
2102 
2103 
2104 
2105  free(polyn1);
2106 
2107  }
2108 
2109 
2110 
2111  return;
2112 
2113  }
2114 
2115 
2116 
2117 
2118 
2119  /**
2120 
2121  \brief Calculate the derivative of Jacobi polynomials
2122 
2123 
2124 
2125  \li Generates a vector \a poly of values of the derivative of the
2126 
2127  \a n th order Jacobi polynomial \f$ P^(\alpha,\beta)_n(z)\f$ at the
2128 
2129  \a np points \a z.
2130 
2131 
2132 
2133  \li To do this we have used the relation
2134 
2135  \n
2136 
2137  \f$ \frac{d}{dz} P^{\alpha,\beta}_n(z)
2138 
2139  = \frac{1}{2} (\alpha + \beta + n + 1) P^{\alpha,\beta}_n(z) \f$
2140 
2141 
2142 
2143  \li This formulation is valid for \f$ -1 \leq z \leq 1 \f$
2144 
2145 
2146 
2147  */
2148 
2149 
2150 
2151  void jacobd(const int np, const double *z, double *polyd, const int n,
2152 
2153  const double alpha, const double beta)
2154 
2155  {
2156 
2157  register int i;
2158 
2159  double one = 1.0;
2160 
2161  if(n == 0)
2162 
2163  for(i = 0; i < np; ++i) polyd[i] = 0.0;
2164 
2165  else{
2166 
2167  //jacobf(np,z,polyd,n-1,alpha+one,beta+one);
2168 
2169  jacobfd(np,z,polyd,NULL,n-1,alpha+one,beta+one);
2170 
2171  for(i = 0; i < np; ++i) polyd[i] *= 0.5*(alpha + beta + (double)n + one);
2172 
2173  }
2174 
2175  return;
2176 
2177  }
2178 
2179 
2180 
2181 
2182 
2183  /**
2184 
2185  \brief Calculate the Gamma function , \f$ \Gamma(n)\f$, for integer
2186 
2187  values and halves.
2188 
2189 
2190 
2191  Determine the value of \f$\Gamma(n)\f$ using:
2192 
2193 
2194 
2195  \f$ \Gamma(n) = (n-1)! \mbox{ or } \Gamma(n+1/2) = (n-1/2)\Gamma(n-1/2)\f$
2196 
2197 
2198 
2199  where \f$ \Gamma(1/2) = \sqrt(\pi)\f$
2200 
2201  */
2202 
2203 
2204 
2205  double gammaF(const double x){
2206 
2207  double gamma = 1.0;
2208 
2209 
2210 
2211  if (x == -0.5) gamma = -2.0*sqrt(M_PI);
2212 
2213  else if (!x) return gamma;
2214 
2215  else if ((x-(int)x) == 0.5){
2216 
2217  int n = (int) x;
2218 
2219  double tmp = x;
2220 
2221 
2222 
2223  gamma = sqrt(M_PI);
2224 
2225  while(n--){
2226 
2227  tmp -= 1.0;
2228 
2229  gamma *= tmp;
2230 
2231  }
2232 
2233  }
2234 
2235  else if ((x-(int)x) == 0.0){
2236 
2237  int n = (int) x;
2238 
2239  double tmp = x;
2240 
2241 
2242 
2243  while(--n){
2244 
2245  tmp -= 1.0;
2246 
2247  gamma *= tmp;
2248 
2249  }
2250 
2251  }
2252 
2253  else
2254 
2255  fprintf(stderr,"%lf is not of integer or half order\n",x);
2256 
2257  return gamma;
2258 
2259  }
2260 
2261 
2262 
2263  /**
2264 
2265  \brief Calculate the \a n zeros, \a z, of the Jacobi polynomial, i.e.
2266 
2267  \f$ P_n^{\alpha,\beta}(z) = 0 \f$
2268 
2269 
2270 
2271  This routine is only value for \f$( \alpha > -1, \beta > -1)\f$
2272 
2273  and uses polynomial deflation in a Newton iteration
2274 
2275  */
2276 
2277 
2278 
2279  static void Jacobz(const int n, double *z, const double alpha,
2280 
2281  const double beta){
2282 
2283  register int i,j,k;
2284 
2285  double dth = M_PI/(2.0*(double)n);
2286 
2287  double poly,pder,rlast=0.0;
2288 
2289  double sum,delr,r;
2290 
2291  double one = 1.0, two = 2.0;
2292 
2293 
2294 
2295  if(!n)
2296 
2297  return;
2298 
2299 
2300 
2301  for(k = 0; k < n; ++k){
2302 
2303  r = -cos((two*(double)k + one) * dth);
2304 
2305  if(k) r = 0.5*(r + rlast);
2306 
2307 
2308 
2309  for(j = 1; j < STOP; ++j){
2310 
2311  jacobfd(1,&r,&poly, &pder, n, alpha, beta);
2312 
2313 
2314 
2315  for(i = 0, sum = 0.0; i < k; ++i) sum += one/(r - z[i]);
2316 
2317 
2318 
2319  delr = -poly / (pder - sum * poly);
2320 
2321  r += delr;
2322 
2323  if( fabs(delr) < EPS ) break;
2324 
2325  }
2326 
2327  z[k] = r;
2328 
2329  rlast = r;
2330 
2331  }
2332 
2333  return;
2334 
2335  }
2336 
2337 
2338 
2339 
2340 
2341  /**
2342 
2343  \brief Zero and Weight determination through the eigenvalues and eigenvectors of a tridiagonal
2344 
2345  matrix from the three term recurrence relationship.
2346 
2347 
2348 
2349  Set up a symmetric tridiagonal matrix
2350 
2351 
2352 
2353  \f$ \left [ \begin{array}{ccccc}
2354 
2355  a[0] & b[0] & & & \\
2356 
2357  b[0] & a[1] & b[1] & & \\
2358 
2359  0 & \ddots & \ddots & \ddots & \\
2360 
2361  & & \ddots & \ddots & b[n-2] \\
2362 
2363  & & & b[n-2] & a[n-1] \end{array} \right ] \f$
2364 
2365 
2366 
2367  Where the coefficients a[n], b[n] come from the recurrence relation
2368 
2369 
2370 
2371  \f$ b_j p_j(z) = (z - a_j ) p_{j-1}(z) - b_{j-1} p_{j-2}(z) \f$
2372 
2373 
2374 
2375  where \f$ j=n+1\f$ and \f$p_j(z)\f$ are the Jacobi (normalized)
2376 
2377  orthogonal polynomials \f$ \alpha,\beta > -1\f$( integer values and
2378 
2379  halves). Since the polynomials are orthonormalized, the tridiagonal
2380 
2381  matrix is guaranteed to be symmetric. The eigenvalues of this
2382 
2383  matrix are the zeros of the Jacobi polynomial.
2384 
2385  */
2386 
2387 
2388 
2389  void JacZeros(const int n, double *a, double*b, const double alpha,
2390 
2391  const double beta){
2392 
2393 
2394 
2395  int i,j;
2396 
2397  RecCoeff(n,a,b,alpha,beta);
2398 
2399 
2400 
2401  double **z = new double*[n];
2402 
2403  for(i = 0; i < n; i++)
2404 
2405  {
2406 
2407  z[i] = new double[n];
2408 
2409  for(j = 0; j < n; j++)
2410 
2411  {
2412 
2413  z[i][j] = 0.0;
2414 
2415  }
2416 
2417  }
2418 
2419  for(i = 0; i < n; i++)
2420 
2421  {
2422 
2423  z[i][i] = 1.0;
2424 
2425  }
2426 
2427 
2428 
2429  // find eigenvalues and eigenvectors
2430 
2431  TriQL(n, a, b,z);
2432 
2433 
2434 
2435  return;
2436 
2437  }
2438 
2439 
2440 
2441  /**
2442 
2443  \brief The routine finds the recurrence coefficients \a a and
2444 
2445  \a b of the orthogonal polynomials
2446 
2447  */
2448 
2449  static void RecCoeff(const int n, double *a, double *b,const double alpha,
2450 
2451  const double beta){
2452 
2453 
2454 
2455  int i;
2456 
2457  double apb, apbi,a2b2;
2458 
2459 
2460 
2461  if(!n)
2462 
2463  return;
2464 
2465 
2466 
2467  // generate normalised terms
2468 
2469  apb = alpha + beta;
2470 
2471  apbi = 2.0 + apb;
2472 
2473 
2474 
2475  b[0] = pow(2.0,apb+1.0)*gammaF(alpha+1.0)*gammaF(beta+1.0)/gammaF(apbi); //MuZero
2476 
2477  a[0] = (beta-alpha)/apbi;
2478 
2479  b[1] = (4.0*(1.0+alpha)*(1.0+beta)/((apbi+1.0)*apbi*apbi));
2480 
2481 
2482 
2483  a2b2 = beta*beta-alpha*alpha;
2484 
2485 
2486 
2487  for(i = 1; i < n-1; i++){
2488 
2489  apbi = 2.0*(i+1) + apb;
2490 
2491  a[i] = a2b2/((apbi-2.0)*apbi);
2492 
2493  b[i+1] = (4.0*(i+1)*(i+1+alpha)*(i+1+beta)*(i+1+apb)/
2494 
2495  ((apbi*apbi-1)*apbi*apbi));
2496 
2497  }
2498 
2499 
2500 
2501  apbi = 2.0*n + apb;
2502 
2503  a[n-1] = a2b2/((apbi-2.0)*apbi);
2504 
2505 
2506 
2507  }
2508 
2509 
2510 
2511 
2512 
2513  /** \brief QL algorithm for symmetric tridiagonal matrix
2514 
2515 
2516 
2517  This subroutine is a translation of an algol procedure,
2518 
2519  num. math. \b 12, 377-383(1968) by martin and wilkinson, as modified
2520 
2521  in num. math. \b 15, 450(1970) by dubrulle. Handbook for
2522 
2523  auto. comp., vol.ii-linear algebra, 241-248(1971). This is a
2524 
2525  modified version from numerical recipes.
2526 
2527 
2528 
2529  This subroutine finds the eigenvalues and first components of the
2530 
2531  eigenvectors of a symmetric tridiagonal matrix by the implicit QL
2532 
2533  method.
2534 
2535 
2536 
2537  on input:
2538 
2539  - n is the order of the matrix;
2540 
2541  - d contains the diagonal elements of the input matrix;
2542 
2543  - e contains the subdiagonal elements of the input matrix
2544 
2545  in its first n-2 positions.
2546 
2547  - z is the n by n identity matrix
2548 
2549 
2550 
2551  on output:
2552 
2553 
2554 
2555  - d contains the eigenvalues in ascending order.
2556 
2557  - e contains the weight values - modifications of the first component
2558 
2559  of normalised eigenvectors
2560 
2561  */
2562 
2563 
2564 
2565  static void TriQL(const int n, double *d,double *e, double **z){
2566 
2567  int m,l,iter,i,k;
2568 
2569  double s,r,p,g,f,dd,c,b;
2570 
2571 
2572 
2573  double MuZero = e[0];
2574 
2575 
2576 
2577  // Renumber the elements of e
2578 
2579  for(i = 0; i < n-1; i++)
2580 
2581  {
2582 
2583  e[i] = sqrt(e[i+1]);
2584 
2585  }
2586 
2587  e[n-1] = 0.0;
2588 
2589 
2590 
2591 
2592 
2593  for (l=0;l<n;l++) {
2594 
2595  iter=0;
2596 
2597  do {
2598 
2599  for (m=l;m<n-1;m++) {
2600 
2601  dd=fabs(d[m])+fabs(d[m+1]);
2602 
2603  if (fabs(e[m])+dd == dd) break;
2604 
2605  }
2606 
2607  if (m != l) {
2608 
2609  if (iter++ == STOP){
2610 
2611  fprintf(stderr,"triQL: Too many iterations in TQLI");
2612 
2613  exit(1);
2614 
2615  }
2616 
2617  g=(d[l+1]-d[l])/(2.0*e[l]);
2618 
2619  r=sqrt((g*g)+1.0);
2620 
2621  g=d[m]-d[l]+e[l]/(g+sign(r,g));
2622 
2623  s=c=1.0;
2624 
2625  p=0.0;
2626 
2627  for (i=m-1;i>=l;i--) {
2628 
2629  f=s*e[i];
2630 
2631  b=c*e[i];
2632 
2633  if (fabs(f) >= fabs(g)) {
2634 
2635  c=g/f;
2636 
2637  r=sqrt((c*c)+1.0);
2638 
2639  e[i+1]=f*r;
2640 
2641  c *= (s=1.0/r);
2642 
2643  } else {
2644 
2645  s=f/g;
2646 
2647  r=sqrt((s*s)+1.0);
2648 
2649  e[i+1]=g*r;
2650 
2651  s *= (c=1.0/r);
2652 
2653  }
2654 
2655  g=d[i+1]-p;
2656 
2657  r=(d[i]-g)*s+2.0*c*b;
2658 
2659  p=s*r;
2660 
2661  d[i+1]=g+p;
2662 
2663  g=c*r-b;
2664 
2665 
2666 
2667  // Calculate the eigenvectors
2668 
2669  for(k = 0; k < n; k++)
2670 
2671  {
2672 
2673  f = z[k][i+1];
2674 
2675  z[k][i+1] = s*z[k][i] + c*f;
2676 
2677  z[k][i] = c*z[k][i] - s*f;
2678 
2679  }
2680 
2681 
2682 
2683  }
2684 
2685  d[l]=d[l]-p;
2686 
2687  e[l]=g;
2688 
2689  e[m]=0.0;
2690 
2691  }
2692 
2693  } while (m != l);
2694 
2695  }
2696 
2697 
2698 
2699  // order eigenvalues
2700 
2701  // Since we only need the first component of the eigenvectors
2702 
2703  // to calcualte the weight, we only swap the first components
2704 
2705  for(i = 0; i < n-1; ++i){
2706 
2707  k = i;
2708 
2709  p = d[i];
2710 
2711  for(l = i+1; l < n; ++l)
2712 
2713  if (d[l] < p) {
2714 
2715  k = l;
2716 
2717  p = d[l];
2718 
2719  }
2720 
2721  d[k] = d[i];
2722 
2723  d[i] = p;
2724 
2725 
2726 
2727  double temp = z[0][k];
2728 
2729  z[0][k] = z[0][i];
2730 
2731  z[0][i] = temp;
2732 
2733  }
2734 
2735 
2736 
2737  // Calculate the weights
2738 
2739  for(i =0 ; i < n; i++)
2740 
2741  {
2742 
2743  e[i] = MuZero*z[0][i]*z[0][i];
2744 
2745  }
2746 
2747  }
2748 
2749 
2750 
2751  /**
2752 
2753  \brief Calcualtes the Jacobi-kronrod matrix by determining the
2754 
2755  \a a and \b coefficients.
2756 
2757 
2758 
2759  The first \a 3n+1 coefficients are already known
2760 
2761 
2762 
2763  For more information refer to:
2764 
2765  "Dirk P. Laurie, Calcualtion of Gauss-Kronrod quadrature rules"
2766 
2767  */
2768 
2769  void JKMatrix(int n, double *a, double *b)
2770 
2771  {
2772 
2773  int i,j,k,m;
2774 
2775  // Working storage
2776 
2777  int size = (int)floor(n/2.0)+2;
2778 
2779  double *s = new double[size];
2780 
2781  double *t = new double[size];
2782 
2783 
2784 
2785  // Initialize s and t to zero
2786 
2787  for(i = 0; i < size; i++)
2788 
2789  {
2790 
2791  s[i] = 0.0;
2792 
2793  t[i] = 0.0;
2794 
2795  }
2796 
2797 
2798 
2799  t[1] = b[n+1];
2800 
2801  for(m = 0; m <= n-2; m++)
2802 
2803  {
2804 
2805  double u = 0.0;
2806 
2807  for(k = (int)floor((m+1)/2.0); k >= 0; k--)
2808 
2809  {
2810 
2811  int l = m-k;
2812 
2813  u = u+(a[k+n+1]-a[l])*t[k+1] + b[k+n+1]*s[k] - b[l]*s[k+1];
2814 
2815  s[k+1] = u;
2816 
2817  }
2818 
2819 
2820 
2821  // Swap the contents of s and t
2822 
2823  double *hold = s;
2824 
2825  s = t;
2826 
2827  t = hold;
2828 
2829  }
2830 
2831 
2832 
2833 
2834 
2835  for(j = (int)floor(n/2.0); j >= 0; j--)
2836 
2837  {
2838 
2839  s[j+1] = s[j];
2840 
2841  }
2842 
2843 
2844 
2845  for(m = n-1; m <= 2*n-3; m++)
2846 
2847  {
2848 
2849  double u = 0;
2850 
2851  for(k = m+1-n; k <= floor((m-1)/2.0); k++)
2852 
2853  {
2854 
2855  int l = m-k;
2856 
2857  j = n-1-l;
2858 
2859  u = u-(a[k+n+1]-a[l])*t[j+1] - b[k+n+1]*s[j+1] + b[l]*s[j+2];
2860 
2861  s[j+1] = u;
2862 
2863  }
2864 
2865 
2866 
2867  if(m%2 == 0)
2868 
2869  {
2870 
2871  k = m/2;
2872 
2873  a[k+n+1] = a[k] + (s[j+1]-b[k+n+1]*s[j+2])/t[j+2];
2874 
2875 
2876 
2877  }else
2878 
2879  {
2880 
2881  k = (m+1)/2;
2882 
2883  b[k+n+1] = s[j+1]/s[j+2];
2884 
2885  }
2886 
2887 
2888 
2889 
2890 
2891  // Swap the contents of s and t
2892 
2893  double *hold = s;
2894 
2895  s = t;
2896 
2897  t = hold;
2898 
2899  }
2900 
2901 
2902 
2903  a[2*n ] = a[n-1]-b[2*n]*s[1]/t[1];
2904 
2905 
2906 
2907  }
2908 
2909 
2910 
2911  /**
2912 
2913  \brief
2914 
2915  Given a weight function \f$w(t)\f$ through the first \a n+1
2916 
2917  coefficients \a a and \a b of its orthogonal polynomials
2918 
2919  this routine generates the first \a n recurrence coefficients for the orthogonal
2920 
2921  polynomials relative to the modified weight function \f$(t-z)w(t)\f$.
2922 
2923 
2924 
2925  The result will be placed in the array \a a0 and \a b0.
2926 
2927  */
2928 
2929 
2930 
2931  void chri1(int n, double* a, double* b, double* a0,
2932 
2933  double* b0,double z)
2934 
2935  {
2936 
2937 
2938 
2939  double q = ceil(3.0*n/2);
2940 
2941  int size = (int)q+1;
2942 
2943  if(size < n+1)
2944 
2945  {
2946 
2947  fprintf(stderr,"input arrays a and b are too short\n");
2948 
2949  }
2950 
2951  double* r = new double[n+1];
2952 
2953  r[0] = z - a[0];
2954 
2955  r[1] = z - a[1] - b[1]/r[0];
2956 
2957  a0[0] = a[1] + r[1] - r[0];
2958 
2959  b0[0] = -r[0]*b[0];
2960 
2961 
2962 
2963  if(n == 1)
2964 
2965  {
2966 
2967  return;
2968 
2969  }
2970 
2971  int k = 0;
2972 
2973  for(k = 1; k < n; k++)
2974 
2975  {
2976 
2977  r[k+1] = z - a[k+1] - b[k+1]/r[k];
2978 
2979  a0[k] = a[k+1] + r[k+1] - r[k];
2980 
2981  b0[k] = b[k] * r[k]/r[k-1];
2982 
2983  }
2984 
2985 
2986  }
2987 
2988  /**
2989 
2990  \brief
2991 
2992  Calcualte the bessel function of the first kind with complex double input y.
2993  Taken from Numerical Recipies in C
2994 
2995  Returns a complex double
2996  */
2997 
2998 
2999  std::complex<Nektar::NekDouble> ImagBesselComp(int n,std::complex<Nektar::NekDouble> y)
3000  {
3001  std::complex<Nektar::NekDouble> z (1.0,0.0);
3002  std::complex<Nektar::NekDouble> zbes (1.0,0.0);
3003  std::complex<Nektar::NekDouble> zarg;
3004  Nektar::NekDouble tol = 1e-15;
3005  int maxit = 10000;
3006  int i = 1;
3007 
3008  zarg = -0.25*y*y;
3009 
3010  while (abs(z) > tol && i <= maxit){
3011  z = z*(1.0/i/(i+n)*zarg);
3012  if (abs(z) <= tol) break;
3013  zbes = zbes + z;
3014  i++;
3015  }
3016  zarg = 0.5*y;
3017  for (i=1;i<=n;i++){
3018  zbes = zbes*zarg;
3019  }
3020  return zbes;
3021 
3022  }
3023 
3024 
3025 
3026 } // end of namespace
3027 
3028 
3029 
void Dgj(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix and its transpose associated.
Definition: Polylib.cpp:937
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:27
static void RecCoeff(const int, double *, double *, const double, const double)
The routine finds the recurrence coefficients a and.
Definition: Polylib.cpp:2449
void zwgrjp(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=1.
Definition: Polylib.cpp:239
double hgj(const int i, const double z, const double *zgj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through.
Definition: Polylib.cpp:1345
void Imgrjp(double *im, const double *zgrj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Radau-Jacobi points.
Definition: Polylib.cpp:1775
void zwglj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Lobatto-Jacobi zeros and weights with end point at z=-1,1.
Definition: Polylib.cpp:311
double hgrjp(const int i, const double z, const double *zgrj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the.
Definition: Polylib.cpp:1505
void zwlk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Lobatto-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:711
void JacZeros(const int n, double *a, double *b, const double alpha, const double beta)
Zero and Weight determination through the eigenvalues and eigenvectors of a tridiagonal.
Definition: Polylib.cpp:2389
void chri1(int, double *, double *, double *, double *, double)
Given a weight function through the first n+1.
Definition: Polylib.cpp:2931
void jacobd(const int np, const double *z, double *polyd, const int n, const double alpha, const double beta)
Calculate the derivative of Jacobi polynomials.
Definition: Polylib.cpp:2151
#define jacobz(n, z, alpha, beta)
zero determination using Newton iteration with polynomial deflation
Definition: Polylib.cpp:45
#define STOP
Maximum number of iterations in polynomial defalation routine Jacobz.
Definition: Polylib.cpp:19
void Imglj(double *im, const double *zglj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Lobatto-Jacobi points.
Definition: Polylib.cpp:1837
void JKMatrix(int, double *, double *)
Calcualtes the Jacobi-kronrod matrix by determining the.
Definition: Polylib.cpp:2769
void Dgrjm(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix and its transpose associated.
Definition: Polylib.cpp:1023
void zwgk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:385
void Imgrjm(double *im, const double *zgrj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Radau-Jacobi points.
Definition: Polylib.cpp:1715
void Dgrjp(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix associated with the.
Definition: Polylib.cpp:1121
double NekDouble
void zwgrjm(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=-1.
Definition: Polylib.cpp:163
The namespace associated with the the Polylib library (Polylib introduction)
Definition: Polylib.cpp:31
static void TriQL(const int, double *, double *, double **)
QL algorithm for symmetric tridiagonal matrix.
Definition: Polylib.cpp:2565
static void Jacobz(const int n, double *z, const double alpha, const double beta)
Calculate the n zeros, z, of the Jacobi polynomial, i.e.
Definition: Polylib.cpp:2279
void zwrk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Radau-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:519
double hglj(const int i, const double z, const double *zglj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the.
Definition: Polylib.cpp:1589
std::complex< Nektar::NekDouble > ImagBesselComp(int n, std::complex< Nektar::NekDouble > y)
Calcualte the bessel function of the first kind with complex double input y. Taken from Numerical Rec...
Definition: Polylib.cpp:2999
double hgrjm(const int i, const double z, const double *zgrj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the.
Definition: Polylib.cpp:1421
#define EPS
Precision tolerance for two points to be similar.
Definition: Polylib.cpp:23
void Dglj(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix associated with the.
Definition: Polylib.cpp:1221
void Imgj(double *im, const double *zgj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Jacobi points to an.
Definition: Polylib.cpp:1657
void jacobfd(const int np, const double *z, double *poly_in, double *polyd, const int n, const double alpha, const double beta)
Routine to calculate Jacobi polynomials, , and their first derivative, .
Definition: Polylib.cpp:1951
double gammaF(const double)
Calculate the Gamma function , , for integer.
Definition: Polylib.cpp:2205
void zwgj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Jacobi zeros and weights.
Definition: Polylib.cpp:107