Intrepid
Intrepid_BurkardtRulesDef.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid Package
5 // Copyright (2007) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Pavel Bochev (pbboche@sandia.gov)
38 // Denis Ridzal (dridzal@sandia.gov), or
39 // Kara Peterson (kjpeter@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
51 #ifdef _MSC_VER
52 #include "winmath.h"
53 #endif
54 
55 //# include <cstdlib>
56 //# include <iomanip>
57 //# include <iostream>
58 //# include <cmath>
59 //# include <ctime>
60 //# include <Teuchos_ScalarTraitsDecl.hpp>
61 
62 namespace Intrepid {
63 
64 //****************************************************************************
65 template<class Scalar>
66 void IntrepidBurkardtRules::chebyshev1_compute ( int n, Scalar x[], Scalar w[] )
67 //****************************************************************************
68 //
69 // Purpose:
70 //
71 // CHEBYSHEV1_COMPUTE computes a Chebyshev type 1 quadrature rule.
72 //
73 // Discussion:
74 //
75 // The integral:
76 //
77 // Integral ( -1 <= X <= 1 ) F(X) / sqrt ( 1 - x^2 ) dX
78 //
79 // The quadrature rule:
80 //
81 // Sum ( 1 <= I <= N ) W(I) * F ( X(I) )
82 //
83 // Licensing:
84 //
85 // This code is distributed under the GNU LGPL license.
86 //
87 // Modified:
88 //
89 // 13 June 2009
90 //
91 // Author:
92 //
93 // John Burkardt
94 //
95 // Reference:
96 //
97 // Philip Davis, Philip Rabinowitz,
98 // Methods of Numerical Integration,
99 // Second Edition,
100 // Dover, 2007,
101 // ISBN: 0486453391,
102 // LC: QA299.3.D28.
103 //
104 // Parameters:
105 //
106 // Input, int N, the order.
107 // 1 <= N.
108 //
109 // Output, Scalar X[N], the abscissas.
110 //
111 // Output, Scalar W[N], the weights.
112 //
113 {
114  if (n<1) {
115  std::cerr << "\n";
116  std::cerr << "CHEBYSHEV1_COMPUTE - Fatal error!\n";
117  std::cerr << " Illegal value of N = " << n << "\n";
118  std::exit (1);
119  }
120 
121  for (int i=0;i<n;i++) {
122  w[i] = M_PI/(Scalar)(n);
123  }
124  for (int i=0;i<n;i++) {
125  x[i] = std::cos(M_PI*(Scalar)(2*n-1-2*i)/(Scalar)(2*n));
126  }
127  if ((n%2)==1) {
128  x[(n-1)/2] = 0.0;
129  }
130 
131  return;
132 }
133 
134 //****************************************************************************
135 template<class Scalar>
137 //****************************************************************************
138 //
139 // Purpose:
140 //
141 // CHEBYSHEV1_COMPUTE_POINTS computes Chebyshev type 1 quadrature points.
142 //
143 // Licensing:
144 //
145 // This code is distributed under the GNU LGPL license.
146 //
147 // Modified:
148 //
149 // 13 June 2009
150 //
151 // Author:
152 //
153 // John Burkardt
154 //
155 // Reference:
156 //
157 // Philip Davis, Philip Rabinowitz,
158 // Methods of Numerical Integration,
159 // Second Edition,
160 // Dover, 2007,
161 // ISBN: 0486453391,
162 // LC: QA299.3.D28.
163 //
164 // Parameters:
165 //
166 // Input, int N, the order.
167 // 1 <= N.
168 //
169 // Output, Scalar X[N], the abscissas.
170 //
171 {
172  if (n<1) {
173  std::cerr << "\n";
174  std::cerr << "CHEBYSHEV1_COMPUTE_POINTS - Fatal error!\n";
175  std::cerr << " Illegal value of N = " << n << "\n";
176  std::exit(1);
177  }
178 
179  for (int i=0;i<n;i++) {
180  x[i] = std::cos(M_PI*(Scalar)(2*n-1-2*i)/(Scalar)(2*n));
181  }
182  if ((n%2)==1) {
183  x[(n-1)/2] = 0.0;
184  }
185 
186  return;
187 }
188 
189 //****************************************************************************
190 template<class Scalar>
192 //****************************************************************************
193 //
194 // Purpose:
195 //
196 // CHEBYSHEV1_COMPUTE_WEIGHTS computes Chebyshev type 1 quadrature weights.
197 //
198 // Licensing:
199 //
200 // This code is distributed under the GNU LGPL license.
201 //
202 // Modified:
203 //
204 // 13 June 2009
205 //
206 // Author:
207 //
208 // John Burkardt
209 //
210 // Reference:
211 //
212 // Philip Davis, Philip Rabinowitz,
213 // Methods of Numerical Integration,
214 // Second Edition,
215 // Dover, 2007,
216 // ISBN: 0486453391,
217 // LC: QA299.3.D28.
218 //
219 // Parameters:
220 //
221 // Input, int N, the order.
222 // 1 <= N.
223 //
224 // Output, Scalar W[N], the weights.
225 //
226 {
227  if (n<1) {
228  std::cerr << "\n";
229  std::cerr << "CHEBYSHEV1_COMPUTE_WEIGHTS - Fatal error!\n";
230  std::cerr << " Illegal value of N = " << n << "\n";
231  std::exit(1);
232  }
233 
234  for (int i=0;i<n;i++) {
235  w[i] = M_PI/(Scalar)n;
236  }
237 
238  return;
239 }
240 
241 //****************************************************************************
242 template<class Scalar>
243 void IntrepidBurkardtRules::chebyshev2_compute ( int n, Scalar x[], Scalar w[] )
244 //****************************************************************************
245 //
246 // Purpose:
247 //
248 // CHEBYSHEV2_COMPUTE computes a Chebyshev type 2 quadrature rule.
249 //
250 // Discussion:
251 //
252 // The integral:
253 //
254 // integral ( -1 <= x <= 1 ) f(x) sqrt ( 1 - x^2 ) dx
255 //
256 // The quadrature rule:
257 //
258 // sum ( 1 <= i <= n ) w(i) * f ( x(i) )
259 //
260 // Licensing:
261 //
262 // This code is distributed under the GNU LGPL license.
263 //
264 // Modified:
265 //
266 // 13 June 2009
267 //
268 // Author:
269 //
270 // John Burkardt
271 //
272 // Reference:
273 //
274 // Philip Davis, Philip Rabinowitz,
275 // Methods of Numerical Integration,
276 // Second Edition,
277 // Dover, 2007,
278 // ISBN: 0486453391,
279 // LC: QA299.3.D28.
280 //
281 // Parameters:
282 //
283 // Input, int N, the order.
284 // 1 <= N.
285 //
286 // Output, Scalar X[N], the abscissas.
287 //
288 // Output, Scalar W[N], the weights.
289 //
290 {
291  Scalar angle;
292 
293  if (n<1) {
294  std::cerr << "\n";
295  std::cerr << "CHEBYSHEV2_COMPUTE - Fatal error!\n";
296  std::cerr << " Illegal value of N = " << n << "\n";
297  std::exit(1);
298  }
299 
300  for (int i=0;i<n;i++) {
301  angle = M_PI*(Scalar)(n-i)/(Scalar)(n+1);
302  w[i] = M_PI/(Scalar)(n+1)*std::pow(std::sin(angle),2);
303  x[i] = std::cos(angle);
304  }
305 
306  if ((n%2)==1) {
307  x[(n-1)/2] = 0.0;
308  }
309 
310  return;
311 }
312 
313 //****************************************************************************
314 template<class Scalar>
316 //****************************************************************************
317 //
318 // Purpose:
319 //
320 // CHEBYSHEV2_COMPUTE_POINTS computes Chebyshev type 2 quadrature points.
321 //
322 // Licensing:
323 //
324 // This code is distributed under the GNU LGPL license.
325 //
326 // Modified:
327 //
328 // 13 June 2009
329 //
330 // Author:
331 //
332 // John Burkardt
333 //
334 // Reference:
335 //
336 // Philip Davis, Philip Rabinowitz,
337 // Methods of Numerical Integration,
338 // Second Edition,
339 // Dover, 2007,
340 // ISBN: 0486453391,
341 // LC: QA299.3.D28.
342 //
343 // Parameters:
344 //
345 // Input, int N, the order.
346 // 1 <= N.
347 //
348 // Output, Scalar X[N], the abscissas.
349 //
350 {
351  Scalar angle;
352 
353  if (n<1) {
354  std::cerr << "\n";
355  std::cerr << "CHEBYSHEV2_COMPUTE_POINTS - Fatal error!\n";
356  std::cerr << " Illegal value of N = " << n << "\n";
357  std::exit(1);
358  }
359 
360  for (int i=0;i<n;i++) {
361  angle = M_PI*(Scalar)(n-i)/(Scalar)(n+1);
362  x[i] = std::cos(angle);
363  }
364 
365  if ((n%2)==1) {
366  x[(n-1)/2] = 0.0;
367  }
368 
369  return;
370 }
371 
372 //****************************************************************************
373 template<class Scalar>
375 //****************************************************************************
376 //
377 // Purpose:
378 //
379 // CHEBYSHEV2_COMPUTE_WEIGHTS computes Chebyshev type 2 quadrature weights.
380 //
381 // Licensing:
382 //
383 // This code is distributed under the GNU LGPL license.
384 //
385 // Modified:
386 //
387 // 13 June 2009
388 //
389 // Author:
390 //
391 // John Burkardt
392 //
393 // Reference:
394 //
395 // Philip Davis, Philip Rabinowitz,
396 // Methods of Numerical Integration,
397 // Second Edition,
398 // Dover, 2007,
399 // ISBN: 0486453391,
400 // LC: QA299.3.D28.
401 //
402 // Parameters:
403 //
404 // Input, int N, the order.
405 // 1 <= N.
406 //
407 // Output, Scalar W[N], the weights.
408 //
409 {
410  Scalar angle;
411 
412  if (n<1) {
413  std::cerr << "\n";
414  std::cerr << "CHEBYSHEV2_COMPUTE_WEIGHTS - Fatal error!\n";
415  std::cerr << " Illegal value of N = " << n << "\n";
416  std::exit(1);
417  }
418 
419  for (int i=0;i<n;i++) {
420  angle = M_PI*(Scalar)(n-i)/(Scalar)(n+1);
421  w[i] = M_PI/(Scalar)(n+1)*std::pow(std::sin(angle),2);
422  }
423 
424  return;
425 }
426 
427 //****************************************************************************
428 template<class Scalar>
429 void IntrepidBurkardtRules::clenshaw_curtis_compute ( int n, Scalar x[], Scalar w[] )
430 //****************************************************************************
431 //
432 // Purpose:
433 //
434 // CLENSHAW_CURTIS_COMPUTE computes a Clenshaw Curtis quadrature rule.
435 //
436 // Discussion:
437 //
438 // The integral:
439 //
440 // Integral ( -1 <= X <= 1 ) F(X) dX
441 //
442 // The quadrature rule:
443 //
444 // Sum ( 1 <= I <= N ) W(I) * F ( X(I) )
445 //
446 // Licensing:
447 //
448 // This code is distributed under the GNU LGPL license.
449 //
450 // Modified:
451 //
452 // 19 March 2009
453 //
454 // Author:
455 //
456 // John Burkardt
457 //
458 // Parameters:
459 //
460 // Input, int N, the order.
461 // 1 <= N.
462 //
463 // Output, Scalar X[N], the abscissas.
464 //
465 // Output, Scalar W[N], the weights.
466 //
467 {
468  Scalar b, theta;
469  int i, j;
470 
471  if (n<1) {
472  std::cerr << "\n";
473  std::cerr << "CLENSHAW_CURTIS_COMPUTE - Fatal error!\n";
474  std::cerr << " Illegal value of N = " << n << "\n";
475  std::exit(1);
476  }
477  else if (n==1) {
478  x[0] = 0.0;
479  w[0] = 2.0;
480  }
481  else {
482  for (i=0;i<n;i++) {
483  x[i] = std::cos((Scalar)(n-1-i)*M_PI/(Scalar)(n-1));
484  }
485  x[0] = -1.0;
486  if ((n%2)==1) {
487  x[(n-1)/2] = 0.0;
488  }
489  x[n-1] = +1.0;
490 
491  for (i=0;i<n;i++) {
492  theta = (Scalar)i*M_PI/(Scalar)(n-1);
493 
494  w[i] = 1.0;
495 
496  for (j=1;j<=(n-1)/2;j++) {
497  if (2*j==(n-1)) {
498  b = 1.0;
499  }
500  else {
501  b = 2.0;
502  }
503 
504  w[i] = w[i]-b*std::cos(2.0*(Scalar)(j)*theta)/(Scalar)(4*j*j-1);
505  }
506  }
507 
508  w[0] = w[0]/(Scalar)(n-1);
509  for (i=1;i<n-1;i++) {
510  w[i] = 2.0*w[i]/(Scalar)(n-1);
511  }
512  w[n-1] = w[n-1]/(Scalar)(n-1);
513  }
514 
515  return;
516 }
517 
518 //****************************************************************************
519 template<class Scalar>
521 //****************************************************************************
522 //
523 // Purpose:
524 //
525 // CLENSHAW_CURTIS_COMPUTE_POINTS computes Clenshaw Curtis quadrature points.
526 //
527 // Discussion:
528 //
529 // Our convention is that the abscissas are numbered from left to right.
530 //
531 // This rule is defined on [-1,1].
532 //
533 // Licensing:
534 //
535 // This code is distributed under the GNU LGPL license.
536 //
537 // Modified:
538 //
539 // 13 June 2009
540 //
541 // Author:
542 //
543 // John Burkardt
544 //
545 // Parameters:
546 //
547 // Input, int N, the order.
548 //
549 // Output, Scalar X[N], the abscissas.
550 //
551 {
552  int index;
553 
554  if (n<1) {
555  std::cerr << "\n";
556  std::cerr << "CLENSHAW_CURTIS_COMPUTE_POINTS - Fatal error!\n";
557  std::cerr << " N < 1.\n";
558  std::exit(1);
559  }
560  else if (n==1) {
561  x[0] = 0.0;
562  }
563  else {
564  for (index=1;index<=n;index++) {
565  x[index-1] = std::cos((Scalar)(n-index)*M_PI/(Scalar)(n-1));
566  }
567  x[0] = -1.0;
568  if ((n%2)==1) {
569  x[(n-1)/2] = 0.0;
570  }
571  x[n-1] = +1.0;
572  }
573  return;
574 }
575 
576 //****************************************************************************
577 template<class Scalar>
579 //****************************************************************************
580 //
581 // Purpose:
582 //
583 // CLENSHAW_CURTIS_COMPUTE_WEIGHTS computes Clenshaw Curtis quadrature weights.
584 //
585 // Discussion:
586 //
587 // The user must preallocate space for the output array W.
588 //
589 // Licensing:
590 //
591 // This code is distributed under the GNU LGPL license.
592 //
593 // Modified:
594 //
595 // 13 June 2009
596 //
597 // Author:
598 //
599 // John Burkardt
600 //
601 // Reference:
602 //
603 // Charles Clenshaw, Alan Curtis,
604 // A Method for Numerical Integration on an Automatic Computer,
605 // Numerische Mathematik,
606 // Volume 2, Number 1, December 1960, pages 197-205.
607 //
608 // Parameters:
609 //
610 // Input, int N, the order.
611 //
612 // Output, Scalar W[N], the weights.
613 //
614 {
615  Scalar b, theta;
616  int i, j;
617 
618  if (n<1) {
619  std::cerr << "\n";
620  std::cerr << "CLENSHAW_CURTIS_COMPUTE_WEIGHTS - Fatal error!\n";
621  std::cerr << " N < 1.\n";
622  std::exit(1);
623  }
624  else if (n==1) {
625  w[0] = 2.0;
626  return;
627  }
628 
629  for (i=1;i<=n;i++) {
630  theta = (Scalar)(i-1)*M_PI/(Scalar)(n-1);
631 
632  w[i-1] = 1.0;
633 
634  for (j=1;j<=(n-1)/2;j++) {
635  if (2*j==(n-1)) {
636  b = 1.0;
637  }
638  else {
639  b = 2.0;
640  }
641 
642  w[i-1] = w[i-1]-b*std::cos(2.0*(Scalar)j*theta)/(Scalar)(4*j*j-1);
643  }
644  }
645 
646  w[0] = w[0]/(Scalar)(n-1);
647  for (i=1;i<n-1;i++) {
648  w[i] = 2.0*w[i]/(Scalar)(n-1);
649  }
650  w[n-1] = w[n-1]/(Scalar)(n-1);
651 
652  return;
653 }
654 
655 //****************************************************************************
656 template<class Scalar>
657 void IntrepidBurkardtRules::fejer2_compute ( int n, Scalar x[], Scalar w[] )
658 //****************************************************************************
659 //
660 // Purpose:
661 //
662 // FEJER2_COMPUTE computes a Fejer type 2 rule.
663 //
664 // Discussion:
665 //
666 // Our convention is that the abscissas are numbered from left to right.
667 //
668 // The rule is defined on [-1,1].
669 //
670 // Licensing:
671 //
672 // This code is distributed under the GNU LGPL license.
673 //
674 // Modified:
675 //
676 // 13 June 2009
677 //
678 // Author:
679 //
680 // John Burkardt
681 //
682 // Parameters:
683 //
684 // Input, int N, the order.
685 // 1 <= N.
686 //
687 // Output, Scalar X[N], the abscissas.
688 //
689 // Output, Scalar W[N], the weights.
690 //
691 {
692  Scalar p, theta;
693 
694  if (n<1) {
695  std::cerr << "\n";
696  std::cerr << "FEJER2_COMPUTE - Fatal error!\n";
697  std::cerr << " Illegal value of N = " << n << "\n";
698  std::exit(1);
699  }
700  else if (n==1) {
701  x[0] = 0.0;
702  w[0] = 2.0;
703  return;
704  }
705 
706  for (int i=0;i<n;i++) {
707  x[i] = std::cos((Scalar)(n-i)*M_PI/(Scalar)(n+1));
708  }
709  if ((n%2)==1) {
710  x[(n-1)/2] = 0.0;
711  }
712 
713  if (n==2) {
714  w[0] = 1.0;
715  w[1] = 1.0;
716  }
717  else {
718  for (int i=0;i<n;i++) {
719  theta = (Scalar)(n-i)*M_PI/(Scalar)(n+1);
720 
721  w[i] = 1.0;
722 
723  for (int j=1;j<=((n-1)/2);j++) {
724  w[i] = w[i]-2.0*std::cos(2.0*(Scalar)j*theta)/(Scalar)(4*j*j-1);
725  }
726  p = 2.0*(Scalar)(((n+1)/2))-1.0;
727  w[i] = w[i]-std::cos((p+1.0)*theta)/p;
728  }
729  for (int i=0;i<n;i++) {
730  w[i] = 2.0*w[i]/(Scalar)(n+1);
731  }
732  }
733  return;
734 }
735 
736 //****************************************************************************
737 template<class Scalar>
739 //****************************************************************************
740 //
741 // Purpose:
742 //
743 // FEJER2_COMPUTE_POINTS computes Fejer type 2 quadrature points.
744 //
745 // Discussion:
746 //
747 // Our convention is that the abscissas are numbered from left to right.
748 //
749 // The rule is defined on [-1,1].
750 //
751 // Licensing:
752 //
753 // This code is distributed under the GNU LGPL license.
754 //
755 // Modified:
756 //
757 // 13 June 2009
758 //
759 // Author:
760 //
761 // John Burkardt
762 //
763 // Parameters:
764 //
765 // Input, int N, the order.
766 // 1 <= N.
767 //
768 // Output, Scalar X[N], the abscissas.
769 //
770 {
771  int i;
772 
773  if (n<1) {
774  std::cerr << "\n";
775  std::cerr << "FEJER2_COMPUTE_POINTS - Fatal error!\n";
776  std::cerr << " N < 1.\n";
777  std::exit(1);
778  }
779  else if (n==1) {
780  x[0] = 0.0;
781  }
782  else {
783  for (i=1;i<=n;i++) {
784  x[i-1] = std::cos((Scalar)(n+1-i)*M_PI/(Scalar)(n+1));
785  }
786  if ((n%2)==1) {
787  x[(n-1)/2] = 0.0;
788  }
789  }
790  return;
791 }
792 
793 //****************************************************************************
794 template<class Scalar>
796 //****************************************************************************
797 //
798 // Purpose:
799 //
800 // FEJER2_COMPUTE_WEIGHTS computes Fejer type 2 quadrature weights.
801 //
802 // Discussion:
803 //
804 // The user must preallocate space for the output array W.
805 //
806 // Licensing:
807 //
808 // This code is distributed under the GNU LGPL license.
809 //
810 // Modified:
811 //
812 // 13 June 2009
813 //
814 // Author:
815 //
816 // John Burkardt
817 //
818 // Reference:
819 //
820 // Philip Davis, Philip Rabinowitz,
821 // Methods of Numerical Integration,
822 // Second Edition,
823 // Dover, 2007,
824 // ISBN: 0486453391,
825 // LC: QA299.3.D28.
826 //
827 // Walter Gautschi,
828 // Numerical Quadrature in the Presence of a Singularity,
829 // SIAM Journal on Numerical Analysis,
830 // Volume 4, Number 3, 1967, pages 357-362.
831 //
832 // Joerg Waldvogel,
833 // Fast Construction of the Fejer and Clenshaw-Curtis Quadrature Rules,
834 // BIT Numerical Mathematics,
835 // Volume 43, Number 1, 2003, pages 1-18.
836 //
837 // Parameters:
838 //
839 // Input, int N, the order.
840 //
841 // Output, Scalar W[N], the weights.
842 //
843 {
844  int i, j;
845  Scalar p, theta;
846 
847  if (n<1) {
848  std::cerr << "\n";
849  std::cerr << "FEJER2_COMPUTE_WEIGHTS - Fatal error!\n";
850  std::cerr << " N < 1.\n";
851  std::exit(1);
852  }
853  else if (n==1) {
854  w[0] = 2.0;
855  }
856  else if (n==2) {
857  w[0] = 1.0;
858  w[1] = 1.0;
859  }
860  else {
861  for (i=1;i<=n;i++) {
862  theta = (Scalar)(n+1-i)*M_PI/(Scalar)(n+1);
863 
864  w[i-1] = 1.0;
865 
866  for (j=1;j<=((n-1)/2);j++) {
867  w[i-1] = w[i-1]-2.0*std::cos(2.0*(Scalar)j*theta)/(Scalar)(4*j*j-1);
868  }
869  p = 2.0*(Scalar)(((n+1)/2))-1.0;
870  w[i-1] = w[i-1]-std::cos((p+1.0)*theta)/p;
871  }
872  for (i=0;i<n;i++) {
873  w[i] = 2.0*w[i]/(Scalar)(n+1);
874  }
875  }
876  return;
877 }
878 
879 //****************************************************************************
880 template<class Scalar>
881 void IntrepidBurkardtRules::hermite_compute ( int n, Scalar x[], Scalar w[] )
882 //****************************************************************************
883 //
884 // Purpose:
885 //
886 // HERMITE_COMPUTE computes a Gauss-Hermite quadrature rule.
887 //
888 // Discussion:
889 //
890 // The code uses an algorithm by Elhay and Kautsky.
891 //
892 // The abscissas are the zeros of the N-th order Hermite polynomial.
893 //
894 // The integral:
895 //
896 // integral ( -oo < x < +oo ) exp ( - x * x ) * f(x) dx
897 //
898 // The quadrature rule:
899 //
900 // sum ( 1 <= i <= n ) w(i) * f ( x(i) )
901 //
902 // Licensing:
903 //
904 // This code is distributed under the GNU LGPL license.
905 //
906 // Modified:
907 //
908 // 19 April 2011
909 //
910 // Author:
911 //
912 // Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky.
913 // C++ version by John Burkardt.
914 //
915 // Reference:
916 //
917 // Sylvan Elhay, Jaroslav Kautsky,
918 // Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of
919 // Interpolatory Quadrature,
920 // ACM Transactions on Mathematical Software,
921 // Volume 13, Number 4, December 1987, pages 399-415.
922 //
923 // Parameters:
924 //
925 // Input, int N, the number of abscissas.
926 //
927 // Output, Scalar X[N], the abscissas.
928 //
929 // Output, Scalar W[N], the weights.
930 //
931 {
932  Scalar* bj;
933 //
934 // Define the zero-th moment.
935 //
936  Scalar zemu = std::sqrt(M_PI);
937 //
938 // Define the Jacobi matrix.
939 //
940  bj = new Scalar[n];
941 
942  for (int i=0;i<n;i++) {
943  bj[i] = std::sqrt((Scalar)(i+1)/2.0);
944  }
945 
946  for (int i=0;i<n;i++) {
947  x[i] = 0.0;
948  }
949 
950  w[0] = std::sqrt (zemu);
951  for (int i=1;i<n;i++) {
952  w[i] = 0.0;
953  }
954 //
955 // Diagonalize the Jacobi matrix.
956 //
957  IntrepidBurkardtRules::imtqlx ( n, x, bj, w );
958 
959  for (int i=0;i<n;i++) {
960  w[i] = w[i]*w[i];
961  }
962 
963  // Ensure that zero is actually zero.
964  if (n%2) {
965  int ind = (int)((Scalar)n/2.0);
966  x[ind] = 0.0;
967  }
968 
969  delete [] bj;
970 
971  return;
972 }
973 
974 //****************************************************************************
975 template<class Scalar>
976 void IntrepidBurkardtRules::hermite_compute_points ( int order, Scalar x[] )
977 //****************************************************************************
978 //
979 // Purpose:
980 //
981 // HERMITE_COMPUTE_POINTS computes Hermite quadrature points.
982 //
983 // Licensing:
984 //
985 // This code is distributed under the GNU LGPL license.
986 //
987 // Modified:
988 //
989 // 13 June 2009
990 //
991 // Author:
992 //
993 // John Burkardt
994 //
995 // Parameters:
996 //
997 // Input, int ORDER, the order.
998 //
999 // Output, Scalar X[ORDER], the abscissas.
1000 //
1001 {
1002  Scalar *w; w = new Scalar[order];
1004  delete [] w;
1005 
1006  return;
1007 }
1008 
1009 //****************************************************************************
1010 template<class Scalar>
1012 //****************************************************************************
1013 //
1014 // Purpose:
1015 //
1016 // HERMITE_COMPUTE_WEIGHTS computes Hermite quadrature weights.
1017 //
1018 // Licensing:
1019 //
1020 // This code is distributed under the GNU LGPL license.
1021 //
1022 // Modified:
1023 //
1024 // 13 June 2009
1025 //
1026 // Author:
1027 //
1028 // John Burkardt
1029 //
1030 // Parameters:
1031 //
1032 // Input, int ORDER, the order.
1033 //
1034 // Output, Scalar W[ORDER], the weights.
1035 //
1036 {
1037  Scalar *x; x = new Scalar[order];
1039  delete [] x;
1040 
1041  return;
1042 }
1043 
1044 //****************************************************************************
1045 template<class Scalar>
1046 void IntrepidBurkardtRules::hermite_genz_keister_lookup ( int n, Scalar x[], Scalar w[] )
1047 //****************************************************************************
1048 //
1049 // Purpose:
1050 //
1051 // HERMITE_GENZ_KEISTER_LOOKUP looks up a Genz-Keister Hermite rule.
1052 //
1053 // Discussion:
1054 //
1055 // The integral:
1056 //
1057 // integral ( -oo <= x <= +oo ) f(x) exp ( - x * x ) dx
1058 //
1059 // The quadrature rule:
1060 //
1061 // sum ( 1 <= i <= n ) w(i) * f ( x(i) )
1062 //
1063 // A nested family of rules for the Hermite integration problem
1064 // was produced by Genz and Keister. The structure of the nested
1065 // family was denoted by 1+2+6+10+16, that is, it comprised rules
1066 // of successive orders O = 1, 3, 9, 19, and 35.
1067 //
1068 // The precisions of these rules are P = 1, 5, 15, 29, and 51.
1069 //
1070 // Licensing:
1071 //
1072 // This code is distributed under the GNU LGPL license.
1073 //
1074 // Modified:
1075 //
1076 // 07 June 2010
1077 //
1078 // Author:
1079 //
1080 // John Burkardt
1081 //
1082 // Reference:
1083 //
1084 // Alan Genz, Bradley Keister,
1085 // Fully symmetric interpolatory rules for multiple integrals
1086 // over infinite regions with Gaussian weight,
1087 // Journal of Computational and Applied Mathematics,
1088 // Volume 71, 1996, pages 299-309
1089 //
1090 // Florian Heiss, Viktor Winschel,
1091 // Likelihood approximation by numerical integration on sparse grids,
1092 // Journal of Econometrics,
1093 // Volume 144, 2008, pages 62-80.
1094 //
1095 // Thomas Patterson,
1096 // The Optimal Addition of Points to Quadrature Formulae,
1097 // Mathematics of Computation,
1098 // Volume 22, Number 104, October 1968, pages 847-856.
1099 //
1100 // Parameters:
1101 //
1102 // Input, int N, the order.
1103 // N must be 1, 3, 9, 19, or 35.
1104 //
1105 // Output, Scalar X[N], the abscissas.
1106 //
1107 // Output, Scalar W[N], the weights.
1108 //
1109 {
1112 
1113  return;
1114 }
1115 
1116 //****************************************************************************
1117 template<class Scalar>
1119 //****************************************************************************
1120 //
1121 // Purpose:
1122 //
1123 // HERMITE_GENZ_KEISTER_LOOKUP_POINTS looks up Genz-Keister Hermite abscissas.
1124 //
1125 // Discussion:
1126 //
1127 // The integral:
1128 //
1129 // integral ( -oo <= x <= +oo ) f(x) exp ( - x * x ) dx
1130 //
1131 // The quadrature rule:
1132 //
1133 // sum ( 1 <= i <= n ) w(i) * f ( x(i) )
1134 //
1135 // A nested family of rules for the Hermite integration problem
1136 // was produced by Genz and Keister. The structure of the nested
1137 // family was denoted by 1+2+6+10+16, that is, it comprised rules
1138 // of successive orders O = 1, 3, 9, 19, and 35.
1139 //
1140 // The precisions of these rules are P = 1, 5, 15, 29, and 51.
1141 //
1142 // Three related families begin the same way, but end with a different final
1143 // rule. As a convenience, this function includes these final rules as well:
1144 //
1145 // Designation Orders Precisions
1146 //
1147 // 1+2+6+10+16, 1,3,9,19,35 1,5,15,29,51
1148 // 1+2+6+10+18 1,3,9,19,37 1,5,15,29,55
1149 // 1+2+6+10+22 1,3,9,19,41 1,5,15,29,63
1150 // 1+2+6+10+24 1,3,9,19,43 1,5,15,29,67
1151 //
1152 // Some of the data in this function was kindly supplied directly by
1153 // Alan Genz on 24 April 2011.
1154 //
1155 // Licensing:
1156 //
1157 // This code is distributed under the GNU LGPL license.
1158 //
1159 // Modified:
1160 //
1161 // 18 May 2011
1162 //
1163 // Author:
1164 //
1165 // John Burkardt
1166 //
1167 // Reference:
1168 //
1169 // Alan Genz, Bradley Keister,
1170 // Fully symmetric interpolatory rules for multiple integrals
1171 // over infinite regions with Gaussian weight,
1172 // Journal of Computational and Applied Mathematics,
1173 // Volume 71, 1996, pages 299-309
1174 //
1175 // Florian Heiss, Viktor Winschel,
1176 // Likelihood approximation by numerical integration on sparse grids,
1177 // Journal of Econometrics,
1178 // Volume 144, 2008, pages 62-80.
1179 //
1180 // Thomas Patterson,
1181 // The Optimal Addition of Points to Quadrature Formulae,
1182 // Mathematics of Computation,
1183 // Volume 22, Number 104, October 1968, pages 847-856.
1184 //
1185 // Parameters:
1186 //
1187 // Input, int N, the order.
1188 // N must be 1, 3, 9, 19, 35, 27, 41, or 43.
1189 //
1190 // Output, Scalar X[N], the abscissas.
1191 //
1192 {
1193  if (n==1) {
1194  x[ 0] = 0.0000000000000000E+00;
1195  }
1196  else if (n==3) {
1197  x[ 0] = -1.2247448713915889E+00;
1198  x[ 1] = 0.0000000000000000E+00;
1199  x[ 2] = 1.2247448713915889E+00;
1200  }
1201  else if (n==9) {
1202  x[ 0] = -2.9592107790638380E+00;
1203  x[ 1] = -2.0232301911005157E+00;
1204  x[ 2] = -1.2247448713915889E+00;
1205  x[ 3] = -5.2403354748695763E-01;
1206  x[ 4] = 0.0000000000000000E+00;
1207  x[ 5] = 5.2403354748695763E-01;
1208  x[ 6] = 1.2247448713915889E+00;
1209  x[ 7] = 2.0232301911005157E+00;
1210  x[ 8] = 2.9592107790638380E+00;
1211  }
1212  else if (n==19) {
1213  x[ 0] = -4.4995993983103881E+00;
1214  x[ 1] = -3.6677742159463378E+00;
1215  x[ 2] = -2.9592107790638380E+00;
1216  x[ 3] = -2.2665132620567876E+00;
1217  x[ 4] = -2.0232301911005157E+00;
1218  x[ 5] = -1.8357079751751868E+00;
1219  x[ 6] = -1.2247448713915889E+00;
1220  x[ 7] = -8.7004089535290285E-01;
1221  x[ 8] = -5.2403354748695763E-01;
1222  x[ 9] = 0.0000000000000000E+00;
1223  x[10] = 5.2403354748695763E-01;
1224  x[11] = 8.7004089535290285E-01;
1225  x[12] = 1.2247448713915889E+00;
1226  x[13] = 1.8357079751751868E+00;
1227  x[14] = 2.0232301911005157E+00;
1228  x[15] = 2.2665132620567876E+00;
1229  x[16] = 2.9592107790638380E+00;
1230  x[17] = 3.6677742159463378E+00;
1231  x[18] = 4.4995993983103881E+00;
1232  }
1233  else if (n==35) {
1234  x[ 0] = -6.3759392709822356E+00;
1235  x[ 1] = -5.6432578578857449E+00;
1236  x[ 2] = -5.0360899444730940E+00;
1237  x[ 3] = -4.4995993983103881E+00;
1238  x[ 4] = -4.0292201405043713E+00;
1239  x[ 5] = -3.6677742159463378E+00;
1240  x[ 6] = -3.3491639537131945E+00;
1241  x[ 7] = -2.9592107790638380E+00;
1242  x[ 8] = -2.5705583765842968E+00;
1243  x[ 9] = -2.2665132620567876E+00;
1244  x[10] = -2.0232301911005157E+00;
1245  x[11] = -1.8357079751751868E+00;
1246  x[12] = -1.5794121348467671E+00;
1247  x[13] = -1.2247448713915889E+00;
1248  x[14] = -8.7004089535290285E-01;
1249  x[15] = -5.2403354748695763E-01;
1250  x[16] = -1.7606414208200893E-01;
1251  x[17] = 0.0000000000000000E+00;
1252  x[18] = 1.7606414208200893E-01;
1253  x[19] = 5.2403354748695763E-01;
1254  x[20] = 8.7004089535290285E-01;
1255  x[21] = 1.2247448713915889E+00;
1256  x[22] = 1.5794121348467671E+00;
1257  x[23] = 1.8357079751751868E+00;
1258  x[24] = 2.0232301911005157E+00;
1259  x[25] = 2.2665132620567876E+00;
1260  x[26] = 2.5705583765842968E+00;
1261  x[27] = 2.9592107790638380E+00;
1262  x[28] = 3.3491639537131945E+00;
1263  x[29] = 3.6677742159463378E+00;
1264  x[30] = 4.0292201405043713E+00;
1265  x[31] = 4.4995993983103881E+00;
1266  x[32] = 5.0360899444730940E+00;
1267  x[33] = 5.6432578578857449E+00;
1268  x[34] = 6.3759392709822356E+00;
1269  }
1270  else if (n==37) {
1271  x[ 0] = -6.853200069757519;
1272  x[ 1] = -6.124527854622158;
1273  x[ 2] = -5.521865209868350;
1274  x[ 3] = -4.986551454150765;
1275  x[ 4] = -4.499599398310388;
1276  x[ 5] = -4.057956316089741;
1277  x[ 6] = -3.667774215946338;
1278  x[ 7] = -3.315584617593290;
1279  x[ 8] = -2.959210779063838;
1280  x[ 9] = -2.597288631188366;
1281  x[10] = -2.266513262056788;
1282  x[11] = -2.023230191100516;
1283  x[12] = -1.835707975175187;
1284  x[13] = -1.561553427651873;
1285  x[14] = -1.224744871391589;
1286  x[15] = -0.870040895352903;
1287  x[16] = -0.524033547486958;
1288  x[17] = -0.214618180588171;
1289  x[18] = 0.000000000000000;
1290  x[19] = 0.214618180588171;
1291  x[20] = 0.524033547486958;
1292  x[21] = 0.870040895352903;
1293  x[22] = 1.224744871391589;
1294  x[23] = 1.561553427651873;
1295  x[24] = 1.835707975175187;
1296  x[25] = 2.023230191100516;
1297  x[26] = 2.266513262056788;
1298  x[27] = 2.597288631188366;
1299  x[28] = 2.959210779063838;
1300  x[29] = 3.315584617593290;
1301  x[30] = 3.667774215946338;
1302  x[31] = 4.057956316089741;
1303  x[32] = 4.499599398310388;
1304  x[33] = 4.986551454150765;
1305  x[34] = 5.521865209868350;
1306  x[35] = 6.124527854622158;
1307  x[36] = 6.853200069757519;
1308  }
1309  else if (n==41) {
1310  x[ 0] = -7.251792998192644;
1311  x[ 1] = -6.547083258397540;
1312  x[ 2] = -5.961461043404500;
1313  x[ 3] = -5.437443360177798;
1314  x[ 4] = -4.953574342912980;
1315  x[ 5] = -4.4995993983103881;
1316  x[ 6] = -4.070919267883068;
1317  x[ 7] = -3.6677742159463378;
1318  x[ 8] = -3.296114596212218;
1319  x[ 9] = -2.9592107790638380;
1320  x[10] = -2.630415236459871;
1321  x[11] = -2.2665132620567876;
1322  x[12] = -2.043834754429505;
1323  x[13] = -2.0232301911005157;
1324  x[14] = -1.8357079751751868;
1325  x[15] = -1.585873011819188;
1326  x[16] = -1.2247448713915889;
1327  x[17] = -0.87004089535290285;
1328  x[18] = -0.52403354748695763;
1329  x[19] = -0.195324784415805;
1330  x[20] = 0.0000000000000000;
1331  x[21] = 0.195324784415805;
1332  x[22] = 0.52403354748695763;
1333  x[23] = 0.87004089535290285;
1334  x[24] = 1.2247448713915889;
1335  x[25] = 1.585873011819188;
1336  x[26] = 1.8357079751751868;
1337  x[27] = 2.0232301911005157;
1338  x[28] = 2.043834754429505;
1339  x[29] = 2.2665132620567876;
1340  x[30] = 2.630415236459871;
1341  x[31] = 2.9592107790638380;
1342  x[32] = 3.296114596212218;
1343  x[33] = 3.6677742159463378;
1344  x[34] = 4.070919267883068;
1345  x[35] = 4.4995993983103881;
1346  x[36] = 4.953574342912980;
1347  x[37] = 5.437443360177798;
1348  x[38] = 5.961461043404500;
1349  x[39] = 6.547083258397540;
1350  x[40] = 7.251792998192644;
1351  }
1352  else if (n==43) {
1353  x[ 0] = -10.167574994881873;
1354  x[ 1] = -7.231746029072501;
1355  x[ 2] = -6.535398426382995;
1356  x[ 3] = -5.954781975039809;
1357  x[ 4] = -5.434053000365068;
1358  x[ 5] = -4.952329763008589;
1359  x[ 6] = -4.4995993983103881;
1360  x[ 7] = -4.071335874253583;
1361  x[ 8] = -3.6677742159463378;
1362  x[ 9] = -3.295265921534226;
1363  x[10] = -2.9592107790638380;
1364  x[11] = -2.633356763661946;
1365  x[12] = -2.2665132620567876;
1366  x[13] = -2.089340389294661;
1367  x[14] = -2.0232301911005157;
1368  x[15] = -1.8357079751751868;
1369  x[16] = -1.583643465293944;
1370  x[17] = -1.2247448713915889;
1371  x[18] = -0.87004089535290285;
1372  x[19] = -0.52403354748695763;
1373  x[20] = -0.196029453662011;
1374  x[21] = 0.0000000000000000;
1375  x[22] = 0.196029453662011;
1376  x[23] = 0.52403354748695763;
1377  x[24] = 0.87004089535290285;
1378  x[25] = 1.2247448713915889;
1379  x[26] = 1.583643465293944;
1380  x[27] = 1.8357079751751868;
1381  x[28] = 2.0232301911005157;
1382  x[29] = 2.089340389294661;
1383  x[30] = 2.2665132620567876;
1384  x[31] = 2.633356763661946;
1385  x[32] = 2.9592107790638380;
1386  x[33] = 3.295265921534226;
1387  x[34] = 3.6677742159463378;
1388  x[35] = 4.071335874253583;
1389  x[36] = 4.4995993983103881;
1390  x[37] = 4.952329763008589;
1391  x[38] = 5.434053000365068;
1392  x[39] = 5.954781975039809;
1393  x[40] = 6.535398426382995;
1394  x[41] = 7.231746029072501;
1395  x[42] = 10.167574994881873;
1396  }
1397  else {
1398  std::cerr << "\n";
1399  std::cerr << "HERMITE_GENZ_KEISTER_LOOKUP_POINTS - Fatal error!\n";
1400  std::cerr << " Illegal input value of N.\n";
1401  std::cerr << " N must be 1, 3, 9, 19, 35, 37, 41 or 43.\n";
1402  std::exit(1);
1403  }
1404  return;
1405 }
1406 
1407 //****************************************************************************
1408 template<class Scalar>
1410 //****************************************************************************
1411 //
1412 // Purpose:
1413 //
1414 // HERMITE_GENZ_KEISTER_LOOKUP_WEIGHTS looks up Genz-Keister Hermite weights.
1415 //
1416 // Discussion:
1417 //
1418 // The integral:
1419 //
1420 // integral ( -oo <= x <= +oo ) f(x) exp ( - x * x ) dx
1421 //
1422 // The quadrature rule:
1423 //
1424 // sum ( 1 <= i <= n ) w(i) * f ( x(i) )
1425 //
1426 // A nested family of rules for the Hermite integration problem
1427 // was produced by Genz and Keister. The structure of the nested
1428 // family was denoted by 1+2+6+10+16, that is, it comprised rules
1429 // of successive orders O = 1, 3, 9, 19, and 35.
1430 //
1431 // The precisions of these rules are P = 1, 5, 15, 29, and 51.
1432 //
1433 // Three related families begin the same way, but end with a different final
1434 // rule. As a convenience, this function includes these final rules as well:
1435 //
1436 // Designation Orders Precisions
1437 //
1438 // 1+2+6+10+16, 1,3,9,19,35 1,5,15,29,51
1439 // 1+2+6+10+18 1,3,9,19,37 1,5,15,29,55
1440 // 1+2+6+10+22 1,3,9,19,41 1,5,15,29,63
1441 // 1+2+6+10+24 1,3,9,19,43 1,5,15,29,67
1442 //
1443 // Some of the data in this function was kindly supplied directly by
1444 // Alan Genz on 24 April 2011.
1445 //
1446 // Licensing:
1447 //
1448 // This code is distributed under the GNU LGPL license.
1449 //
1450 // Modified:
1451 //
1452 // 18 May 2011
1453 //
1454 // Author:
1455 //
1456 // John Burkardt
1457 //
1458 // Reference:
1459 //
1460 // Alan Genz, Bradley Keister,
1461 // Fully symmetric interpolatory rules for multiple integrals
1462 // over infinite regions with Gaussian weight,
1463 // Journal of Computational and Applied Mathematics,
1464 // Volume 71, 1996, pages 299-309
1465 //
1466 // Florian Heiss, Viktor Winschel,
1467 // Likelihood approximation by numerical integration on sparse grids,
1468 // Journal of Econometrics,
1469 // Volume 144, 2008, pages 62-80.
1470 //
1471 // Thomas Patterson,
1472 // The Optimal Addition of Points to Quadrature Formulae,
1473 // Mathematics of Computation,
1474 // Volume 22, Number 104, October 1968, pages 847-856.
1475 //
1476 // Parameters:
1477 //
1478 // Input, int N, the order.
1479 // N must be 1, 3, 9, 19, 35, 37, 41, or 43.
1480 //
1481 // Output, Scalar W[N], the weights.
1482 //
1483 {
1484  if (n==1) {
1485  w[ 0] = 1.7724538509055159E+00;
1486  }
1487  else if (n==3) {
1488  w[ 0] = 2.9540897515091930E-01;
1489  w[ 1] = 1.1816359006036772E+00;
1490  w[ 2] = 2.9540897515091930E-01;
1491  }
1492  else if (n==9) {
1493  w[ 0] = 1.6708826306882348E-04;
1494  w[ 1] = 1.4173117873979098E-02;
1495  w[ 2] = 1.6811892894767771E-01;
1496  w[ 3] = 4.7869428549114124E-01;
1497  w[ 4] = 4.5014700975378197E-01;
1498  w[ 5] = 4.7869428549114124E-01;
1499  w[ 6] = 1.6811892894767771E-01;
1500  w[ 7] = 1.4173117873979098E-02;
1501  w[ 8] = 1.6708826306882348E-04;
1502  }
1503  else if (n==19) {
1504  w[ 0] = 1.5295717705322357E-09;
1505  w[ 1] = 1.0802767206624762E-06;
1506  w[ 2] = 1.0656589772852267E-04;
1507  w[ 3] = 5.1133174390883855E-03;
1508  w[ 4] = -1.1232438489069229E-02;
1509  w[ 5] = 3.2055243099445879E-02;
1510  w[ 6] = 1.1360729895748269E-01;
1511  w[ 7] = 1.0838861955003017E-01;
1512  w[ 8] = 3.6924643368920851E-01;
1513  w[ 9] = 5.3788160700510168E-01;
1514  w[10] = 3.6924643368920851E-01;
1515  w[11] = 1.0838861955003017E-01;
1516  w[12] = 1.1360729895748269E-01;
1517  w[13] = 3.2055243099445879E-02;
1518  w[14] = -1.1232438489069229E-02;
1519  w[15] = 5.1133174390883855E-03;
1520  w[16] = 1.0656589772852267E-04;
1521  w[17] = 1.0802767206624762E-06;
1522  w[18] = 1.5295717705322357E-09;
1523  }
1524  else if (n==35) {
1525  w[ 0] = 1.8684014894510604E-18;
1526  w[ 1] = 9.6599466278563243E-15;
1527  w[ 2] = 5.4896836948499462E-12;
1528  w[ 3] = 8.1553721816916897E-10;
1529  w[ 4] = 3.7920222392319532E-08;
1530  w[ 5] = 4.3737818040926989E-07;
1531  w[ 6] = 4.8462799737020461E-06;
1532  w[ 7] = 6.3328620805617891E-05;
1533  w[ 8] = 4.8785399304443770E-04;
1534  w[ 9] = 1.4515580425155904E-03;
1535  w[10] = 4.0967527720344047E-03;
1536  w[11] = 5.5928828911469180E-03;
1537  w[12] = 2.7780508908535097E-02;
1538  w[13] = 8.0245518147390893E-02;
1539  w[14] = 1.6371221555735804E-01;
1540  w[15] = 2.6244871488784277E-01;
1541  w[16] = 3.3988595585585218E-01;
1542  w[17] = 9.1262675363737921E-04;
1543  w[18] = 3.3988595585585218E-01;
1544  w[19] = 2.6244871488784277E-01;
1545  w[20] = 1.6371221555735804E-01;
1546  w[21] = 8.0245518147390893E-02;
1547  w[22] = 2.7780508908535097E-02;
1548  w[23] = 5.5928828911469180E-03;
1549  w[24] = 4.0967527720344047E-03;
1550  w[25] = 1.4515580425155904E-03;
1551  w[26] = 4.8785399304443770E-04;
1552  w[27] = 6.3328620805617891E-05;
1553  w[28] = 4.8462799737020461E-06;
1554  w[29] = 4.3737818040926989E-07;
1555  w[30] = 3.7920222392319532E-08;
1556  w[31] = 8.1553721816916897E-10;
1557  w[32] = 5.4896836948499462E-12;
1558  w[33] = 9.6599466278563243E-15;
1559  w[34] = 1.8684014894510604E-18;
1560  }
1561  else if (n==37) {
1562  w[ 0] = 0.19030350940130498E-20;
1563  w[ 1] = 0.187781893143728947E-16;
1564  w[ 2] = 0.182242751549129356E-13;
1565  w[ 3] = 0.45661763676186859E-11;
1566  w[ 4] = 0.422525843963111041E-09;
1567  w[ 5] = 0.16595448809389819E-07;
1568  w[ 6] = 0.295907520230744049E-06;
1569  w[ 7] = 0.330975870979203419E-05;
1570  w[ 8] = 0.32265185983739747E-04;
1571  w[ 9] = 0.234940366465975222E-03;
1572  w[10] = 0.985827582996483824E-03;
1573  w[11] = 0.176802225818295443E-02;
1574  w[12] = 0.43334988122723492E-02;
1575  w[13] = 0.15513109874859354E-01;
1576  w[14] = 0.442116442189845444E-01;
1577  w[15] = 0.937208280655245902E-01;
1578  w[16] = 0.143099302896833389E+00;
1579  w[17] = 0.147655710402686249E+00;
1580  w[18] = 0.968824552928425499E-01;
1581  w[19] = 0.147655710402686249E+00;
1582  w[20] = 0.143099302896833389E+00;
1583  w[21] = 0.937208280655245902E-01;
1584  w[22] = 0.442116442189845444E-01;
1585  w[23] = 0.15513109874859354E-01;
1586  w[24] = 0.43334988122723492E-02;
1587  w[25] = 0.176802225818295443E-02;
1588  w[26] = 0.985827582996483824E-03;
1589  w[27] = 0.234940366465975222E-03;
1590  w[28] = 0.32265185983739747E-04;
1591  w[29] = 0.330975870979203419E-05;
1592  w[30] = 0.295907520230744049E-06;
1593  w[31] = 0.16595448809389819E-07;
1594  w[32] = 0.422525843963111041E-09;
1595  w[33] = 0.45661763676186859E-11;
1596  w[34] = 0.182242751549129356E-13;
1597  w[35] = 0.187781893143728947E-16;
1598  w[36] = 0.19030350940130498E-20;
1599  }
1600  else if (n==41) {
1601  w[ 0] = 0.664195893812757801E-23;
1602  w[ 1] = 0.860427172512207236E-19;
1603  w[ 2] = 0.1140700785308509E-15;
1604  w[ 3] = 0.408820161202505983E-13;
1605  w[ 4] = 0.581803393170320419E-11;
1606  w[ 5] = 0.400784141604834759E-09;
1607  w[ 6] = 0.149158210417831408E-07;
1608  w[ 7] = 0.315372265852264871E-06;
1609  w[ 8] = 0.381182791749177506E-05;
1610  w[ 9] = 0.288976780274478689E-04;
1611  w[10] = 0.189010909805097887E-03;
1612  w[11] = 0.140697424065246825E-02;
1613  w[12] = - 0.144528422206988237E-01;
1614  w[13] = 0.178852543033699732E-01;
1615  w[14] = 0.705471110122962612E-03;
1616  w[15] = 0.165445526705860772E-01;
1617  w[16] = 0.45109010335859128E-01;
1618  w[17] = 0.928338228510111845E-01;
1619  w[18] = 0.145966293895926429E+00;
1620  w[19] = 0.165639740400529554E+00;
1621  w[20] = 0.562793426043218877E-01;
1622  w[21] = 0.165639740400529554E+00;
1623  w[22] = 0.145966293895926429E+00;
1624  w[23] = 0.928338228510111845E-01;
1625  w[24] = 0.45109010335859128E-01;
1626  w[25] = 0.165445526705860772E-01;
1627  w[26] = 0.705471110122962612E-03;
1628  w[27] = 0.178852543033699732E-01;
1629  w[28] = - 0.144528422206988237E-01;
1630  w[29] = 0.140697424065246825E-02;
1631  w[30] = 0.189010909805097887E-03;
1632  w[31] = 0.288976780274478689E-04;
1633  w[32] = 0.381182791749177506E-05;
1634  w[33] = 0.315372265852264871E-06;
1635  w[34] = 0.149158210417831408E-07;
1636  w[35] = 0.400784141604834759E-09;
1637  w[36] = 0.581803393170320419E-11;
1638  w[37] = 0.408820161202505983E-13;
1639  w[38] = 0.1140700785308509E-15;
1640  w[39] = 0.860427172512207236E-19;
1641  w[40] = 0.664195893812757801E-23;
1642  }
1643  else if (n==43) {
1644  w[ 0] = 0.546191947478318097E-37;
1645  w[ 1] = 0.87544909871323873E-23;
1646  w[ 2] = 0.992619971560149097E-19;
1647  w[ 3] = 0.122619614947864357E-15;
1648  w[ 4] = 0.421921851448196032E-13;
1649  w[ 5] = 0.586915885251734856E-11;
1650  w[ 6] = 0.400030575425776948E-09;
1651  w[ 7] = 0.148653643571796457E-07;
1652  w[ 8] = 0.316018363221289247E-06;
1653  w[ 9] = 0.383880761947398577E-05;
1654  w[10] = 0.286802318064777813E-04;
1655  w[11] = 0.184789465688357423E-03;
1656  w[12] = 0.150909333211638847E-02;
1657  w[13] = - 0.38799558623877157E-02;
1658  w[14] = 0.67354758901013295E-02;
1659  w[15] = 0.139966252291568061E-02;
1660  w[16] = 0.163616873493832402E-01;
1661  w[17] = 0.450612329041864976E-01;
1662  w[18] = 0.928711584442575456E-01;
1663  w[19] = 0.145863292632147353E+00;
1664  w[20] = 0.164880913687436689E+00;
1665  w[21] = 0.579595986101181095E-01;
1666  w[22] = 0.164880913687436689E+00;
1667  w[23] = 0.145863292632147353E+00;
1668  w[24] = 0.928711584442575456E-01;
1669  w[25] = 0.450612329041864976E-01;
1670  w[26] = 0.163616873493832402E-01;
1671  w[27] = 0.139966252291568061E-02;
1672  w[28] = 0.67354758901013295E-02;
1673  w[29] = - 0.38799558623877157E-02;
1674  w[30] = 0.150909333211638847E-02;
1675  w[31] = 0.184789465688357423E-03;
1676  w[32] = 0.286802318064777813E-04;
1677  w[33] = 0.383880761947398577E-05;
1678  w[34] = 0.316018363221289247E-06;
1679  w[35] = 0.148653643571796457E-07;
1680  w[36] = 0.400030575425776948E-09;
1681  w[37] = 0.586915885251734856E-11;
1682  w[38] = 0.421921851448196032E-13;
1683  w[39] = 0.122619614947864357E-15;
1684  w[40] = 0.992619971560149097E-19;
1685  w[41] = 0.87544909871323873E-23;
1686  w[42] = 0.546191947478318097E-37;
1687  }
1688  else {
1689  std::cerr << "\n";
1690  std::cerr << "HERMITE_GENZ_KEISTER_LOOKUP_WEIGHTS - Fatal error!\n";
1691  std::cerr << " Illegal input value of N.\n";
1692  std::cerr << " N must be 1, 3, 9, 19, 35, 37, 41 or 43.\n";
1693  std::exit(1);
1694  }
1695  return;
1696 }
1697 
1698 //****************************************************************************
1699 template<class Scalar>
1700 void IntrepidBurkardtRules::hermite_lookup ( int n, Scalar x[], Scalar w[] )
1701 //****************************************************************************
1702 //
1703 // Purpose:
1704 //
1705 // HERMITE_LOOKUP looks up abscissas and weights for Gauss-Hermite quadrature.
1706 //
1707 // Licensing:
1708 //
1709 // This code is distributed under the GNU LGPL license.
1710 //
1711 // Modified:
1712 //
1713 // 27 April 2010
1714 //
1715 // Author:
1716 //
1717 // John Burkardt
1718 //
1719 // Reference:
1720 //
1721 // Milton Abramowitz, Irene Stegun,
1722 // Handbook of Mathematical Functions,
1723 // National Bureau of Standards, 1964,
1724 // ISBN: 0-486-61272-4,
1725 // LC: QA47.A34.
1726 //
1727 // Vladimir Krylov,
1728 // Approximate Calculation of Integrals,
1729 // Dover, 2006,
1730 // ISBN: 0486445798.
1731 // LC: QA311.K713.
1732 //
1733 // Arthur Stroud, Don Secrest,
1734 // Gaussian Quadrature Formulas,
1735 // Prentice Hall, 1966,
1736 // LC: QA299.4G3S7.
1737 //
1738 // Stephen Wolfram,
1739 // The Mathematica Book,
1740 // Fourth Edition,
1741 // Cambridge University Press, 1999,
1742 // ISBN: 0-521-64314-7,
1743 // LC: QA76.95.W65.
1744 //
1745 // Daniel Zwillinger, editor,
1746 // CRC Standard Mathematical Tables and Formulae,
1747 // 30th Edition,
1748 // CRC Press, 1996,
1749 // ISBN: 0-8493-2479-3,
1750 // LC: QA47.M315.
1751 //
1752 // Parameters:
1753 //
1754 // Input, int N, the order.
1755 // N must be between 1 and 20.
1756 //
1757 // Output, Scalar X[N], the abscissas.
1758 //
1759 // Output, Scalar W[N], the weights.
1760 //
1761 {
1764 
1765  return;
1766 }
1767 
1768 //****************************************************************************
1769 template<class Scalar>
1771 //****************************************************************************
1772 //
1773 // Purpose:
1774 //
1775 // HERMITE_LOOKUP_POINTS looks up abscissas for Hermite quadrature.
1776 //
1777 // Discussion:
1778 //
1779 // The integral:
1780 //
1781 // integral ( -oo < x < +oo ) exp ( - x * x ) * f(x) dx
1782 //
1783 // The quadrature rule:
1784 //
1785 // sum ( 1 <= i <= n ) w(i) * f ( x(i) ).
1786 //
1787 // Mathematica can numerically estimate the abscissas
1788 // of order N to P digits by the command:
1789 //
1790 // NSolve [ HermiteH [ n, x ] == 0, x, p ]
1791 //
1792 // Licensing:
1793 //
1794 // This code is distributed under the GNU LGPL license.
1795 //
1796 // Modified:
1797 //
1798 // 27 April 2010
1799 //
1800 // Author:
1801 //
1802 // John Burkardt
1803 //
1804 // Reference:
1805 //
1806 // Milton Abramowitz, Irene Stegun,
1807 // Handbook of Mathematical Functions,
1808 // National Bureau of Standards, 1964,
1809 // ISBN: 0-486-61272-4,
1810 // LC: QA47.A34.
1811 //
1812 // Vladimir Krylov,
1813 // Approximate Calculation of Integrals,
1814 // Dover, 2006,
1815 // ISBN: 0486445798,
1816 // LC: QA311.K713.
1817 //
1818 // Arthur Stroud, Don Secrest,
1819 // Gaussian Quadrature Formulas,
1820 // Prentice Hall, 1966,
1821 // LC: QA299.4G3S7.
1822 //
1823 // Stephen Wolfram,
1824 // The Mathematica Book,
1825 // Fourth Edition,
1826 // Cambridge University Press, 1999,
1827 // ISBN: 0-521-64314-7,
1828 // LC: QA76.95.W65.
1829 //
1830 // Daniel Zwillinger, editor,
1831 // CRC Standard Mathematical Tables and Formulae,
1832 // 30th Edition,
1833 // CRC Press, 1996,
1834 // ISBN: 0-8493-2479-3,
1835 // LC: QA47.M315.
1836 //
1837 // Parameters:
1838 //
1839 // Input, int N, the order.
1840 // N must be between 1 and 20.
1841 //
1842 // Output, Scalar X[N], the abscissas.
1843 //
1844 {
1845  if (n==1) {
1846  x[ 0] = 0.0;
1847  }
1848  else if(n==2) {
1849  x[ 0] = - 0.707106781186547524400844362105E+00;
1850  x[ 1] = 0.707106781186547524400844362105E+00;
1851  }
1852  else if (n==3) {
1853  x[ 0] = - 0.122474487139158904909864203735E+01;
1854  x[ 1] = 0.0E+00;
1855  x[ 2] = 0.122474487139158904909864203735E+01;
1856  }
1857  else if (n==4) {
1858  x[ 0] = - 0.165068012388578455588334111112E+01;
1859  x[ 1] = - 0.524647623275290317884060253835E+00;
1860  x[ 2] = 0.524647623275290317884060253835E+00;
1861  x[ 3] = 0.165068012388578455588334111112E+01;
1862  }
1863  else if (n==5) {
1864  x[ 0] = - 0.202018287045608563292872408814E+01;
1865  x[ 1] = - 0.958572464613818507112770593893E+00;
1866  x[ 2] = 0.0E+00;
1867  x[ 3] = 0.958572464613818507112770593893E+00;
1868  x[ 4] = 0.202018287045608563292872408814E+01;
1869  }
1870  else if (n==6) {
1871  x[ 0] = - 0.235060497367449222283392198706E+01;
1872  x[ 1] = - 0.133584907401369694971489528297E+01;
1873  x[ 2] = - 0.436077411927616508679215948251E+00;
1874  x[ 3] = 0.436077411927616508679215948251E+00;
1875  x[ 4] = 0.133584907401369694971489528297E+01;
1876  x[ 5] = 0.235060497367449222283392198706E+01;
1877  }
1878  else if (n==7) {
1879  x[ 0] = - 0.265196135683523349244708200652E+01;
1880  x[ 1] = - 0.167355162876747144503180139830E+01;
1881  x[ 2] = - 0.816287882858964663038710959027E+00;
1882  x[ 3] = 0.0E+00;
1883  x[ 4] = 0.816287882858964663038710959027E+00;
1884  x[ 5] = 0.167355162876747144503180139830E+01;
1885  x[ 6] = 0.265196135683523349244708200652E+01;
1886  }
1887  else if (n==8) {
1888  x[ 0] = - 0.293063742025724401922350270524E+01;
1889  x[ 1] = - 0.198165675669584292585463063977E+01;
1890  x[ 2] = - 0.115719371244678019472076577906E+01;
1891  x[ 3] = - 0.381186990207322116854718885584E+00;
1892  x[ 4] = 0.381186990207322116854718885584E+00;
1893  x[ 5] = 0.115719371244678019472076577906E+01;
1894  x[ 6] = 0.198165675669584292585463063977E+01;
1895  x[ 7] = 0.293063742025724401922350270524E+01;
1896  }
1897  else if (n==9) {
1898  x[ 0] = - 0.319099320178152760723004779538E+01;
1899  x[ 1] = - 0.226658058453184311180209693284E+01;
1900  x[ 2] = - 0.146855328921666793166701573925E+01;
1901  x[ 3] = - 0.723551018752837573322639864579E+00;
1902  x[ 4] = 0.0E+00;
1903  x[ 5] = 0.723551018752837573322639864579E+00;
1904  x[ 6] = 0.146855328921666793166701573925E+01;
1905  x[ 7] = 0.226658058453184311180209693284E+01;
1906  x[ 8] = 0.319099320178152760723004779538E+01;
1907  }
1908  else if (n==10) {
1909  x[ 0] = - 0.343615911883773760332672549432E+01;
1910  x[ 1] = - 0.253273167423278979640896079775E+01;
1911  x[ 2] = - 0.175668364929988177345140122011E+01;
1912  x[ 3] = - 0.103661082978951365417749191676E+01;
1913  x[ 4] = - 0.342901327223704608789165025557E+00;
1914  x[ 5] = 0.342901327223704608789165025557E+00;
1915  x[ 6] = 0.103661082978951365417749191676E+01;
1916  x[ 7] = 0.175668364929988177345140122011E+01;
1917  x[ 8] = 0.253273167423278979640896079775E+01;
1918  x[ 9] = 0.343615911883773760332672549432E+01;
1919  }
1920  else if (n==11) {
1921  x[ 0] = - 0.366847084655958251845837146485E+01;
1922  x[ 1] = - 0.278329009978165177083671870152E+01;
1923  x[ 2] = - 0.202594801582575533516591283121E+01;
1924  x[ 3] = - 0.132655708449493285594973473558E+01;
1925  x[ 4] = - 0.656809566882099765024611575383E+00;
1926  x[ 5] = 0.0E+00;
1927  x[ 6] = 0.656809566882099765024611575383E+00;
1928  x[ 7] = 0.132655708449493285594973473558E+01;
1929  x[ 8] = 0.202594801582575533516591283121E+01;
1930  x[ 9] = 0.278329009978165177083671870152E+01;
1931  x[10] = 0.366847084655958251845837146485E+01;
1932  }
1933  else if (n==12) {
1934  x[ 0] = - 0.388972489786978191927164274724E+01;
1935  x[ 1] = - 0.302063702512088977171067937518E+01;
1936  x[ 2] = - 0.227950708050105990018772856942E+01;
1937  x[ 3] = - 0.159768263515260479670966277090E+01;
1938  x[ 4] = - 0.947788391240163743704578131060E+00;
1939  x[ 5] = - 0.314240376254359111276611634095E+00;
1940  x[ 6] = 0.314240376254359111276611634095E+00;
1941  x[ 7] = 0.947788391240163743704578131060E+00;
1942  x[ 8] = 0.159768263515260479670966277090E+01;
1943  x[ 9] = 0.227950708050105990018772856942E+01;
1944  x[10] = 0.302063702512088977171067937518E+01;
1945  x[11] = 0.388972489786978191927164274724E+01;
1946  }
1947  else if (n==13) {
1948  x[ 0] = - 0.410133759617863964117891508007E+01;
1949  x[ 1] = - 0.324660897837240998812205115236E+01;
1950  x[ 2] = - 0.251973568567823788343040913628E+01;
1951  x[ 3] = - 0.185310765160151214200350644316E+01;
1952  x[ 4] = - 0.122005503659074842622205526637E+01;
1953  x[ 5] = - 0.605763879171060113080537108602E+00;
1954  x[ 6] = 0.0E+00;
1955  x[ 7] = 0.605763879171060113080537108602E+00;
1956  x[ 8] = 0.122005503659074842622205526637E+01;
1957  x[ 9] = 0.185310765160151214200350644316E+01;
1958  x[10] = 0.251973568567823788343040913628E+01;
1959  x[11] = 0.324660897837240998812205115236E+01;
1960  x[12] = 0.410133759617863964117891508007E+01;
1961  }
1962  else if (n==14) {
1963  x[ 0] = - 0.430444857047363181262129810037E+01;
1964  x[ 1] = - 0.346265693360227055020891736115E+01;
1965  x[ 2] = - 0.274847072498540256862499852415E+01;
1966  x[ 3] = - 0.209518325850771681573497272630E+01;
1967  x[ 4] = - 0.147668273114114087058350654421E+01;
1968  x[ 5] = - 0.878713787329399416114679311861E+00;
1969  x[ 6] = - 0.291745510672562078446113075799E+00;
1970  x[ 7] = 0.291745510672562078446113075799E+00;
1971  x[ 8] = 0.878713787329399416114679311861E+00;
1972  x[ 9] = 0.147668273114114087058350654421E+01;
1973  x[10] = 0.209518325850771681573497272630E+01;
1974  x[11] = 0.274847072498540256862499852415E+01;
1975  x[12] = 0.346265693360227055020891736115E+01;
1976  x[13] = 0.430444857047363181262129810037E+01;
1977  }
1978  else if (n==15) {
1979  x[ 0] = - 0.449999070730939155366438053053E+01;
1980  x[ 1] = - 0.366995037340445253472922383312E+01;
1981  x[ 2] = - 0.296716692790560324848896036355E+01;
1982  x[ 3] = - 0.232573248617385774545404479449E+01;
1983  x[ 4] = - 0.171999257518648893241583152515E+01;
1984  x[ 5] = - 0.113611558521092066631913490556E+01;
1985  x[ 6] = - 0.565069583255575748526020337198E+00;
1986  x[ 7] = 0.0E+00;
1987  x[ 8] = 0.565069583255575748526020337198E+00;
1988  x[ 9] = 0.113611558521092066631913490556E+01;
1989  x[10] = 0.171999257518648893241583152515E+01;
1990  x[11] = 0.232573248617385774545404479449E+01;
1991  x[12] = 0.296716692790560324848896036355E+01;
1992  x[13] = 0.366995037340445253472922383312E+01;
1993  x[14] = 0.449999070730939155366438053053E+01;
1994  }
1995  else if (n==16) {
1996  x[ 0] = - 0.468873893930581836468849864875E+01;
1997  x[ 1] = - 0.386944790486012269871942409801E+01;
1998  x[ 2] = - 0.317699916197995602681399455926E+01;
1999  x[ 3] = - 0.254620215784748136215932870545E+01;
2000  x[ 4] = - 0.195178799091625397743465541496E+01;
2001  x[ 5] = - 0.138025853919888079637208966969E+01;
2002  x[ 6] = - 0.822951449144655892582454496734E+00;
2003  x[ 7] = - 0.273481046138152452158280401965E+00;
2004  x[ 8] = 0.273481046138152452158280401965E+00;
2005  x[ 9] = 0.822951449144655892582454496734E+00;
2006  x[10] = 0.138025853919888079637208966969E+01;
2007  x[11] = 0.195178799091625397743465541496E+01;
2008  x[12] = 0.254620215784748136215932870545E+01;
2009  x[13] = 0.317699916197995602681399455926E+01;
2010  x[14] = 0.386944790486012269871942409801E+01;
2011  x[15] = 0.468873893930581836468849864875E+01;
2012  }
2013  else if (n==17) {
2014  x[ 0] = - 0.487134519367440308834927655662E+01;
2015  x[ 1] = - 0.406194667587547430689245559698E+01;
2016  x[ 2] = - 0.337893209114149408338327069289E+01;
2017  x[ 3] = - 0.275776291570388873092640349574E+01;
2018  x[ 4] = - 0.217350282666662081927537907149E+01;
2019  x[ 5] = - 0.161292431422123133311288254454E+01;
2020  x[ 6] = - 0.106764872574345055363045773799E+01;
2021  x[ 7] = - 0.531633001342654731349086553718E+00;
2022  x[ 8] = 0.0E+00;
2023  x[ 9] = 0.531633001342654731349086553718E+00;
2024  x[10] = 0.106764872574345055363045773799E+01;
2025  x[11] = 0.161292431422123133311288254454E+01;
2026  x[12] = 0.217350282666662081927537907149E+01;
2027  x[13] = 0.275776291570388873092640349574E+01;
2028  x[14] = 0.337893209114149408338327069289E+01;
2029  x[15] = 0.406194667587547430689245559698E+01;
2030  x[16] = 0.487134519367440308834927655662E+01;
2031  }
2032  else if (n==18) {
2033  x[ 0] = - 0.504836400887446676837203757885E+01;
2034  x[ 1] = - 0.424811787356812646302342016090E+01;
2035  x[ 2] = - 0.357376906848626607950067599377E+01;
2036  x[ 3] = - 0.296137750553160684477863254906E+01;
2037  x[ 4] = - 0.238629908916668600026459301424E+01;
2038  x[ 5] = - 0.183553160426162889225383944409E+01;
2039  x[ 6] = - 0.130092085838961736566626555439E+01;
2040  x[ 7] = - 0.776682919267411661316659462284E+00;
2041  x[ 8] = - 0.258267750519096759258116098711E+00;
2042  x[ 9] = 0.258267750519096759258116098711E+00;
2043  x[10] = 0.776682919267411661316659462284E+00;
2044  x[11] = 0.130092085838961736566626555439E+01;
2045  x[12] = 0.183553160426162889225383944409E+01;
2046  x[13] = 0.238629908916668600026459301424E+01;
2047  x[14] = 0.296137750553160684477863254906E+01;
2048  x[15] = 0.357376906848626607950067599377E+01;
2049  x[16] = 0.424811787356812646302342016090E+01;
2050  x[17] = 0.504836400887446676837203757885E+01;
2051  }
2052  else if (n==19) {
2053  x[ 0] = - 0.522027169053748216460967142500E+01;
2054  x[ 1] = - 0.442853280660377943723498532226E+01;
2055  x[ 2] = - 0.376218735196402009751489394104E+01;
2056  x[ 3] = - 0.315784881834760228184318034120E+01;
2057  x[ 4] = - 0.259113378979454256492128084112E+01;
2058  x[ 5] = - 0.204923170985061937575050838669E+01;
2059  x[ 6] = - 0.152417061939353303183354859367E+01;
2060  x[ 7] = - 0.101036838713431135136859873726E+01;
2061  x[ 8] = - 0.503520163423888209373811765050E+00;
2062  x[ 9] = 0.0E+00;
2063  x[10] = 0.503520163423888209373811765050E+00;
2064  x[11] = 0.101036838713431135136859873726E+01;
2065  x[12] = 0.152417061939353303183354859367E+01;
2066  x[13] = 0.204923170985061937575050838669E+01;
2067  x[14] = 0.259113378979454256492128084112E+01;
2068  x[15] = 0.315784881834760228184318034120E+01;
2069  x[16] = 0.376218735196402009751489394104E+01;
2070  x[17] = 0.442853280660377943723498532226E+01;
2071  x[18] = 0.522027169053748216460967142500E+01;
2072  }
2073  else if (n==20) {
2074  x[ 0] = - 0.538748089001123286201690041068E+01;
2075  x[ 1] = - 0.460368244955074427307767524898E+01;
2076  x[ 2] = - 0.394476404011562521037562880052E+01;
2077  x[ 3] = - 0.334785456738321632691492452300E+01;
2078  x[ 4] = - 0.278880605842813048052503375640E+01;
2079  x[ 5] = - 0.225497400208927552308233334473E+01;
2080  x[ 6] = - 0.173853771211658620678086566214E+01;
2081  x[ 7] = - 0.123407621539532300788581834696E+01;
2082  x[ 8] = - 0.737473728545394358705605144252E+00;
2083  x[ 9] = - 0.245340708300901249903836530634E+00;
2084  x[10] = 0.245340708300901249903836530634E+00;
2085  x[11] = 0.737473728545394358705605144252E+00;
2086  x[12] = 0.123407621539532300788581834696E+01;
2087  x[13] = 0.173853771211658620678086566214E+01;
2088  x[14] = 0.225497400208927552308233334473E+01;
2089  x[15] = 0.278880605842813048052503375640E+01;
2090  x[16] = 0.334785456738321632691492452300E+01;
2091  x[17] = 0.394476404011562521037562880052E+01;
2092  x[18] = 0.460368244955074427307767524898E+01;
2093  x[19] = 0.538748089001123286201690041068E+01;
2094  }
2095  else {
2096  std::cerr << "\n";
2097  std::cerr << "HERMITE_LOOKUP_POINTS - Fatal error!\n";
2098  std::cerr << " Illegal value of N = " << n << "\n";
2099  std::cerr << " Legal values are 1 through 20.\n";
2100  std::exit(1);
2101  }
2102 
2103  return;
2104 }
2105 
2106 //****************************************************************************
2107 template<class Scalar>
2109 //****************************************************************************
2110 //
2111 // Purpose:
2112 //
2113 // HERMITE_LOOKUP_WEIGHTS looks up weights for Hermite quadrature.
2114 //
2115 // Discussion:
2116 //
2117 // The integral:
2118 //
2119 // integral ( -oo < x < +oo ) exp ( - x * x ) * f(x) dx
2120 //
2121 // The quadrature rule:
2122 //
2123 // sum ( 1 <= i <= n ) w(i) * f ( x(i) ).
2124 //
2125 // Mathematica can numerically estimate the abscissas
2126 // of order N to P digits by the command:
2127 //
2128 // NSolve [ HermiteH [ n, x ] == 0, x, p ]
2129 //
2130 // Licensing:
2131 //
2132 // This code is distributed under the GNU LGPL license.
2133 //
2134 // Modified:
2135 //
2136 // 27 April 2010
2137 //
2138 // Author:
2139 //
2140 // John Burkardt
2141 //
2142 // Reference:
2143 //
2144 // Milton Abramowitz, Irene Stegun,
2145 // Handbook of Mathematical Functions,
2146 // National Bureau of Standards, 1964,
2147 // ISBN: 0-486-61272-4,
2148 // LC: QA47.A34.
2149 //
2150 // Vladimir Krylov,
2151 // Approximate Calculation of Integrals,
2152 // Dover, 2006,
2153 // ISBN: 0486445798,
2154 // LC: QA311.K713.
2155 //
2156 // Arthur Stroud, Don Secrest,
2157 // Gaussian Quadrature Formulas,
2158 // Prentice Hall, 1966,
2159 // LC: QA299.4G3S7.
2160 //
2161 // Stephen Wolfram,
2162 // The Mathematica Book,
2163 // Fourth Edition,
2164 // Cambridge University Press, 1999,
2165 // ISBN: 0-521-64314-7,
2166 // LC: QA76.95.W65.
2167 //
2168 // Daniel Zwillinger, editor,
2169 // CRC Standard Mathematical Tables and Formulae,
2170 // 30th Edition,
2171 // CRC Press, 1996,
2172 // ISBN: 0-8493-2479-3,
2173 // LC: QA47.M315.
2174 //
2175 // Parameters:
2176 //
2177 // Input, int N, the order.
2178 // N must be between 1 and 20.
2179 //
2180 // Output, Scalar W[N], the weights.
2181 //
2182 {
2183  if (n==1) {
2184  w[ 0] = 1.77245385090551602729816748334;
2185  }
2186  else if (n==2) {
2187  w[ 0] = 0.886226925452758013649083741671E+00;
2188  w[ 1] = 0.886226925452758013649083741671E+00;
2189  }
2190  else if (n==3) {
2191  w[ 0] = 0.295408975150919337883027913890E+00;
2192  w[ 1] = 0.118163590060367735153211165556E+01;
2193  w[ 2] = 0.295408975150919337883027913890E+00;
2194  }
2195  else if (n==4) {
2196  w[ 0] = 0.813128354472451771430345571899E-01;
2197  w[ 1] = 0.804914090005512836506049184481E+00;
2198  w[ 2] = 0.804914090005512836506049184481E+00;
2199  w[ 3] = 0.813128354472451771430345571899E-01;
2200  }
2201  else if (n==5) {
2202  w[ 0] = 0.199532420590459132077434585942E-01;
2203  w[ 1] = 0.393619323152241159828495620852E+00;
2204  w[ 2] = 0.945308720482941881225689324449E+00;
2205  w[ 3] = 0.393619323152241159828495620852E+00;
2206  w[ 4] = 0.199532420590459132077434585942E-01;
2207  }
2208  else if (n==6) {
2209  w[ 0] = 0.453000990550884564085747256463E-02;
2210  w[ 1] = 0.157067320322856643916311563508E+00;
2211  w[ 2] = 0.724629595224392524091914705598E+00;
2212  w[ 3] = 0.724629595224392524091914705598E+00;
2213  w[ 4] = 0.157067320322856643916311563508E+00;
2214  w[ 5] = 0.453000990550884564085747256463E-02;
2215  }
2216  else if (n==7) {
2217  w[ 0] = 0.971781245099519154149424255939E-03;
2218  w[ 1] = 0.545155828191270305921785688417E-01;
2219  w[ 2] = 0.425607252610127800520317466666E+00;
2220  w[ 3] = 0.810264617556807326764876563813E+00;
2221  w[ 4] = 0.425607252610127800520317466666E+00;
2222  w[ 5] = 0.545155828191270305921785688417E-01;
2223  w[ 6] = 0.971781245099519154149424255939E-03;
2224  }
2225  else if (n==8) {
2226  w[ 0] = 0.199604072211367619206090452544E-03;
2227  w[ 1] = 0.170779830074134754562030564364E-01;
2228  w[ 2] = 0.207802325814891879543258620286E+00;
2229  w[ 3] = 0.661147012558241291030415974496E+00;
2230  w[ 4] = 0.661147012558241291030415974496E+00;
2231  w[ 5] = 0.207802325814891879543258620286E+00;
2232  w[ 6] = 0.170779830074134754562030564364E-01;
2233  w[ 7] = 0.199604072211367619206090452544E-03;
2234  }
2235  else if (n==9) {
2236  w[ 0] = 0.396069772632643819045862946425E-04;
2237  w[ 1] = 0.494362427553694721722456597763E-02;
2238  w[ 2] = 0.884745273943765732879751147476E-01;
2239  w[ 3] = 0.432651559002555750199812112956E+00;
2240  w[ 4] = 0.720235215606050957124334723389E+00;
2241  w[ 5] = 0.432651559002555750199812112956E+00;
2242  w[ 6] = 0.884745273943765732879751147476E-01;
2243  w[ 7] = 0.494362427553694721722456597763E-02;
2244  w[ 8] = 0.396069772632643819045862946425E-04;
2245  }
2246  else if (n==10) {
2247  w[ 0] = 0.764043285523262062915936785960E-05;
2248  w[ 1] = 0.134364574678123269220156558585E-02;
2249  w[ 2] = 0.338743944554810631361647312776E-01;
2250  w[ 3] = 0.240138611082314686416523295006E+00;
2251  w[ 4] = 0.610862633735325798783564990433E+00;
2252  w[ 5] = 0.610862633735325798783564990433E+00;
2253  w[ 6] = 0.240138611082314686416523295006E+00;
2254  w[ 7] = 0.338743944554810631361647312776E-01;
2255  w[ 8] = 0.134364574678123269220156558585E-02;
2256  w[ 9] = 0.764043285523262062915936785960E-05;
2257  }
2258  else if (n==11) {
2259  w[ 0] = 0.143956039371425822033088366032E-05;
2260  w[ 1] = 0.346819466323345510643413772940E-03;
2261  w[ 2] = 0.119113954449115324503874202916E-01;
2262  w[ 3] = 0.117227875167708503381788649308E+00;
2263  w[ 4] = 0.429359752356125028446073598601E+00;
2264  w[ 5] = 0.654759286914591779203940657627E+00;
2265  w[ 6] = 0.429359752356125028446073598601E+00;
2266  w[ 7] = 0.117227875167708503381788649308E+00;
2267  w[ 8] = 0.119113954449115324503874202916E-01;
2268  w[ 9] = 0.346819466323345510643413772940E-03;
2269  w[10] = 0.143956039371425822033088366032E-05;
2270  }
2271  else if (n==12) {
2272  w[ 0] = 0.265855168435630160602311400877E-06;
2273  w[ 1] = 0.857368704358785865456906323153E-04;
2274  w[ 2] = 0.390539058462906185999438432620E-02;
2275  w[ 3] = 0.516079856158839299918734423606E-01;
2276  w[ 4] = 0.260492310264161129233396139765E+00;
2277  w[ 5] = 0.570135236262479578347113482275E+00;
2278  w[ 6] = 0.570135236262479578347113482275E+00;
2279  w[ 7] = 0.260492310264161129233396139765E+00;
2280  w[ 8] = 0.516079856158839299918734423606E-01;
2281  w[ 9] = 0.390539058462906185999438432620E-02;
2282  w[10] = 0.857368704358785865456906323153E-04;
2283  w[11] = 0.265855168435630160602311400877E-06;
2284  }
2285  else if (n==13) {
2286  w[ 0] = 0.482573185007313108834997332342E-07;
2287  w[ 1] = 0.204303604027070731248669432937E-04;
2288  w[ 2] = 0.120745999271938594730924899224E-02;
2289  w[ 3] = 0.208627752961699392166033805050E-01;
2290  w[ 4] = 0.140323320687023437762792268873E+00;
2291  w[ 5] = 0.421616296898543221746893558568E+00;
2292  w[ 6] = 0.604393187921161642342099068579E+00;
2293  w[ 7] = 0.421616296898543221746893558568E+00;
2294  w[ 8] = 0.140323320687023437762792268873E+00;
2295  w[ 9] = 0.208627752961699392166033805050E-01;
2296  w[10] = 0.120745999271938594730924899224E-02;
2297  w[11] = 0.204303604027070731248669432937E-04;
2298  w[12] = 0.482573185007313108834997332342E-07;
2299  }
2300  else if (n==14) {
2301  w[ 0] = 0.862859116812515794532041783429E-08;
2302  w[ 1] = 0.471648435501891674887688950105E-05;
2303  w[ 2] = 0.355092613551923610483661076691E-03;
2304  w[ 3] = 0.785005472645794431048644334608E-02;
2305  w[ 4] = 0.685055342234652055387163312367E-01;
2306  w[ 5] = 0.273105609064246603352569187026E+00;
2307  w[ 6] = 0.536405909712090149794921296776E+00;
2308  w[ 7] = 0.536405909712090149794921296776E+00;
2309  w[ 8] = 0.273105609064246603352569187026E+00;
2310  w[ 9] = 0.685055342234652055387163312367E-01;
2311  w[10] = 0.785005472645794431048644334608E-02;
2312  w[11] = 0.355092613551923610483661076691E-03;
2313  w[12] = 0.471648435501891674887688950105E-05;
2314  w[13] = 0.862859116812515794532041783429E-08;
2315  }
2316  else if (n==15) {
2317  w[ 0] = 0.152247580425351702016062666965E-08;
2318  w[ 1] = 0.105911554771106663577520791055E-05;
2319  w[ 2] = 0.100004441232499868127296736177E-03;
2320  w[ 3] = 0.277806884291277589607887049229E-02;
2321  w[ 4] = 0.307800338725460822286814158758E-01;
2322  w[ 5] = 0.158488915795935746883839384960E+00;
2323  w[ 6] = 0.412028687498898627025891079568E+00;
2324  w[ 7] = 0.564100308726417532852625797340E+00;
2325  w[ 8] = 0.412028687498898627025891079568E+00;
2326  w[ 9] = 0.158488915795935746883839384960E+00;
2327  w[10] = 0.307800338725460822286814158758E-01;
2328  w[11] = 0.277806884291277589607887049229E-02;
2329  w[12] = 0.100004441232499868127296736177E-03;
2330  w[13] = 0.105911554771106663577520791055E-05;
2331  w[14] = 0.152247580425351702016062666965E-08;
2332  }
2333  else if (n==16) {
2334  w[ 0] = 0.265480747401118224470926366050E-09;
2335  w[ 1] = 0.232098084486521065338749423185E-06;
2336  w[ 2] = 0.271186009253788151201891432244E-04;
2337  w[ 3] = 0.932284008624180529914277305537E-03;
2338  w[ 4] = 0.128803115355099736834642999312E-01;
2339  w[ 5] = 0.838100413989858294154207349001E-01;
2340  w[ 6] = 0.280647458528533675369463335380E+00;
2341  w[ 7] = 0.507929479016613741913517341791E+00;
2342  w[ 8] = 0.507929479016613741913517341791E+00;
2343  w[ 9] = 0.280647458528533675369463335380E+00;
2344  w[10] = 0.838100413989858294154207349001E-01;
2345  w[11] = 0.128803115355099736834642999312E-01;
2346  w[12] = 0.932284008624180529914277305537E-03;
2347  w[13] = 0.271186009253788151201891432244E-04;
2348  w[14] = 0.232098084486521065338749423185E-06;
2349  w[15] = 0.265480747401118224470926366050E-09;
2350  }
2351  else if (n==17) {
2352  w[ 0] = 0.458057893079863330580889281222E-10;
2353  w[ 1] = 0.497707898163079405227863353715E-07;
2354  w[ 2] = 0.711228914002130958353327376218E-05;
2355  w[ 3] = 0.298643286697753041151336643059E-03;
2356  w[ 4] = 0.506734995762753791170069495879E-02;
2357  w[ 5] = 0.409200341495762798094994877854E-01;
2358  w[ 6] = 0.172648297670097079217645196219E+00;
2359  w[ 7] = 0.401826469470411956577635085257E+00;
2360  w[ 8] = 0.530917937624863560331883103379E+00;
2361  w[ 9] = 0.401826469470411956577635085257E+00;
2362  w[10] = 0.172648297670097079217645196219E+00;
2363  w[11] = 0.409200341495762798094994877854E-01;
2364  w[12] = 0.506734995762753791170069495879E-02;
2365  w[13] = 0.298643286697753041151336643059E-03;
2366  w[14] = 0.711228914002130958353327376218E-05;
2367  w[15] = 0.497707898163079405227863353715E-07;
2368  w[16] = 0.458057893079863330580889281222E-10;
2369  }
2370  else if (n==18) {
2371  w[ 0] = 0.782819977211589102925147471012E-11;
2372  w[ 1] = 0.104672057957920824443559608435E-07;
2373  w[ 2] = 0.181065448109343040959702385911E-05;
2374  w[ 3] = 0.918112686792940352914675407371E-04;
2375  w[ 4] = 0.188852263026841789438175325426E-02;
2376  w[ 5] = 0.186400423875446519219315221973E-01;
2377  w[ 6] = 0.973017476413154293308537234155E-01;
2378  w[ 7] = 0.284807285669979578595606820713E+00;
2379  w[ 8] = 0.483495694725455552876410522141E+00;
2380  w[ 9] = 0.483495694725455552876410522141E+00;
2381  w[10] = 0.284807285669979578595606820713E+00;
2382  w[11] = 0.973017476413154293308537234155E-01;
2383  w[12] = 0.186400423875446519219315221973E-01;
2384  w[13] = 0.188852263026841789438175325426E-02;
2385  w[14] = 0.918112686792940352914675407371E-04;
2386  w[15] = 0.181065448109343040959702385911E-05;
2387  w[16] = 0.104672057957920824443559608435E-07;
2388  w[17] = 0.782819977211589102925147471012E-11;
2389  }
2390  else if (n==19) {
2391  w[ 0] = 0.132629709449851575185289154385E-11;
2392  w[ 1] = 0.216305100986355475019693077221E-08;
2393  w[ 2] = 0.448824314722312295179447915594E-06;
2394  w[ 3] = 0.272091977631616257711941025214E-04;
2395  w[ 4] = 0.670877521407181106194696282100E-03;
2396  w[ 5] = 0.798886677772299020922211491861E-02;
2397  w[ 6] = 0.508103869090520673569908110358E-01;
2398  w[ 7] = 0.183632701306997074156148485766E+00;
2399  w[ 8] = 0.391608988613030244504042313621E+00;
2400  w[ 9] = 0.502974888276186530840731361096E+00;
2401  w[10] = 0.391608988613030244504042313621E+00;
2402  w[11] = 0.183632701306997074156148485766E+00;
2403  w[12] = 0.508103869090520673569908110358E-01;
2404  w[13] = 0.798886677772299020922211491861E-02;
2405  w[14] = 0.670877521407181106194696282100E-03;
2406  w[15] = 0.272091977631616257711941025214E-04;
2407  w[16] = 0.448824314722312295179447915594E-06;
2408  w[17] = 0.216305100986355475019693077221E-08;
2409  w[18] = 0.132629709449851575185289154385E-11;
2410  }
2411  else if (n==20) {
2412  w[ 0] = 0.222939364553415129252250061603E-12;
2413  w[ 1] = 0.439934099227318055362885145547E-09;
2414  w[ 2] = 0.108606937076928169399952456345E-06;
2415  w[ 3] = 0.780255647853206369414599199965E-05;
2416  w[ 4] = 0.228338636016353967257145917963E-03;
2417  w[ 5] = 0.324377334223786183218324713235E-02;
2418  w[ 6] = 0.248105208874636108821649525589E-01;
2419  w[ 7] = 0.109017206020023320013755033535E+00;
2420  w[ 8] = 0.286675505362834129719659706228E+00;
2421  w[ 9] = 0.462243669600610089650328639861E+00;
2422  w[10] = 0.462243669600610089650328639861E+00;
2423  w[11] = 0.286675505362834129719659706228E+00;
2424  w[12] = 0.109017206020023320013755033535E+00;
2425  w[13] = 0.248105208874636108821649525589E-01;
2426  w[14] = 0.324377334223786183218324713235E-02;
2427  w[15] = 0.228338636016353967257145917963E-03;
2428  w[16] = 0.780255647853206369414599199965E-05;
2429  w[17] = 0.108606937076928169399952456345E-06;
2430  w[18] = 0.439934099227318055362885145547E-09;
2431  w[19] = 0.222939364553415129252250061603E-12;
2432  }
2433  else {
2434  std::cerr << "\n";
2435  std::cerr << "HERMITE_LOOKUP_WEIGHTS - Fatal error!\n";
2436  std::cerr << " Illegal value of N = " << n << "\n";
2437  std::cerr << " Legal values are 1 through 20.\n";
2438  std::exit(1);
2439  }
2440 
2441  return;
2442 }
2443 
2444 //****************************************************************************
2445 template<class Scalar>
2446 void IntrepidBurkardtRules::imtqlx ( int n, Scalar d[], Scalar e[], Scalar z[] )
2447 //****************************************************************************
2448 //
2449 // Purpose:
2450 //
2451 // IMTQLX diagonalizes a symmetric tridiagonal matrix.
2452 //
2453 // Discussion:
2454 //
2455 // This routine is a slightly modified version of the EISPACK routine to
2456 // perform the implicit QL algorithm on a symmetric tridiagonal matrix.
2457 //
2458 // The authors thank the authors of EISPACK for permission to use this
2459 // routine.
2460 //
2461 // It has been modified to produce the product Q' * Z, where Z is an input
2462 // vector and Q is the orthogonal matrix diagonalizing the input matrix.
2463 // The changes consist (essentially) of applying the orthogonal transformations
2464 // directly to Z as they are generated.
2465 //
2466 // Licensing:
2467 //
2468 // This code is distributed under the GNU LGPL license.
2469 //
2470 // Modified:
2471 //
2472 // 08 January 2010
2473 //
2474 // Author:
2475 //
2476 // Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky.
2477 // C++ version by John Burkardt.
2478 //
2479 // Reference:
2480 //
2481 // Sylvan Elhay, Jaroslav Kautsky,
2482 // Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of
2483 // Interpolatory Quadrature,
2484 // ACM Transactions on Mathematical Software,
2485 // Volume 13, Number 4, December 1987, pages 399-415.
2486 //
2487 // Roger Martin, James Wilkinson,
2488 // The Implicit QL Algorithm,
2489 // Numerische Mathematik,
2490 // Volume 12, Number 5, December 1968, pages 377-383.
2491 //
2492 // Parameters:
2493 //
2494 // Input, int N, the order of the matrix.
2495 //
2496 // Input/output, Scalar D(N), the diagonal entries of the matrix.
2497 // On output, the information in D has been overwritten.
2498 //
2499 // Input/output, Scalar E(N), the subdiagonal entries of the
2500 // matrix, in entries E(1) through E(N-1). On output, the information in
2501 // E has been overwritten.
2502 //
2503 // Input/output, Scalar Z(N). On input, a vector. On output,
2504 // the value of Q' * Z, where Q is the matrix that diagonalizes the
2505 // input symmetric tridiagonal matrix.
2506 //
2507 {
2508  Scalar b = 0, c = 0, f = 0, g = 0, p = 0, r = 0, s = 0;
2509  int i = 0, ii = 0, j = 0, k = 0, l = 0, mml = 0, m = 0, itn = 30;
2510  Scalar prec = IntrepidBurkardtRules::r8_epsilon(1.0); //2.22E-16;
2511 
2512  if (n==1) {
2513  return;
2514  }
2515 
2516  e[n-1] = 0.0;
2517 
2518  for (l=1;l<=n;l++) {
2519  j = 0;
2520  for ( ; ; ) {
2521  for (m=l;m<=n;m++) {
2522  if (m==n) {
2523  break;
2524  }
2525  if (std::abs(e[m-1])<=prec*(std::abs(d[m-1])+std::abs(d[m]))) {
2526  break;
2527  }
2528  }
2529  p = d[l-1];
2530  if (m==l) {
2531  break;
2532  }
2533  if (itn<=j) {
2534  std::cerr << "\n";
2535  std::cerr << "IMTQLX - Fatal error!\n";
2536  std::cerr << " Iteration limit exceeded\n";
2537  std::exit(1);
2538  }
2539  j = j+1;
2540  g = (d[l]-p)/(2.0*e[l-1]);
2541  r = std::sqrt(g*g+1.0);
2542  g = d[m-1]-p+e[l-1]/(g+std::abs(r)*IntrepidBurkardtRules::r8_sign(g));
2543  s = 1.0;
2544  c = 1.0;
2545  p = 0.0;
2546  mml = m-l;
2547 
2548  for (ii=1;ii<=mml;ii++) {
2549  i = m-ii;
2550  f = s*e[i-1];
2551  b = c*e[i-1];
2552 
2553  if (std::abs(g)<=std::abs(f)) {
2554  c = g/f;
2555  r = std::sqrt(c*c+1.0);
2556  e[i] = f*r;
2557  s = 1.0/r;
2558  c = c*s;
2559  }
2560  else {
2561  s = f/g;
2562  r = std::sqrt(s*s+1.0);
2563  e[i] = g*r;
2564  c = 1.0/r;
2565  s = s*c;
2566  }
2567  g = d[i]-p;
2568  r = (d[i-1]-g)*s+2.0*c*b;
2569  p = s*r;
2570  d[i] = g+p;
2571  g = c*r-b;
2572  f = z[i];
2573  z[i] = s*z[i-1]+c*f;
2574  z[i-1] = c*z[i-1]-s*f;
2575  }
2576  d[l-1] = d[l-1]-p;
2577  e[l-1] = g;
2578  e[m-1] = 0.0;
2579  }
2580  }
2581 //
2582 // Sorting.
2583 //
2584  for (ii=2;ii<=m;ii++) {
2585  i = ii-1;
2586  k = i;
2587  p = d[i-1];
2588 
2589  for (j=ii;j<=n;j++) {
2590  if (d[j-1]<p) {
2591  k = j;
2592  p = d[j-1];
2593  }
2594  }
2595 
2596  if (k!=i) {
2597  d[k-1] = d[i-1];
2598  d[i-1] = p;
2599  p = z[i-1];
2600  z[i-1] = z[k-1];
2601  z[k-1] = p;
2602  }
2603  }
2604  return;
2605 }
2606 
2607 //****************************************************************************
2608 template<class Scalar>
2609 void IntrepidBurkardtRules::laguerre_compute ( int n, Scalar x[], Scalar w[] )
2610 //****************************************************************************
2611 //
2612 // Purpose:
2613 //
2614 // LAGUERRE_COMPUTE: Laguerre quadrature rule by the Elhay-Kautsky method.
2615 //
2616 // Licensing:
2617 //
2618 // This code is distributed under the GNU LGPL license.
2619 //
2620 // Modified:
2621 //
2622 // 23 April 2011
2623 //
2624 // Author:
2625 //
2626 // Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky.
2627 // C++ version by John Burkardt.
2628 //
2629 // Reference:
2630 //
2631 // Sylvan Elhay, Jaroslav Kautsky,
2632 // Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of
2633 // Interpolatory Quadrature,
2634 // ACM Transactions on Mathematical Software,
2635 // Volume 13, Number 4, December 1987, pages 399-415.
2636 //
2637 // Parameters:
2638 //
2639 // Input, int N, the order.
2640 //
2641 // Output, Scalar X[N], the abscissas.
2642 //
2643 // Output, Scalar W[N], the weights.
2644 //
2645 {
2646  Scalar *bj;
2647  int i;
2648  Scalar zemu;
2649 //
2650 // Define the zero-th moment.
2651 //
2652  zemu = 1.0;
2653 //
2654 // Define the Jacobi matrix.
2655 //
2656  bj = new Scalar[n];
2657 
2658  for (i=0;i<n;i++) {
2659  bj[i] = (Scalar)(i+1);
2660  }
2661 
2662  for (i=0;i<n;i++) {
2663  x[i] = (Scalar)(2*i+1);
2664  }
2665 
2666  w[0] = std::sqrt(zemu);
2667 
2668  for (i=1;i<n;i++) {
2669  w[i] = 0.0;
2670  }
2671 //
2672 // Diagonalize the Jacobi matrix.
2673 //
2674  IntrepidBurkardtRules::imtqlx(n,x,bj,w);
2675 
2676  for (i=0;i<n;i++) {
2677  w[i] = w[i]*w[i];
2678  }
2679 
2680  delete [] bj;
2681 
2682  return;
2683 }
2684 
2685 //****************************************************************************
2686 template<class Scalar>
2688 //****************************************************************************
2689 //
2690 // Purpose:
2691 //
2692 // LAGUERRE_COMPUTE_POINTS computes Laguerre quadrature points.
2693 //
2694 // Licensing:
2695 //
2696 // This code is distributed under the GNU LGPL license.
2697 //
2698 // Modified:
2699 //
2700 // 13 June 2009
2701 //
2702 // Author:
2703 //
2704 // John Burkardt
2705 //
2706 // Parameters:
2707 //
2708 // Input, int ORDER, the order.
2709 //
2710 // Output, Scalar X[ORDER], the abscissas.
2711 //
2712 {
2713  Scalar *w; w = new Scalar[order];
2715  delete [] w;
2716 
2717  return;
2718 }
2719 
2720 //****************************************************************************
2721 template<class Scalar>
2723 //****************************************************************************
2724 //
2725 // Purpose:
2726 //
2727 // LAGUERRE_COMPUTE_WEIGHTS computes Laguerre quadrature weights.
2728 //
2729 // Licensing:
2730 //
2731 // This code is distributed under the GNU LGPL license.
2732 //
2733 // Modified:
2734 //
2735 // 13 June 2009
2736 //
2737 // Author:
2738 //
2739 // John Burkardt
2740 //
2741 // Parameters:
2742 //
2743 // Input, int ORDER, the order.
2744 //
2745 // Output, Scalar W[ORDER], the weights.
2746 //
2747 {
2748  Scalar *x; x = new Scalar[order];
2750  delete [] x;
2751 
2752  return;
2753 }
2754 
2755 //****************************************************************************
2756 template<class Scalar>
2757 void IntrepidBurkardtRules::laguerre_lookup ( int n, Scalar x[], Scalar w[] )
2758 //****************************************************************************
2759 //
2760 // Purpose:
2761 //
2762 // LAGUERRE_LOOKUP looks up abscissas and weights for Laguerre quadrature.
2763 //
2764 // Discussion:
2765 //
2766 // The abscissas are the zeroes of the Laguerre polynomial L(N)(X).
2767 //
2768 // The integral:
2769 //
2770 // Integral ( 0 <= X < +oo ) exp ( -X ) * F(X) dX
2771 //
2772 // The quadrature rule:
2773 //
2774 // Sum ( 1 <= I <= N ) W(I) * f ( X(I) )
2775 //
2776 // The integral:
2777 //
2778 // Integral ( 0 <= X < +oo ) F(X) dX
2779 //
2780 // The quadrature rule:
2781 //
2782 // Sum ( 1 <= I <= N ) W(I) * exp ( X(I) ) * f ( X(I) )
2783 //
2784 // Mathematica can numerically estimate the abscissas for the
2785 // n-th order polynomial to p digits of precision by the command:
2786 //
2787 // NSolve [ LaguerreL[n,x] == 0, x, p ]
2788 //
2789 // Licensing:
2790 //
2791 // This code is distributed under the GNU LGPL license.
2792 //
2793 // Modified:
2794 //
2795 // 27 April 2010
2796 //
2797 // Author:
2798 //
2799 // John Burkardt
2800 //
2801 // Reference:
2802 //
2803 // Milton Abramowitz, Irene Stegun,
2804 // Handbook of Mathematical Functions,
2805 // National Bureau of Standards, 1964,
2806 // ISBN: 0-486-61272-4,
2807 // LC: QA47.A34.
2808 //
2809 // Vladimir Krylov,
2810 // Approximate Calculation of Integrals,
2811 // Dover, 2006,
2812 // ISBN: 0486445798,
2813 // LC: QA311.K713.
2814 //
2815 // Arthur Stroud, Don Secrest,
2816 // Gaussian Quadrature Formulas,
2817 // Prentice Hall, 1966,
2818 // LC: QA299.4G3S7.
2819 //
2820 // Stephen Wolfram,
2821 // The Mathematica Book,
2822 // Fourth Edition,
2823 // Cambridge University Press, 1999,
2824 // ISBN: 0-521-64314-7,
2825 // LC: QA76.95.W65.
2826 //
2827 // Daniel Zwillinger, editor,
2828 // CRC Standard Mathematical Tables and Formulae,
2829 // 30th Edition,
2830 // CRC Press, 1996,
2831 // ISBN: 0-8493-2479-3.
2832 //
2833 // Parameters:
2834 //
2835 // Input, int N, the order.
2836 // N must be between 1 and 20.
2837 //
2838 // Output, Scalar X[N], the abscissas.
2839 //
2840 // Output, Scalar W[N], the weights.
2841 //
2842 {
2845 
2846  return;
2847 }
2848 
2849 //****************************************************************************
2850 template<class Scalar>
2852 //****************************************************************************
2853 //
2854 // Purpose:
2855 //
2856 // LAGUERRE_LOOKUP_POINTS looks up abscissas for Laguerre quadrature.
2857 //
2858 // Licensing:
2859 //
2860 // This code is distributed under the GNU LGPL license.
2861 //
2862 // Modified:
2863 //
2864 // 27 April 2010
2865 //
2866 // Author:
2867 //
2868 // John Burkardt
2869 //
2870 // Reference:
2871 //
2872 // Milton Abramowitz, Irene Stegun,
2873 // Handbook of Mathematical Functions,
2874 // National Bureau of Standards, 1964,
2875 // ISBN: 0-486-61272-4,
2876 // LC: QA47.A34.
2877 //
2878 // Vladimir Krylov,
2879 // Approximate Calculation of Integrals,
2880 // Dover, 2006,
2881 // ISBN: 0486445798,
2882 // LC: QA311.K713.
2883 //
2884 // Arthur Stroud, Don Secrest,
2885 // Gaussian Quadrature Formulas,
2886 // Prentice Hall, 1966,
2887 // LC: QA299.4G3S7.
2888 //
2889 // Stephen Wolfram,
2890 // The Mathematica Book,
2891 // Fourth Edition,
2892 // Cambridge University Press, 1999,
2893 // ISBN: 0-521-64314-7,
2894 // LC: QA76.95.W65.
2895 //
2896 // Daniel Zwillinger, editor,
2897 // CRC Standard Mathematical Tables and Formulae,
2898 // 30th Edition,
2899 // CRC Press, 1996,
2900 // ISBN: 0-8493-2479-3.
2901 //
2902 // Parameters:
2903 //
2904 // Input, int N, the order.
2905 // N must be between 1 and 20.
2906 //
2907 // Output, Scalar X[N], the abscissas.
2908 //
2909 {
2910  if (n==1) {
2911  x[ 0] = 1.00000000000000000000000000000E+00;
2912  }
2913  else if (n==2) {
2914  x[ 0] = 0.585786437626904951198311275790E+00;
2915  x[ 1] = 3.41421356237309504880168872421E+00;
2916  }
2917  else if (n==3) {
2918  x[ 0] = 0.415774556783479083311533873128E+00;
2919  x[ 1] = 2.29428036027904171982205036136E+00;
2920  x[ 2] = 6.28994508293747919686641576551E+00;
2921  }
2922  else if (n==4) {
2923  x[ 0] = 0.322547689619392311800361459104E+00;
2924  x[ 1] = 1.74576110115834657568681671252E+00;
2925  x[ 2] = 4.53662029692112798327928538496E+00;
2926  x[ 3] = 9.39507091230113312923353644342E+00;
2927  }
2928  else if (n==5) {
2929  x[ 0] = 0.263560319718140910203061943361E+00;
2930  x[ 1] = 1.41340305910651679221840798019E+00;
2931  x[ 2] = 3.59642577104072208122318658878E+00;
2932  x[ 3] = 7.08581000585883755692212418111E+00;
2933  x[ 4] = 12.6408008442757826594332193066E+00;
2934  }
2935  else if (n==6) {
2936  x[ 0] = 0.222846604179260689464354826787E+00;
2937  x[ 1] = 1.18893210167262303074315092194E+00;
2938  x[ 2] = 2.99273632605931407769132528451E+00;
2939  x[ 3] = 5.77514356910451050183983036943E+00;
2940  x[ 4] = 9.83746741838258991771554702994E+00;
2941  x[ 5] = 15.9828739806017017825457915674E+00;
2942  }
2943  else if (n==7) {
2944  x[ 0] = 0.193043676560362413838247885004E+00;
2945  x[ 1] = 1.02666489533919195034519944317E+00;
2946  x[ 2] = 2.56787674495074620690778622666E+00;
2947  x[ 3] = 4.90035308452648456810171437810E+00;
2948  x[ 4] = 8.18215344456286079108182755123E+00;
2949  x[ 5] = 12.7341802917978137580126424582E+00;
2950  x[ 6] = 19.3957278622625403117125820576E+00;
2951  }
2952  else if (n==8) {
2953  x[ 0] = 0.170279632305100999788861856608E+00;
2954  x[ 1] = 0.903701776799379912186020223555E+00;
2955  x[ 2] = 2.25108662986613068930711836697E+00;
2956  x[ 3] = 4.26670017028765879364942182690E+00;
2957  x[ 4] = 7.04590540239346569727932548212E+00;
2958  x[ 5] = 10.7585160101809952240599567880E+00;
2959  x[ 6] = 15.7406786412780045780287611584E+00;
2960  x[ 7] = 22.8631317368892641057005342974E+00;
2961  }
2962  else if (n==9) {
2963  x[ 0] = 0.152322227731808247428107073127E+00;
2964  x[ 1] = 0.807220022742255847741419210952E+00;
2965  x[ 2] = 2.00513515561934712298303324701E+00;
2966  x[ 3] = 3.78347397333123299167540609364E+00;
2967  x[ 4] = 6.20495677787661260697353521006E+00;
2968  x[ 5] = 9.37298525168757620180971073215E+00;
2969  x[ 6] = 13.4662369110920935710978818397E+00;
2970  x[ 7] = 18.8335977889916966141498992996E+00;
2971  x[ 8] = 26.3740718909273767961410072937E+00;
2972  }
2973  else if (n==10) {
2974  x[ 0] = 0.137793470540492430830772505653E+00;
2975  x[ 1] = 0.729454549503170498160373121676E+00;
2976  x[ 2] = 1.80834290174031604823292007575E+00;
2977  x[ 3] = 3.40143369785489951448253222141E+00;
2978  x[ 4] = 5.55249614006380363241755848687E+00;
2979  x[ 5] = 8.33015274676449670023876719727E+00;
2980  x[ 6] = 11.8437858379000655649185389191E+00;
2981  x[ 7] = 16.2792578313781020995326539358E+00;
2982  x[ 8] = 21.9965858119807619512770901956E+00;
2983  x[ 9] = 29.9206970122738915599087933408E+00;
2984  }
2985  else if (n==11) {
2986  x[ 0] = 0.125796442187967522675794577516E+00;
2987  x[ 1] = 0.665418255839227841678127839420E+00;
2988  x[ 2] = 1.64715054587216930958700321365E+00;
2989  x[ 3] = 3.09113814303525495330195934259E+00;
2990  x[ 4] = 5.02928440157983321236999508366E+00;
2991  x[ 5] = 7.50988786380661681941099714450E+00;
2992  x[ 6] = 10.6059509995469677805559216457E+00;
2993  x[ 7] = 14.4316137580641855353200450349E+00;
2994  x[ 8] = 19.1788574032146786478174853989E+00;
2995  x[ 9] = 25.2177093396775611040909447797E+00;
2996  x[10] = 33.4971928471755372731917259395E+00;
2997  }
2998  else if (n==12) {
2999  x[ 0] = 0.115722117358020675267196428240E+00;
3000  x[ 1] = 0.611757484515130665391630053042E+00;
3001  x[ 2] = 1.51261026977641878678173792687E+00;
3002  x[ 3] = 2.83375133774350722862747177657E+00;
3003  x[ 4] = 4.59922763941834848460572922485E+00;
3004  x[ 5] = 6.84452545311517734775433041849E+00;
3005  x[ 6] = 9.62131684245686704391238234923E+00;
3006  x[ 7] = 13.0060549933063477203460524294E+00;
3007  x[ 8] = 17.1168551874622557281840528008E+00;
3008  x[ 9] = 22.1510903793970056699218950837E+00;
3009  x[10] = 28.4879672509840003125686072325E+00;
3010  x[11] = 37.0991210444669203366389142764E+00;
3011  }
3012  else if (n==13) {
3013  x[ 0] = 0.107142388472252310648493376977E+00;
3014  x[ 1] = 0.566131899040401853406036347177E+00;
3015  x[ 2] = 1.39856433645101971792750259921E+00;
3016  x[ 3] = 2.61659710840641129808364008472E+00;
3017  x[ 4] = 4.23884592901703327937303389926E+00;
3018  x[ 5] = 6.29225627114007378039376523025E+00;
3019  x[ 6] = 8.81500194118697804733348868036E+00;
3020  x[ 7] = 11.8614035888112425762212021880E+00;
3021  x[ 8] = 15.5107620377037527818478532958E+00;
3022  x[ 9] = 19.8846356638802283332036594634E+00;
3023  x[10] = 25.1852638646777580842970297823E+00;
3024  x[11] = 31.8003863019472683713663283526E+00;
3025  x[12] = 40.7230086692655795658979667001E+00;
3026  }
3027  else if (n==14) {
3028  x[ 0] = 0.0997475070325975745736829452514E+00;
3029  x[ 1] = 0.526857648851902896404583451502E+00;
3030  x[ 2] = 1.30062912125149648170842022116E+00;
3031  x[ 3] = 2.43080107873084463616999751038E+00;
3032  x[ 4] = 3.93210282229321888213134366778E+00;
3033  x[ 5] = 5.82553621830170841933899983898E+00;
3034  x[ 6] = 8.14024014156514503005978046052E+00;
3035  x[ 7] = 10.9164995073660188408130510904E+00;
3036  x[ 8] = 14.2108050111612886831059780825E+00;
3037  x[ 9] = 18.1048922202180984125546272083E+00;
3038  x[10] = 22.7233816282696248232280886985E+00;
3039  x[11] = 28.2729817232482056954158923218E+00;
3040  x[12] = 35.1494436605924265828643121364E+00;
3041  x[13] = 44.3660817111174230416312423666E+00;
3042  }
3043  else if (n==15) {
3044  x[ 0] = 0.0933078120172818047629030383672E+00;
3045  x[ 1] = 0.492691740301883908960101791412E+00;
3046  x[ 2] = 1.21559541207094946372992716488E+00;
3047  x[ 3] = 2.26994952620374320247421741375E+00;
3048  x[ 4] = 3.66762272175143727724905959436E+00;
3049  x[ 5] = 5.42533662741355316534358132596E+00;
3050  x[ 6] = 7.56591622661306786049739555812E+00;
3051  x[ 7] = 10.1202285680191127347927394568E+00;
3052  x[ 8] = 13.1302824821757235640991204176E+00;
3053  x[ 9] = 16.6544077083299578225202408430E+00;
3054  x[10] = 20.7764788994487667729157175676E+00;
3055  x[11] = 25.6238942267287801445868285977E+00;
3056  x[12] = 31.4075191697539385152432196202E+00;
3057  x[13] = 38.5306833064860094162515167595E+00;
3058  x[14] = 48.0260855726857943465734308508E+00;
3059  }
3060  else if (n==16) {
3061  x[ 0] = 0.0876494104789278403601980973401E+00;
3062  x[ 1] = 0.462696328915080831880838260664E+00;
3063  x[ 2] = 1.14105777483122685687794501811E+00;
3064  x[ 3] = 2.12928364509838061632615907066E+00;
3065  x[ 4] = 3.43708663389320664523510701675E+00;
3066  x[ 5] = 5.07801861454976791292305830814E+00;
3067  x[ 6] = 7.07033853504823413039598947080E+00;
3068  x[ 7] = 9.43831433639193878394724672911E+00;
3069  x[ 8] = 12.2142233688661587369391246088E+00;
3070  x[ 9] = 15.4415273687816170767647741622E+00;
3071  x[10] = 19.1801568567531348546631409497E+00;
3072  x[11] = 23.5159056939919085318231872752E+00;
3073  x[12] = 28.5787297428821403675206137099E+00;
3074  x[13] = 34.5833987022866258145276871778E+00;
3075  x[14] = 41.9404526476883326354722330252E+00;
3076  x[15] = 51.7011603395433183643426971197E+00;
3077  }
3078  else if (n==17) {
3079  x[ 0] = 0.0826382147089476690543986151980E+00;
3080  x[ 1] = 0.436150323558710436375959029847E+00;
3081  x[ 2] = 1.07517657751142857732980316755E+00;
3082  x[ 3] = 2.00519353164923224070293371933E+00;
3083  x[ 4] = 3.23425612404744376157380120696E+00;
3084  x[ 5] = 4.77351351370019726480932076262E+00;
3085  x[ 6] = 6.63782920536495266541643929703E+00;
3086  x[ 7] = 8.84668551116980005369470571184E+00;
3087  x[ 8] = 11.4255293193733525869726151469E+00;
3088  x[ 9] = 14.4078230374813180021982874959E+00;
3089  x[10] = 17.8382847307011409290658752412E+00;
3090  x[11] = 21.7782682577222653261749080522E+00;
3091  x[12] = 26.3153178112487997766149598369E+00;
3092  x[13] = 31.5817716804567331343908517497E+00;
3093  x[14] = 37.7960938374771007286092846663E+00;
3094  x[15] = 45.3757165339889661829258363215E+00;
3095  x[16] = 55.3897517898396106640900199790E+00;
3096  }
3097  else if (n==18) {
3098  x[ 0] = 0.0781691666697054712986747615334E+00;
3099  x[ 1] = 0.412490085259129291039101536536E+00;
3100  x[ 2] = 1.01652017962353968919093686187E+00;
3101  x[ 3] = 1.89488850996976091426727831954E+00;
3102  x[ 4] = 3.05435311320265975115241130719E+00;
3103  x[ 5] = 4.50420553888989282633795571455E+00;
3104  x[ 6] = 6.25672507394911145274209116326E+00;
3105  x[ 7] = 8.32782515660563002170470261564E+00;
3106  x[ 8] = 10.7379900477576093352179033397E+00;
3107  x[ 9] = 13.5136562075550898190863812108E+00;
3108  x[10] = 16.6893062819301059378183984163E+00;
3109  x[11] = 20.3107676262677428561313764553E+00;
3110  x[12] = 24.4406813592837027656442257980E+00;
3111  x[13] = 29.1682086625796161312980677805E+00;
3112  x[14] = 34.6279270656601721454012429438E+00;
3113  x[15] = 41.0418167728087581392948614284E+00;
3114  x[16] = 48.8339227160865227486586093290E+00;
3115  x[17] = 59.0905464359012507037157810181E+00;
3116  }
3117  else if (n==19) {
3118  x[ 0] = 0.0741587837572050877131369916024E+00;
3119  x[ 1] = 0.391268613319994607337648350299E+00;
3120  x[ 2] = 0.963957343997958058624878377130E+00;
3121  x[ 3] = 1.79617558206832812557725825252E+00;
3122  x[ 4] = 2.89365138187378399116494713237E+00;
3123  x[ 5] = 4.26421553962776647436040018167E+00;
3124  x[ 6] = 5.91814156164404855815360191408E+00;
3125  x[ 7] = 7.86861891533473373105668358176E+00;
3126  x[ 8] = 10.1324237168152659251627415800E+00;
3127  x[ 9] = 12.7308814638423980045092979656E+00;
3128  x[10] = 15.6912783398358885454136069861E+00;
3129  x[11] = 19.0489932098235501532136429732E+00;
3130  x[12] = 22.8508497608294829323930586693E+00;
3131  x[13] = 27.1606693274114488789963947149E+00;
3132  x[14] = 32.0691222518622423224362865906E+00;
3133  x[15] = 37.7129058012196494770647508283E+00;
3134  x[16] = 44.3173627958314961196067736013E+00;
3135  x[17] = 52.3129024574043831658644222420E+00;
3136  x[18] = 62.8024231535003758413504690673E+00;
3137  }
3138  else if (n==20) {
3139  x[ 0] = 0.0705398896919887533666890045842E+00;
3140  x[ 1] = 0.372126818001611443794241388761E+00;
3141  x[ 2] = 0.916582102483273564667716277074E+00;
3142  x[ 3] = 1.70730653102834388068768966741E+00;
3143  x[ 4] = 2.74919925530943212964503046049E+00;
3144  x[ 5] = 4.04892531385088692237495336913E+00;
3145  x[ 6] = 5.61517497086161651410453988565E+00;
3146  x[ 7] = 7.45901745367106330976886021837E+00;
3147  x[ 8] = 9.59439286958109677247367273428E+00;
3148  x[ 9] = 12.0388025469643163096234092989E+00;
3149  x[10] = 14.8142934426307399785126797100E+00;
3150  x[11] = 17.9488955205193760173657909926E+00;
3151  x[12] = 21.4787882402850109757351703696E+00;
3152  x[13] = 25.4517027931869055035186774846E+00;
3153  x[14] = 29.9325546317006120067136561352E+00;
3154  x[15] = 35.0134342404790000062849359067E+00;
3155  x[16] = 40.8330570567285710620295677078E+00;
3156  x[17] = 47.6199940473465021399416271529E+00;
3157  x[18] = 55.8107957500638988907507734445E+00;
3158  x[19] = 66.5244165256157538186403187915E+00;
3159  }
3160  else {
3161  std::cerr << "\n";
3162  std::cerr << "LAGUERRE_LOOKUP_POINTS - Fatal error!\n";
3163  std::cerr << " Illegal value of N = " << n << "\n";
3164  std::cerr << " Legal values are 1 through 20.\n";
3165  std::exit(1);
3166  }
3167 
3168  return;
3169 }
3170 
3171 //****************************************************************************
3172 template<class Scalar>
3174 //****************************************************************************
3175 //
3176 // Purpose:
3177 //
3178 // LAGUERRE_LOOKUP_WEIGHTS looks up weights for Laguerre quadrature.
3179 //
3180 // Licensing:
3181 //
3182 // This code is distributed under the GNU LGPL license.
3183 //
3184 // Modified:
3185 //
3186 // 27 April 2010
3187 //
3188 // Author:
3189 //
3190 // John Burkardt
3191 //
3192 // Reference:
3193 //
3194 // Milton Abramowitz, Irene Stegun,
3195 // Handbook of Mathematical Functions,
3196 // National Bureau of Standards, 1964,
3197 // ISBN: 0-486-61272-4,
3198 // LC: QA47.A34.
3199 //
3200 // Vladimir Krylov,
3201 // Approximate Calculation of Integrals,
3202 // Dover, 2006,
3203 // ISBN: 0486445798,
3204 // LC: QA311.K713.
3205 //
3206 // Arthur Stroud, Don Secrest,
3207 // Gaussian Quadrature Formulas,
3208 // Prentice Hall, 1966,
3209 // LC: QA299.4G3S7.
3210 //
3211 // Stephen Wolfram,
3212 // The Mathematica Book,
3213 // Fourth Edition,
3214 // Cambridge University Press, 1999,
3215 // ISBN: 0-521-64314-7,
3216 // LC: QA76.95.W65.
3217 //
3218 // Daniel Zwillinger, editor,
3219 // CRC Standard Mathematical Tables and Formulae,
3220 // 30th Edition,
3221 // CRC Press, 1996,
3222 // ISBN: 0-8493-2479-3.
3223 //
3224 // Parameters:
3225 //
3226 // Input, int N, the order.
3227 // N must be between 1 and 20.
3228 //
3229 // Output, Scalar W[N], the weights.
3230 //
3231 {
3232  if (n==1) {
3233  w[ 0] = 1.00000000000000000000000000000E+00;
3234  }
3235  else if (n==2) {
3236  w[ 0] = 0.85355339059327376220042218105E+00;
3237  w[ 1] = 0.146446609406726237799577818948E+00;
3238  }
3239  else if (n==3) {
3240  w[ 0] = 0.71109300992917301544959019114E+00;
3241  w[ 1] = 0.27851773356924084880144488846E+00;
3242  w[ 2] = 0.010389256501586135748964920401E+00;
3243  }
3244  else if (n==4) {
3245  w[ 0] = 0.60315410434163360163596602382E+00;
3246  w[ 1] = 0.35741869243779968664149201746E+00;
3247  w[ 2] = 0.03888790851500538427243816816E+00;
3248  w[ 3] = 0.0005392947055613274501037905676E+00;
3249  }
3250  else if (n==5) {
3251  w[ 0] = 0.52175561058280865247586092879E+00;
3252  w[ 1] = 0.3986668110831759274541333481E+00;
3253  w[ 2] = 0.0759424496817075953876533114E+00;
3254  w[ 3] = 0.00361175867992204845446126257E+00;
3255  w[ 4] = 0.00002336997238577622789114908455E+00;
3256  }
3257  else if (n==6) {
3258  w[ 0] = 0.45896467394996359356828487771E+00;
3259  w[ 1] = 0.4170008307721209941133775662E+00;
3260  w[ 2] = 0.1133733820740449757387061851E+00;
3261  w[ 3] = 0.01039919745314907489891330285E+00;
3262  w[ 4] = 0.000261017202814932059479242860E+00;
3263  w[ 5] = 8.98547906429621238825292053E-07;
3264  }
3265  else if (n==7) {
3266  w[ 0] = 0.40931895170127390213043288002E+00;
3267  w[ 1] = 0.4218312778617197799292810054E+00;
3268  w[ 2] = 0.1471263486575052783953741846E+00;
3269  w[ 3] = 0.0206335144687169398657056150E+00;
3270  w[ 4] = 0.00107401014328074552213195963E+00;
3271  w[ 5] = 0.0000158654643485642012687326223E+00;
3272  w[ 6] = 3.17031547899558056227132215E-08;
3273  }
3274  else if (n==8) {
3275  w[ 0] = 0.36918858934163752992058283938E+00;
3276  w[ 1] = 0.4187867808143429560769785813E+00;
3277  w[ 2] = 0.175794986637171805699659867E+00;
3278  w[ 3] = 0.033343492261215651522132535E+00;
3279  w[ 4] = 0.0027945362352256725249389241E+00;
3280  w[ 5] = 0.00009076508773358213104238501E+00;
3281  w[ 6] = 8.4857467162725315448680183E-07;
3282  w[ 7] = 1.04800117487151038161508854E-09;
3283  }
3284  else if (n==9) {
3285  w[ 0] = 0.336126421797962519673467717606E+00;
3286  w[ 1] = 0.411213980423984387309146942793E+00;
3287  w[ 2] = 0.199287525370885580860575607212E+00;
3288  w[ 3] = 0.0474605627656515992621163600479E+00;
3289  w[ 4] = 0.00559962661079458317700419900556E+00;
3290  w[ 5] = 0.000305249767093210566305412824291E+00;
3291  w[ 6] = 6.59212302607535239225572284875E-06;
3292  w[ 7] = 4.1107693303495484429024104033E-08;
3293  w[ 8] = 3.29087403035070757646681380323E-11;
3294  }
3295  else if (n==10) {
3296  w[ 0] = 0.30844111576502014154747083468E+00;
3297  w[ 1] = 0.4011199291552735515157803099E+00;
3298  w[ 2] = 0.218068287611809421588648523E+00;
3299  w[ 3] = 0.062087456098677747392902129E+00;
3300  w[ 4] = 0.009501516975181100553839072E+00;
3301  w[ 5] = 0.0007530083885875387754559644E+00;
3302  w[ 6] = 0.00002825923349599565567422564E+00;
3303  w[ 7] = 4.249313984962686372586577E-07;
3304  w[ 8] = 1.839564823979630780921535E-09;
3305  w[ 9] = 9.911827219609008558377547E-13;
3306  }
3307  else if (n==11) {
3308  w[ 0] = 0.28493321289420060505605102472E+00;
3309  w[ 1] = 0.3897208895278493779375535080E+00;
3310  w[ 2] = 0.232781831848991333940223796E+00;
3311  w[ 3] = 0.076564453546196686400854179E+00;
3312  w[ 4] = 0.014393282767350695091863919E+00;
3313  w[ 5] = 0.001518880846484873069847776E+00;
3314  w[ 6] = 0.0000851312243547192259720424E+00;
3315  w[ 7] = 2.29240387957450407857683E-06;
3316  w[ 8] = 2.48635370276779587373391E-08;
3317  w[ 9] = 7.71262693369132047028153E-11;
3318  w[10] = 2.883775868323623861597778E-14;
3319  }
3320  else if (n==12) {
3321  w[ 0] = 0.26473137105544319034973889206E+00;
3322  w[ 1] = 0.3777592758731379820244905567E+00;
3323  w[ 2] = 0.244082011319877564254870818E+00;
3324  w[ 3] = 0.09044922221168093072750549E+00;
3325  w[ 4] = 0.02010238115463409652266129E+00;
3326  w[ 5] = 0.002663973541865315881054158E+00;
3327  w[ 6] = 0.000203231592662999392121433E+00;
3328  w[ 7] = 8.3650558568197987453363E-06;
3329  w[ 8] = 1.66849387654091026116990E-07;
3330  w[ 9] = 1.34239103051500414552392E-09;
3331  w[10] = 3.06160163503502078142408E-12;
3332  w[11] = 8.148077467426241682473119E-16;
3333  }
3334  else if (n==13) {
3335  w[ 0] = 0.24718870842996262134624918596E+00;
3336  w[ 1] = 0.3656888229005219453067175309E+00;
3337  w[ 2] = 0.252562420057658502356824289E+00;
3338  w[ 3] = 0.10347075802418370511421863E+00;
3339  w[ 4] = 0.02643275441556161577815877E+00;
3340  w[ 5] = 0.00422039604025475276555209E+00;
3341  w[ 6] = 0.000411881770472734774892473E+00;
3342  w[ 7] = 0.0000235154739815532386882897E+00;
3343  w[ 8] = 7.3173116202490991040105E-07;
3344  w[ 9] = 1.10884162570398067979151E-08;
3345  w[10] = 6.7708266922058988406462E-11;
3346  w[11] = 1.15997995990507606094507E-13;
3347  w[12] = 2.245093203892758415991872E-17;
3348  }
3349  else if (n==14) {
3350  w[ 0] = 0.23181557714486497784077486110E+00;
3351  w[ 1] = 0.3537846915975431518023313013E+00;
3352  w[ 2] = 0.258734610245428085987320561E+00;
3353  w[ 3] = 0.11548289355692321008730499E+00;
3354  w[ 4] = 0.03319209215933736003874996E+00;
3355  w[ 5] = 0.00619286943700661021678786E+00;
3356  w[ 6] = 0.00073989037786738594242589E+00;
3357  w[ 7] = 0.000054907194668416983785733E+00;
3358  w[ 8] = 2.4095857640853774967578E-06;
3359  w[ 9] = 5.801543981676495180886E-08;
3360  w[10] = 6.819314692484974119616E-10;
3361  w[11] = 3.2212077518948479398089E-12;
3362  w[12] = 4.2213524405165873515980E-15;
3363  w[13] = 6.05237502228918880839871E-19;
3364  }
3365  else if (n==15) {
3366  w[ 0] = 0.21823488594008688985641323645E+00;
3367  w[ 1] = 0.3422101779228833296389489568E+00;
3368  w[ 2] = 0.263027577941680097414812275E+00;
3369  w[ 3] = 0.12642581810593053584303055E+00;
3370  w[ 4] = 0.04020686492100091484158548E+00;
3371  w[ 5] = 0.00856387780361183836391576E+00;
3372  w[ 6] = 0.00121243614721425207621921E+00;
3373  w[ 7] = 0.00011167439234425194199258E+00;
3374  w[ 8] = 6.459926762022900924653E-06;
3375  w[ 9] = 2.226316907096272630332E-07;
3376  w[10] = 4.227430384979365007351E-09;
3377  w[11] = 3.921897267041089290385E-11;
3378  w[12] = 1.4565152640731264063327E-13;
3379  w[13] = 1.4830270511133013354616E-16;
3380  w[14] = 1.60059490621113323104998E-20;
3381  }
3382  else if (n==16) {
3383  w[ 0] = 0.20615171495780099433427363674E+00;
3384  w[ 1] = 0.3310578549508841659929830987E+00;
3385  w[ 2] = 0.265795777644214152599502021E+00;
3386  w[ 3] = 0.13629693429637753997554751E+00;
3387  w[ 4] = 0.0473289286941252189780623E+00;
3388  w[ 5] = 0.0112999000803394532312490E+00;
3389  w[ 6] = 0.0018490709435263108642918E+00;
3390  w[ 7] = 0.00020427191530827846012602E+00;
3391  w[ 8] = 0.00001484458687398129877135E+00;
3392  w[ 9] = 6.828319330871199564396E-07;
3393  w[10] = 1.881024841079673213882E-08;
3394  w[11] = 2.862350242973881619631E-10;
3395  w[12] = 2.127079033224102967390E-12;
3396  w[13] = 6.297967002517867787174E-15;
3397  w[14] = 5.050473700035512820402E-18;
3398  w[15] = 4.1614623703728551904265E-22;
3399  }
3400  else if (n==17) {
3401  w[ 0] = 0.19533220525177083214592729770E+00;
3402  w[ 1] = 0.3203753572745402813366256320E+00;
3403  w[ 2] = 0.267329726357171097238809604E+00;
3404  w[ 3] = 0.14512985435875862540742645E+00;
3405  w[ 4] = 0.0544369432453384577793806E+00;
3406  w[ 5] = 0.0143572977660618672917767E+00;
3407  w[ 6] = 0.0026628247355727725684324E+00;
3408  w[ 7] = 0.0003436797271562999206118E+00;
3409  w[ 8] = 0.00003027551783782870109437E+00;
3410  w[ 9] = 1.768515053231676895381E-06;
3411  w[10] = 6.57627288681043332199E-08;
3412  w[11] = 1.469730932159546790344E-09;
3413  w[12] = 1.81691036255544979555E-11;
3414  w[13] = 1.095401388928687402976E-13;
3415  w[14] = 2.617373882223370421551E-16;
3416  w[15] = 1.6729356931461546908502E-19;
3417  w[16] = 1.06562631627404278815253E-23;
3418  }
3419  else if (n==18) {
3420  w[ 0] = 0.18558860314691880562333775228E+00;
3421  w[ 1] = 0.3101817663702252936495975957E+00;
3422  w[ 2] = 0.267866567148536354820854395E+00;
3423  w[ 3] = 0.15297974746807490655384308E+00;
3424  w[ 4] = 0.0614349178609616527076780E+00;
3425  w[ 5] = 0.0176872130807729312772600E+00;
3426  w[ 6] = 0.0036601797677599177980266E+00;
3427  w[ 7] = 0.0005406227870077353231284E+00;
3428  w[ 8] = 0.0000561696505121423113818E+00;
3429  w[ 9] = 4.01530788370115755859E-06;
3430  w[10] = 1.91466985667567497969E-07;
3431  w[11] = 5.8360952686315941292E-09;
3432  w[12] = 1.07171126695539012773E-10;
3433  w[13] = 1.08909871388883385562E-12;
3434  w[14] = 5.38666474837830887608E-15;
3435  w[15] = 1.049865978035703408779E-17;
3436  w[16] = 5.405398451631053643566E-21;
3437  w[17] = 2.6916532692010286270838E-25;
3438  }
3439  else if (n==19) {
3440  w[ 0] = 0.17676847491591250225103547981E+00;
3441  w[ 1] = 0.3004781436072543794821568077E+00;
3442  w[ 2] = 0.267599547038175030772695441E+00;
3443  w[ 3] = 0.15991337213558021678551215E+00;
3444  w[ 4] = 0.0682493799761491134552355E+00;
3445  w[ 5] = 0.0212393076065443249244062E+00;
3446  w[ 6] = 0.0048416273511483959672501E+00;
3447  w[ 7] = 0.0008049127473813667665946E+00;
3448  w[ 8] = 0.0000965247209315350170843E+00;
3449  w[ 9] = 8.20730525805103054409E-06;
3450  w[10] = 4.8305667247307725394E-07;
3451  w[11] = 1.90499136112328569994E-08;
3452  w[12] = 4.8166846309280615577E-10;
3453  w[13] = 7.3482588395511443768E-12;
3454  w[14] = 6.2022753875726163989E-14;
3455  w[15] = 2.54143084301542272372E-16;
3456  w[16] = 4.07886129682571235007E-19;
3457  w[17] = 1.707750187593837061004E-22;
3458  w[18] = 6.715064649908189959990E-27;
3459  }
3460  else if (n==20) {
3461  w[ 0] = 0.168746801851113862149223899689E+00;
3462  w[ 1] = 0.291254362006068281716795323812E+00;
3463  w[ 2] = 0.266686102867001288549520868998E+00;
3464  w[ 3] = 0.166002453269506840031469127816E+00;
3465  w[ 4] = 0.0748260646687923705400624639615E+00;
3466  w[ 5] = 0.0249644173092832210728227383234E+00;
3467  w[ 6] = 0.00620255084457223684744754785395E+00;
3468  w[ 7] = 0.00114496238647690824203955356969E+00;
3469  w[ 8] = 0.000155741773027811974779809513214E+00;
3470  w[ 9] = 0.0000154014408652249156893806714048E+00;
3471  w[10] = 1.08648636651798235147970004439E-06;
3472  w[11] = 5.33012090955671475092780244305E-08;
3473  w[12] = 1.7579811790505820035778763784E-09;
3474  w[13] = 3.72550240251232087262924585338E-11;
3475  w[14] = 4.76752925157819052449488071613E-13;
3476  w[15] = 3.37284424336243841236506064991E-15;
3477  w[16] = 1.15501433950039883096396247181E-17;
3478  w[17] = 1.53952214058234355346383319667E-20;
3479  w[18] = 5.28644272556915782880273587683E-24;
3480  w[19] = 1.65645661249902329590781908529E-28;
3481  }
3482  else {
3483  std::cerr << "\n";
3484  std::cerr << "LAGUERRE_LOOKUP_WEIGHTS - Fatal error!\n";
3485  std::cerr << " Illegal value of N = " << n << "\n";
3486  std::cerr << " Legal values are 1 through 20.\n";
3487  std::exit(1);
3488  }
3489 
3490  return;
3491 }
3492 
3493 //****************************************************************************
3494 template<class Scalar>
3495 void IntrepidBurkardtRules::legendre_compute ( int n, Scalar x[], Scalar w[] )
3496 //****************************************************************************
3497 //
3498 // Purpose:
3499 //
3500 // LEGENDRE_COMPUTE: Legendre quadrature rule by the Elhay-Kautsky method.
3501 //
3502 // Licensing:
3503 //
3504 // This code is distributed under the GNU LGPL license.
3505 //
3506 // Modified:
3507 //
3508 // 19 April 2011
3509 //
3510 // Author:
3511 //
3512 // Original FORTRAN77 version by Sylvan Elhay, Jaroslav Kautsky.
3513 // C++ version by John Burkardt.
3514 //
3515 // Reference:
3516 //
3517 // Sylvan Elhay, Jaroslav Kautsky,
3518 // Algorithm 655: IQPACK, FORTRAN Subroutines for the Weights of
3519 // Interpolatory Quadrature,
3520 // ACM Transactions on Mathematical Software,
3521 // Volume 13, Number 4, December 1987, pages 399-415.
3522 //
3523 // Parameters:
3524 //
3525 // Input, int N, the order.
3526 //
3527 // Output, Scalar X[N], the abscissas.
3528 //
3529 // Output, Scalar W[N], the weights.
3530 //
3531 {
3532  Scalar *bj;
3533  int i;
3534  Scalar zemu;
3535 //
3536 // Define the zero-th moment.
3537 //
3538  zemu = 2.0;
3539 //
3540 // Define the Jacobi matrix.
3541 //
3542  bj = new Scalar[n];
3543 
3544  for (i=0;i<n;i++) {
3545  bj[i] = (Scalar)((i+1)*(i+1))/(Scalar)(4*(i+1)*(i+1)-1);
3546  bj[i] = std::sqrt(bj[i]);
3547  }
3548 
3549  for (i=0;i<n;i++) {
3550  x[i] = 0.0;
3551  }
3552 
3553  w[0] = std::sqrt(zemu);
3554 
3555  for (i=1;i<n;i++) {
3556  w[i] = 0.0;
3557  }
3558 //
3559 // Diagonalize the Jacobi matrix.
3560 //
3561  IntrepidBurkardtRules::imtqlx(n,x,bj,w);
3562 
3563  for (i=0;i<n;i++) {
3564  w[i] = w[i]*w[i];
3565  }
3566 
3567  // Ensure that zero is actually zero.
3568  if (n%2) {
3569  int ind = (int)((Scalar)n/2.0);
3570  x[ind] = 0.0;
3571  }
3572 
3573  delete [] bj;
3574 
3575  return;
3576 }
3577 
3578 //****************************************************************************
3579 template<class Scalar>
3581 //****************************************************************************
3582 //
3583 // Purpose:
3584 //
3585 // LEGENDRE_COMPUTE_POINTS computes Legendre quadrature points.
3586 //
3587 // Licensing:
3588 //
3589 // This code is distributed under the GNU LGPL license.
3590 //
3591 // Modified:
3592 //
3593 // 13 June 2009
3594 //
3595 // Author:
3596 //
3597 // John Burkardt
3598 //
3599 // Parameters:
3600 //
3601 // Input, int N, the order.
3602 //
3603 // Output, Scalar X[N], the abscissas.
3604 //
3605 {
3606  Scalar *w; w= new Scalar[n];
3608  delete [] w;
3609 
3610  return;
3611 }
3612 
3613 //****************************************************************************
3614 template<class Scalar>
3616 //****************************************************************************
3617 //
3618 // Purpose:
3619 //
3620 // LEGENDRE_COMPUTE_WEIGHTS computes Legendre quadrature weights.
3621 //
3622 // Licensing:
3623 //
3624 // This code is distributed under the GNU LGPL license.
3625 //
3626 // Modified:
3627 //
3628 // 13 June 2009
3629 //
3630 // Author:
3631 //
3632 // John Burkardt
3633 //
3634 // Parameters:
3635 //
3636 // Input, int N, the order.
3637 //
3638 // Output, Scalar W[N], the weights.
3639 //
3640 {
3641  Scalar *x; x = new Scalar[n];
3643  delete [] x;
3644 
3645  return;
3646 }
3647 
3648 //****************************************************************************
3649 template<class Scalar>
3650 void IntrepidBurkardtRules::legendre_lookup ( int n, Scalar x[], Scalar w[] )
3651 //****************************************************************************
3652 //
3653 // Purpose:
3654 //
3655 // LEGENDRE_LOOKUP looks up abscissas and weights for Gauss-Legendre quadrature.
3656 //
3657 // Licensing:
3658 //
3659 // This code is distributed under the GNU LGPL license.
3660 //
3661 // Modified:
3662 //
3663 // 27 April 2010
3664 //
3665 // Author:
3666 //
3667 // John Burkardt
3668 //
3669 // Reference:
3670 //
3671 // Milton Abramowitz, Irene Stegun,
3672 // Handbook of Mathematical Functions,
3673 // National Bureau of Standards, 1964,
3674 // ISBN: 0-486-61272-4,
3675 // LC: QA47.A34.
3676 //
3677 // Vladimir Krylov,
3678 // Approximate Calculation of Integrals,
3679 // Dover, 2006,
3680 // ISBN: 0486445798.
3681 // LC: QA311.K713.
3682 //
3683 // Arthur Stroud, Don Secrest,
3684 // Gaussian Quadrature Formulas,
3685 // Prentice Hall, 1966,
3686 // LC: QA299.4G3S7.
3687 //
3688 // Stephen Wolfram,
3689 // The Mathematica Book,
3690 // Fourth Edition,
3691 // Cambridge University Press, 1999,
3692 // ISBN: 0-521-64314-7,
3693 // LC: QA76.95.W65.
3694 //
3695 // Daniel Zwillinger, editor,
3696 // CRC Standard Mathematical Tables and Formulae,
3697 // 30th Edition,
3698 // CRC Press, 1996,
3699 // ISBN: 0-8493-2479-3,
3700 // LC: QA47.M315.
3701 //
3702 // Parameters:
3703 //
3704 // Input, int N, the order.
3705 // N must be between 1 and 33.
3706 //
3707 // Output, Scalar X[N], the abscissas.
3708 //
3709 // Output, Scalar W[N], the abscissas.
3710 //
3711 {
3714 
3715  return;
3716 }
3717 
3718 //****************************************************************************
3719 template<class Scalar>
3721 //****************************************************************************
3722 //
3723 // Purpose:
3724 //
3725 // LEGENDRE_LOOKUP_POINTS looks up abscissas for Gauss-Legendre quadrature.
3726 //
3727 // Licensing:
3728 //
3729 // This code is distributed under the GNU LGPL license.
3730 //
3731 // Modified:
3732 //
3733 // 27 April 2010
3734 //
3735 // Author:
3736 //
3737 // John Burkardt
3738 //
3739 // Reference:
3740 //
3741 // Milton Abramowitz, Irene Stegun,
3742 // Handbook of Mathematical Functions,
3743 // National Bureau of Standards, 1964,
3744 // ISBN: 0-486-61272-4,
3745 // LC: QA47.A34.
3746 //
3747 // Vladimir Krylov,
3748 // Approximate Calculation of Integrals,
3749 // Dover, 2006,
3750 // ISBN: 0486445798.
3751 // LC: QA311.K713.
3752 //
3753 // Arthur Stroud, Don Secrest,
3754 // Gaussian Quadrature Formulas,
3755 // Prentice Hall, 1966,
3756 // LC: QA299.4G3S7.
3757 //
3758 // Stephen Wolfram,
3759 // The Mathematica Book,
3760 // Fourth Edition,
3761 // Cambridge University Press, 1999,
3762 // ISBN: 0-521-64314-7,
3763 // LC: QA76.95.W65.
3764 //
3765 // Daniel Zwillinger, editor,
3766 // CRC Standard Mathematical Tables and Formulae,
3767 // 30th Edition,
3768 // CRC Press, 1996,
3769 // ISBN: 0-8493-2479-3,
3770 // LC: QA47.M315.
3771 //
3772 // Parameters:
3773 //
3774 // Input, int N, the order.
3775 // N must be between 1 and 33.
3776 //
3777 // Output, Scalar X[N], the abscissas.
3778 //
3779 {
3780  if (n==1) {
3781  x[ 0] = 0.000000000000000000000000000000;
3782  }
3783  else if (n==2) {
3784  x[ 0] = -0.577350269189625764509148780502;
3785  x[ 1] = 0.577350269189625764509148780502;
3786  }
3787  else if (n==3) {
3788  x[ 0] = -0.774596669241483377035853079956;
3789  x[ 1] = 0.000000000000000000000000000000;
3790  x[ 2] = 0.774596669241483377035853079956;
3791  }
3792  else if (n==4) {
3793  x[ 0] = -0.861136311594052575223946488893;
3794  x[ 1] = -0.339981043584856264802665759103;
3795  x[ 2] = 0.339981043584856264802665759103;
3796  x[ 3] = 0.861136311594052575223946488893;
3797  }
3798  else if (n==5) {
3799  x[ 0] = -0.906179845938663992797626878299;
3800  x[ 1] = -0.538469310105683091036314420700;
3801  x[ 2] = 0.000000000000000000000000000000;
3802  x[ 3] = 0.538469310105683091036314420700;
3803  x[ 4] = 0.906179845938663992797626878299;
3804  }
3805  else if (n==6) {
3806  x[ 0] = -0.932469514203152027812301554494;
3807  x[ 1] = -0.661209386466264513661399595020;
3808  x[ 2] = -0.238619186083196908630501721681;
3809  x[ 3] = 0.238619186083196908630501721681;
3810  x[ 4] = 0.661209386466264513661399595020;
3811  x[ 5] = 0.932469514203152027812301554494;
3812  }
3813  else if (n==7) {
3814  x[ 0] = -0.949107912342758524526189684048;
3815  x[ 1] = -0.741531185599394439863864773281;
3816  x[ 2] = -0.405845151377397166906606412077;
3817  x[ 3] = 0.000000000000000000000000000000;
3818  x[ 4] = 0.405845151377397166906606412077;
3819  x[ 5] = 0.741531185599394439863864773281;
3820  x[ 6] = 0.949107912342758524526189684048;
3821  }
3822  else if (n==8) {
3823  x[ 0] = -0.960289856497536231683560868569;
3824  x[ 1] = -0.796666477413626739591553936476;
3825  x[ 2] = -0.525532409916328985817739049189;
3826  x[ 3] = -0.183434642495649804939476142360;
3827  x[ 4] = 0.183434642495649804939476142360;
3828  x[ 5] = 0.525532409916328985817739049189;
3829  x[ 6] = 0.796666477413626739591553936476;
3830  x[ 7] = 0.960289856497536231683560868569;
3831  }
3832  else if (n==9) {
3833  x[ 0] = -0.968160239507626089835576203;
3834  x[ 1] = -0.836031107326635794299429788;
3835  x[ 2] = -0.613371432700590397308702039;
3836  x[ 3] = -0.324253423403808929038538015;
3837  x[ 4] = 0.000000000000000000000000000;
3838  x[ 5] = 0.324253423403808929038538015;
3839  x[ 6] = 0.613371432700590397308702039;
3840  x[ 7] = 0.836031107326635794299429788;
3841  x[ 8] = 0.968160239507626089835576203;
3842  }
3843  else if (n==10) {
3844  x[ 0] = -0.973906528517171720077964012;
3845  x[ 1] = -0.865063366688984510732096688;
3846  x[ 2] = -0.679409568299024406234327365;
3847  x[ 3] = -0.433395394129247190799265943;
3848  x[ 4] = -0.148874338981631210884826001;
3849  x[ 5] = 0.148874338981631210884826001;
3850  x[ 6] = 0.433395394129247190799265943;
3851  x[ 7] = 0.679409568299024406234327365;
3852  x[ 8] = 0.865063366688984510732096688;
3853  x[ 9] = 0.973906528517171720077964012;
3854  }
3855  else if (n==11) {
3856  x[ 0] = -0.978228658146056992803938001;
3857  x[ 1] = -0.887062599768095299075157769;
3858  x[ 2] = -0.730152005574049324093416252;
3859  x[ 3] = -0.519096129206811815925725669;
3860  x[ 4] = -0.269543155952344972331531985;
3861  x[ 5] = 0.000000000000000000000000000;
3862  x[ 6] = 0.269543155952344972331531985;
3863  x[ 7] = 0.519096129206811815925725669;
3864  x[ 8] = 0.730152005574049324093416252;
3865  x[ 9] = 0.887062599768095299075157769;
3866  x[10] = 0.978228658146056992803938001;
3867  }
3868  else if (n==12) {
3869  x[ 0] = -0.981560634246719250690549090;
3870  x[ 1] = -0.904117256370474856678465866;
3871  x[ 2] = -0.769902674194304687036893833;
3872  x[ 3] = -0.587317954286617447296702419;
3873  x[ 4] = -0.367831498998180193752691537;
3874  x[ 5] = -0.125233408511468915472441369;
3875  x[ 6] = 0.125233408511468915472441369;
3876  x[ 7] = 0.367831498998180193752691537;
3877  x[ 8] = 0.587317954286617447296702419;
3878  x[ 9] = 0.769902674194304687036893833;
3879  x[10] = 0.904117256370474856678465866;
3880  x[11] = 0.981560634246719250690549090;
3881  }
3882  else if (n==13) {
3883  x[ 0] = -0.984183054718588149472829449;
3884  x[ 1] = -0.917598399222977965206547837;
3885  x[ 2] = -0.801578090733309912794206490;
3886  x[ 3] = -0.642349339440340220643984607;
3887  x[ 4] = -0.448492751036446852877912852;
3888  x[ 5] = -0.230458315955134794065528121;
3889  x[ 6] = 0.000000000000000000000000000;
3890  x[ 7] = 0.230458315955134794065528121;
3891  x[ 8] = 0.448492751036446852877912852;
3892  x[ 9] = 0.642349339440340220643984607;
3893  x[10] = 0.80157809073330991279420649;
3894  x[11] = 0.91759839922297796520654784;
3895  x[12] = 0.98418305471858814947282945;
3896  }
3897  else if (n==14) {
3898  x[ 0] = -0.986283808696812338841597267;
3899  x[ 1] = -0.928434883663573517336391139;
3900  x[ 2] = -0.827201315069764993189794743;
3901  x[ 3] = -0.687292904811685470148019803;
3902  x[ 4] = -0.515248636358154091965290719;
3903  x[ 5] = -0.319112368927889760435671824;
3904  x[ 6] = -0.108054948707343662066244650;
3905  x[ 7] = 0.108054948707343662066244650;
3906  x[ 8] = 0.31911236892788976043567182;
3907  x[ 9] = 0.51524863635815409196529072;
3908  x[10] = 0.68729290481168547014801980;
3909  x[11] = 0.82720131506976499318979474;
3910  x[12] = 0.92843488366357351733639114;
3911  x[13] = 0.98628380869681233884159727;
3912  }
3913  else if (n==15) {
3914  x[ 0] = -0.987992518020485428489565719;
3915  x[ 1] = -0.937273392400705904307758948;
3916  x[ 2] = -0.848206583410427216200648321;
3917  x[ 3] = -0.724417731360170047416186055;
3918  x[ 4] = -0.570972172608538847537226737;
3919  x[ 5] = -0.394151347077563369897207371;
3920  x[ 6] = -0.201194093997434522300628303;
3921  x[ 7] = 0.00000000000000000000000000;
3922  x[ 8] = 0.20119409399743452230062830;
3923  x[ 9] = 0.39415134707756336989720737;
3924  x[10] = 0.57097217260853884753722674;
3925  x[11] = 0.72441773136017004741618605;
3926  x[12] = 0.84820658341042721620064832;
3927  x[13] = 0.93727339240070590430775895;
3928  x[14] = 0.98799251802048542848956572;
3929  }
3930  else if (n==16) {
3931  x[ 0] = -0.989400934991649932596154173;
3932  x[ 1] = -0.944575023073232576077988416;
3933  x[ 2] = -0.865631202387831743880467898;
3934  x[ 3] = -0.755404408355003033895101195;
3935  x[ 4] = -0.617876244402643748446671764;
3936  x[ 5] = -0.458016777657227386342419443;
3937  x[ 6] = -0.281603550779258913230460501;
3938  x[ 7] = -0.09501250983763744018531934;
3939  x[ 8] = 0.09501250983763744018531934;
3940  x[ 9] = 0.28160355077925891323046050;
3941  x[10] = 0.45801677765722738634241944;
3942  x[11] = 0.61787624440264374844667176;
3943  x[12] = 0.75540440835500303389510119;
3944  x[13] = 0.86563120238783174388046790;
3945  x[14] = 0.94457502307323257607798842;
3946  x[15] = 0.98940093499164993259615417;
3947  }
3948  else if (n==17) {
3949  x[ 0] = -0.990575475314417335675434020;
3950  x[ 1] = -0.950675521768767761222716958;
3951  x[ 2] = -0.880239153726985902122955694;
3952  x[ 3] = -0.781514003896801406925230056;
3953  x[ 4] = -0.657671159216690765850302217;
3954  x[ 5] = -0.512690537086476967886246569;
3955  x[ 6] = -0.35123176345387631529718552;
3956  x[ 7] = -0.17848418149584785585067749;
3957  x[ 8] = 0.00000000000000000000000000;
3958  x[ 9] = 0.17848418149584785585067749;
3959  x[10] = 0.35123176345387631529718552;
3960  x[11] = 0.51269053708647696788624657;
3961  x[12] = 0.65767115921669076585030222;
3962  x[13] = 0.78151400389680140692523006;
3963  x[14] = 0.88023915372698590212295569;
3964  x[15] = 0.95067552176876776122271696;
3965  x[16] = 0.99057547531441733567543402;
3966  }
3967  else if (n==18) {
3968  x[ 0] = -0.991565168420930946730016005;
3969  x[ 1] = -0.955823949571397755181195893;
3970  x[ 2] = -0.892602466497555739206060591;
3971  x[ 3] = -0.803704958972523115682417455;
3972  x[ 4] = -0.691687043060353207874891081;
3973  x[ 5] = -0.55977083107394753460787155;
3974  x[ 6] = -0.41175116146284264603593179;
3975  x[ 7] = -0.25188622569150550958897285;
3976  x[ 8] = -0.08477501304173530124226185;
3977  x[ 9] = 0.08477501304173530124226185;
3978  x[10] = 0.25188622569150550958897285;
3979  x[11] = 0.41175116146284264603593179;
3980  x[12] = 0.55977083107394753460787155;
3981  x[13] = 0.69168704306035320787489108;
3982  x[14] = 0.80370495897252311568241746;
3983  x[15] = 0.89260246649755573920606059;
3984  x[16] = 0.95582394957139775518119589;
3985  x[17] = 0.99156516842093094673001600;
3986  }
3987  else if (n==19) {
3988  x[ 0] = -0.992406843843584403189017670;
3989  x[ 1] = -0.960208152134830030852778841;
3990  x[ 2] = -0.903155903614817901642660929;
3991  x[ 3] = -0.822714656537142824978922487;
3992  x[ 4] = -0.72096617733522937861709586;
3993  x[ 5] = -0.60054530466168102346963816;
3994  x[ 6] = -0.46457074137596094571726715;
3995  x[ 7] = -0.31656409996362983199011733;
3996  x[ 8] = -0.16035864564022537586809612;
3997  x[ 9] = 0.00000000000000000000000000;
3998  x[10] = 0.16035864564022537586809612;
3999  x[11] = 0.31656409996362983199011733;
4000  x[12] = 0.46457074137596094571726715;
4001  x[13] = 0.60054530466168102346963816;
4002  x[14] = 0.72096617733522937861709586;
4003  x[15] = 0.82271465653714282497892249;
4004  x[16] = 0.90315590361481790164266093;
4005  x[17] = 0.96020815213483003085277884;
4006  x[18] = 0.99240684384358440318901767;
4007  }
4008  else if (n==20) {
4009  x[ 0] = -0.993128599185094924786122388;
4010  x[ 1] = -0.963971927277913791267666131;
4011  x[ 2] = -0.912234428251325905867752441;
4012  x[ 3] = -0.83911697182221882339452906;
4013  x[ 4] = -0.74633190646015079261430507;
4014  x[ 5] = -0.63605368072651502545283670;
4015  x[ 6] = -0.51086700195082709800436405;
4016  x[ 7] = -0.37370608871541956067254818;
4017  x[ 8] = -0.22778585114164507808049620;
4018  x[ 9] = -0.07652652113349733375464041;
4019  x[10] = 0.07652652113349733375464041;
4020  x[11] = 0.22778585114164507808049620;
4021  x[12] = 0.37370608871541956067254818;
4022  x[13] = 0.51086700195082709800436405;
4023  x[14] = 0.63605368072651502545283670;
4024  x[15] = 0.74633190646015079261430507;
4025  x[16] = 0.83911697182221882339452906;
4026  x[17] = 0.91223442825132590586775244;
4027  x[18] = 0.96397192727791379126766613;
4028  x[19] = 0.99312859918509492478612239;
4029  }
4030  else if (n==21) {
4031  x[ 0] = -0.99375217062038950026024204;
4032  x[ 1] = -0.96722683856630629431662221;
4033  x[ 2] = -0.92009933415040082879018713;
4034  x[ 3] = -0.85336336458331728364725064;
4035  x[ 4] = -0.76843996347567790861587785;
4036  x[ 5] = -0.66713880419741231930596667;
4037  x[ 6] = -0.55161883588721980705901880;
4038  x[ 7] = -0.42434212020743878357366889;
4039  x[ 8] = -0.28802131680240109660079252;
4040  x[ 9] = -0.14556185416089509093703098;
4041  x[10] = 0.00000000000000000000000000;
4042  x[11] = +0.14556185416089509093703098;
4043  x[12] = +0.28802131680240109660079252;
4044  x[13] = +0.42434212020743878357366889;
4045  x[14] = +0.55161883588721980705901880;
4046  x[15] = +0.66713880419741231930596667;
4047  x[16] = +0.76843996347567790861587785;
4048  x[17] = +0.85336336458331728364725064;
4049  x[18] = +0.92009933415040082879018713;
4050  x[19] = +0.96722683856630629431662221;
4051  x[20] = +0.99375217062038950026024204;
4052  }
4053  else if (n==22) {
4054  x[ 0] = -0.99429458548239929207303142;
4055  x[ 1] = -0.97006049783542872712395099;
4056  x[ 2] = -0.92695677218717400052069294;
4057  x[ 3] = -0.86581257772030013653642564;
4058  x[ 4] = -0.78781680597920816200427796;
4059  x[ 5] = -0.69448726318668278005068984;
4060  x[ 6] = -0.58764040350691159295887693;
4061  x[ 7] = -0.46935583798675702640633071;
4062  x[ 8] = -0.34193582089208422515814742;
4063  x[ 9] = -0.20786042668822128547884653;
4064  x[10] = -0.06973927331972222121384180;
4065  x[11] = 0.06973927331972222121384180;
4066  x[12] = 0.20786042668822128547884653;
4067  x[13] = 0.34193582089208422515814742;
4068  x[14] = 0.46935583798675702640633071;
4069  x[15] = 0.58764040350691159295887693;
4070  x[16] = 0.69448726318668278005068984;
4071  x[17] = 0.78781680597920816200427796;
4072  x[18] = 0.86581257772030013653642564;
4073  x[19] = 0.92695677218717400052069294;
4074  x[20] = 0.97006049783542872712395099;
4075  x[21] = 0.99429458548239929207303142;
4076  }
4077  else if (n==23) {
4078  x[ 0] = -0.99476933499755212352392572;
4079  x[ 1] = -0.97254247121811523195602408;
4080  x[ 2] = -0.93297108682601610234919699;
4081  x[ 3] = -0.87675235827044166737815689;
4082  x[ 4] = -0.80488840161883989215111841;
4083  x[ 5] = -0.71866136313195019446162448;
4084  x[ 6] = -0.61960987576364615638509731;
4085  x[ 7] = -0.50950147784600754968979305;
4086  x[ 8] = -0.39030103803029083142148887;
4087  x[ 9] = -0.26413568097034493053386954;
4088  x[10] = -0.13325682429846611093174268;
4089  x[11] = 0.00000000000000000000000000;
4090  x[12] = 0.13325682429846611093174268;
4091  x[13] = 0.26413568097034493053386954;
4092  x[14] = 0.39030103803029083142148887;
4093  x[15] = 0.50950147784600754968979305;
4094  x[16] = 0.61960987576364615638509731;
4095  x[17] = 0.71866136313195019446162448;
4096  x[18] = 0.80488840161883989215111841;
4097  x[19] = 0.87675235827044166737815689;
4098  x[20] = 0.93297108682601610234919699;
4099  x[21] = 0.97254247121811523195602408;
4100  x[22] = 0.99476933499755212352392572;
4101  }
4102  else if (n==24) {
4103  x[ 0] = -0.99518721999702136017999741;
4104  x[ 1] = -0.97472855597130949819839199;
4105  x[ 2] = -0.93827455200273275852364900;
4106  x[ 3] = -0.88641552700440103421315434;
4107  x[ 4] = -0.82000198597390292195394987;
4108  x[ 5] = -0.74012419157855436424382810;
4109  x[ 6] = -0.64809365193697556925249579;
4110  x[ 7] = -0.54542147138883953565837562;
4111  x[ 8] = -0.43379350762604513848708423;
4112  x[ 9] = -0.31504267969616337438679329;
4113  x[10] = -0.19111886747361630915863982;
4114  x[11] = -0.06405689286260562608504308;
4115  x[12] = 0.06405689286260562608504308;
4116  x[13] = 0.19111886747361630915863982;
4117  x[14] = 0.31504267969616337438679329;
4118  x[15] = 0.43379350762604513848708423;
4119  x[16] = 0.54542147138883953565837562;
4120  x[17] = 0.64809365193697556925249579;
4121  x[18] = 0.74012419157855436424382810;
4122  x[19] = 0.82000198597390292195394987;
4123  x[20] = 0.88641552700440103421315434;
4124  x[21] = 0.93827455200273275852364900;
4125  x[22] = 0.97472855597130949819839199;
4126  x[23] = 0.99518721999702136017999741;
4127  }
4128  else if (n==25) {
4129  x[ 0] = -0.99555696979049809790878495;
4130  x[ 1] = -0.97666392145951751149831539;
4131  x[ 2] = -0.94297457122897433941401117;
4132  x[ 3] = -0.89499199787827536885104201;
4133  x[ 4] = -0.83344262876083400142102111;
4134  x[ 5] = -0.75925926303735763057728287;
4135  x[ 6] = -0.67356636847346836448512063;
4136  x[ 7] = -0.57766293024122296772368984;
4137  x[ 8] = -0.47300273144571496052218212;
4138  x[ 9] = -0.36117230580938783773582173;
4139  x[10] = -0.24386688372098843204519036;
4140  x[11] = -0.12286469261071039638735982;
4141  x[12] = 0.00000000000000000000000000;
4142  x[13] = 0.12286469261071039638735982;
4143  x[14] = 0.24386688372098843204519036;
4144  x[15] = 0.36117230580938783773582173;
4145  x[16] = 0.47300273144571496052218212;
4146  x[17] = 0.57766293024122296772368984;
4147  x[18] = 0.67356636847346836448512063;
4148  x[19] = 0.75925926303735763057728287;
4149  x[20] = 0.83344262876083400142102111;
4150  x[21] = 0.89499199787827536885104201;
4151  x[22] = 0.94297457122897433941401117;
4152  x[23] = 0.97666392145951751149831539;
4153  x[24] = 0.99555696979049809790878495;
4154  }
4155  else if (n==26) {
4156  x[ 0] = -0.99588570114561692900321696;
4157  x[ 1] = -0.97838544595647099110058035;
4158  x[ 2] = -0.94715906666171425013591528;
4159  x[ 3] = -0.90263786198430707421766560;
4160  x[ 4] = -0.84544594278849801879750706;
4161  x[ 5] = -0.77638594882067885619296725;
4162  x[ 6] = -0.69642726041995726486381391;
4163  x[ 7] = -0.60669229301761806323197875;
4164  x[ 8] = -0.50844071482450571769570306;
4165  x[ 9] = -0.40305175512348630648107738;
4166  x[10] = -0.29200483948595689514283538;
4167  x[11] = -0.17685882035689018396905775;
4168  x[12] = -0.05923009342931320709371858;
4169  x[13] = 0.05923009342931320709371858;
4170  x[14] = 0.17685882035689018396905775;
4171  x[15] = 0.29200483948595689514283538;
4172  x[16] = 0.40305175512348630648107738;
4173  x[17] = 0.50844071482450571769570306;
4174  x[18] = 0.60669229301761806323197875;
4175  x[19] = 0.69642726041995726486381391;
4176  x[20] = 0.77638594882067885619296725;
4177  x[21] = 0.84544594278849801879750706;
4178  x[22] = 0.90263786198430707421766560;
4179  x[23] = 0.94715906666171425013591528;
4180  x[24] = 0.97838544595647099110058035;
4181  x[25] = 0.99588570114561692900321696;
4182  }
4183  else if (n==27) {
4184  x[ 0] = -0.99617926288898856693888721;
4185  x[ 1] = -0.97992347596150122285587336;
4186  x[ 2] = -0.95090055781470500685190803;
4187  x[ 3] = -0.90948232067749110430064502;
4188  x[ 4] = -0.85620790801829449030273722;
4189  x[ 5] = -0.79177163907050822714439734;
4190  x[ 6] = -0.71701347373942369929481621;
4191  x[ 7] = -0.63290797194649514092773464;
4192  x[ 8] = -0.54055156457945689490030094;
4193  x[ 9] = -0.44114825175002688058597416;
4194  x[10] = -0.33599390363850889973031903;
4195  x[11] = -0.22645936543953685885723911;
4196  x[12] = -0.11397258560952996693289498;
4197  x[13] = 0.00000000000000000000000000;
4198  x[14] = 0.11397258560952996693289498;
4199  x[15] = 0.22645936543953685885723911;
4200  x[16] = 0.33599390363850889973031903;
4201  x[17] = 0.44114825175002688058597416;
4202  x[18] = 0.54055156457945689490030094;
4203  x[19] = 0.63290797194649514092773464;
4204  x[20] = 0.71701347373942369929481621;
4205  x[21] = 0.79177163907050822714439734;
4206  x[22] = 0.85620790801829449030273722;
4207  x[23] = 0.90948232067749110430064502;
4208  x[24] = 0.95090055781470500685190803;
4209  x[25] = 0.97992347596150122285587336;
4210  x[26] = 0.99617926288898856693888721;
4211  }
4212  else if (n==28) {
4213  x[ 0] = -0.99644249757395444995043639;
4214  x[ 1] = -0.98130316537087275369455995;
4215  x[ 2] = -0.95425928062893819725410184;
4216  x[ 3] = -0.91563302639213207386968942;
4217  x[ 4] = -0.86589252257439504894225457;
4218  x[ 5] = -0.80564137091717917144788596;
4219  x[ 6] = -0.73561087801363177202814451;
4220  x[ 7] = -0.65665109403886496121989818;
4221  x[ 8] = -0.56972047181140171930800328;
4222  x[ 9] = -0.47587422495511826103441185;
4223  x[10] = -0.37625151608907871022135721;
4224  x[11] = -0.27206162763517807767682636;
4225  x[12] = -0.16456928213338077128147178;
4226  x[13] = -0.05507928988403427042651653;
4227  x[14] = 0.05507928988403427042651653;
4228  x[15] = 0.16456928213338077128147178;
4229  x[16] = 0.27206162763517807767682636;
4230  x[17] = 0.37625151608907871022135721;
4231  x[18] = 0.47587422495511826103441185;
4232  x[19] = 0.56972047181140171930800328;
4233  x[20] = 0.65665109403886496121989818;
4234  x[21] = 0.73561087801363177202814451;
4235  x[22] = 0.80564137091717917144788596;
4236  x[23] = 0.86589252257439504894225457;
4237  x[24] = 0.91563302639213207386968942;
4238  x[25] = 0.95425928062893819725410184;
4239  x[26] = 0.98130316537087275369455995;
4240  x[27] = 0.99644249757395444995043639;
4241  }
4242  else if (n==29) {
4243  x[ 0] = -0.99667944226059658616319153;
4244  x[ 1] = -0.98254550526141317487092602;
4245  x[ 2] = -0.95728559577808772579820804;
4246  x[ 3] = -0.92118023295305878509375344;
4247  x[ 4] = -0.87463780492010279041779342;
4248  x[ 5] = -0.81818548761525244498957221;
4249  x[ 6] = -0.75246285173447713391261008;
4250  x[ 7] = -0.67821453760268651515618501;
4251  x[ 8] = -0.59628179713822782037958621;
4252  x[ 9] = -0.50759295512422764210262792;
4253  x[10] = -0.41315288817400866389070659;
4254  x[11] = -0.31403163786763993494819592;
4255  x[12] = -0.21135228616600107450637573;
4256  x[13] = -0.10627823013267923017098239;
4257  x[14] = 0.00000000000000000000000000;
4258  x[15] = 0.10627823013267923017098239;
4259  x[16] = 0.21135228616600107450637573;
4260  x[17] = 0.31403163786763993494819592;
4261  x[18] = 0.41315288817400866389070659;
4262  x[19] = 0.50759295512422764210262792;
4263  x[20] = 0.59628179713822782037958621;
4264  x[21] = 0.67821453760268651515618501;
4265  x[22] = 0.75246285173447713391261008;
4266  x[23] = 0.81818548761525244498957221;
4267  x[24] = 0.87463780492010279041779342;
4268  x[25] = 0.92118023295305878509375344;
4269  x[26] = 0.95728559577808772579820804;
4270  x[27] = 0.98254550526141317487092602;
4271  x[28] = 0.99667944226059658616319153;
4272  }
4273  else if (n==30) {
4274  x[ 0] = -0.99689348407464954027163005;
4275  x[ 1] = -0.98366812327974720997003258;
4276  x[ 2] = -0.96002186496830751221687103;
4277  x[ 3] = -0.92620004742927432587932428;
4278  x[ 4] = -0.88256053579205268154311646;
4279  x[ 5] = -0.82956576238276839744289812;
4280  x[ 6] = -0.76777743210482619491797734;
4281  x[ 7] = -0.69785049479331579693229239;
4282  x[ 8] = -0.62052618298924286114047756;
4283  x[ 9] = -0.53662414814201989926416979;
4284  x[10] = -0.44703376953808917678060990;
4285  x[11] = -0.35270472553087811347103721;
4286  x[12] = -0.25463692616788984643980513;
4287  x[13] = -0.15386991360858354696379467;
4288  x[14] = -0.05147184255531769583302521;
4289  x[15] = 0.05147184255531769583302521;
4290  x[16] = 0.15386991360858354696379467;
4291  x[17] = 0.25463692616788984643980513;
4292  x[18] = 0.35270472553087811347103721;
4293  x[19] = 0.44703376953808917678060990;
4294  x[20] = 0.53662414814201989926416979;
4295  x[21] = 0.62052618298924286114047756;
4296  x[22] = 0.69785049479331579693229239;
4297  x[23] = 0.76777743210482619491797734;
4298  x[24] = 0.82956576238276839744289812;
4299  x[25] = 0.88256053579205268154311646;
4300  x[26] = 0.92620004742927432587932428;
4301  x[27] = 0.96002186496830751221687103;
4302  x[28] = 0.98366812327974720997003258;
4303  x[29] = 0.99689348407464954027163005;
4304  }
4305  else if (n==31) {
4306  x[ 0] = -0.99708748181947707405562655;
4307  x[ 1] = -0.98468590966515248400246517;
4308  x[ 2] = -0.96250392509294966178905240;
4309  x[ 3] = -0.93075699789664816495694576;
4310  x[ 4] = -0.88976002994827104337419201;
4311  x[ 5] = -0.83992032014626734008690454;
4312  x[ 6] = -0.78173314841662494040636002;
4313  x[ 7] = -0.71577678458685328390597087;
4314  x[ 8] = -0.64270672292426034618441820;
4315  x[ 9] = -0.56324916140714926272094492;
4316  x[10] = -0.47819378204490248044059404;
4317  x[11] = -0.38838590160823294306135146;
4318  x[12] = -0.29471806998170161661790390;
4319  x[13] = -0.19812119933557062877241300;
4320  x[14] = -0.09955531215234152032517479;
4321  x[15] = 0.00000000000000000000000000;
4322  x[16] = 0.09955531215234152032517479;
4323  x[17] = 0.19812119933557062877241300;
4324  x[18] = 0.29471806998170161661790390;
4325  x[19] = 0.38838590160823294306135146;
4326  x[20] = 0.47819378204490248044059404;
4327  x[21] = 0.56324916140714926272094492;
4328  x[22] = 0.64270672292426034618441820;
4329  x[23] = 0.71577678458685328390597087;
4330  x[24] = 0.78173314841662494040636002;
4331  x[25] = 0.83992032014626734008690454;
4332  x[26] = 0.88976002994827104337419201;
4333  x[27] = 0.93075699789664816495694576;
4334  x[28] = 0.96250392509294966178905240;
4335  x[29] = 0.98468590966515248400246517;
4336  x[30] = 0.99708748181947707405562655;
4337  }
4338  else if (n==32) {
4339  x[ 0] = -0.99726386184948156354498113;
4340  x[ 1] = -0.98561151154526833540017504;
4341  x[ 2] = -0.96476225558750643077381193;
4342  x[ 3] = -0.93490607593773968917091913;
4343  x[ 4] = -0.89632115576605212396530724;
4344  x[ 5] = -0.84936761373256997013369300;
4345  x[ 6] = -0.79448379596794240696309730;
4346  x[ 7] = -0.73218211874028968038742667;
4347  x[ 8] = -0.66304426693021520097511517;
4348  x[ 9] = -0.58771575724076232904074548;
4349  x[10] = -0.50689990893222939002374747;
4350  x[11] = -0.42135127613063534536411944;
4351  x[12] = -0.33186860228212764977991681;
4352  x[13] = -0.23928736225213707454460321;
4353  x[14] = -0.14447196158279649348518637;
4354  x[15] = -0.04830766568773831623481257;
4355  x[16] = 0.04830766568773831623481257;
4356  x[17] = 0.14447196158279649348518637;
4357  x[18] = 0.23928736225213707454460321;
4358  x[19] = 0.33186860228212764977991681;
4359  x[20] = 0.42135127613063534536411944;
4360  x[21] = 0.50689990893222939002374747;
4361  x[22] = 0.58771575724076232904074548;
4362  x[23] = 0.66304426693021520097511517;
4363  x[24] = 0.73218211874028968038742667;
4364  x[25] = 0.79448379596794240696309730;
4365  x[26] = 0.84936761373256997013369300;
4366  x[27] = 0.89632115576605212396530724;
4367  x[28] = 0.93490607593773968917091913;
4368  x[29] = 0.96476225558750643077381193;
4369  x[30] = 0.98561151154526833540017504;
4370  x[31] = 0.99726386184948156354498113;
4371  }
4372  else if (n==33) {
4373  x[ 0] = -0.99742469424645521726616802;
4374  x[ 1] = -0.98645572623064248811037570;
4375  x[ 2] = -0.96682290968999276892837771;
4376  x[ 3] = -0.93869437261116835035583512;
4377  x[ 4] = -0.90231676774343358304053133;
4378  x[ 5] = -0.85800965267650406464306148;
4379  x[ 6] = -0.80616235627416658979620087;
4380  x[ 7] = -0.74723049644956215785905512;
4381  x[ 8] = -0.68173195996974278626821595;
4382  x[ 9] = -0.61024234583637902730728751;
4383  x[10] = -0.53338990478634764354889426;
4384  x[11] = -0.45185001727245069572599328;
4385  x[12] = -0.36633925774807334107022062;
4386  x[13] = -0.27760909715249702940324807;
4387  x[14] = -0.18643929882799157233579876;
4388  x[15] = -0.09363106585473338567074292;
4389  x[16] = 0.00000000000000000000000000;
4390  x[17] = 0.09363106585473338567074292;
4391  x[18] = 0.18643929882799157233579876;
4392  x[19] = 0.27760909715249702940324807;
4393  x[20] = 0.36633925774807334107022062;
4394  x[21] = 0.45185001727245069572599328;
4395  x[22] = 0.53338990478634764354889426;
4396  x[23] = 0.61024234583637902730728751;
4397  x[24] = 0.68173195996974278626821595;
4398  x[25] = 0.74723049644956215785905512;
4399  x[26] = 0.80616235627416658979620087;
4400  x[27] = 0.85800965267650406464306148;
4401  x[28] = 0.90231676774343358304053133;
4402  x[29] = 0.93869437261116835035583512;
4403  x[30] = 0.96682290968999276892837771;
4404  x[31] = 0.98645572623064248811037570;
4405  x[32] = 0.99742469424645521726616802;
4406  }
4407  else {
4408  std::cerr << "\n";
4409  std::cerr << "LEGENDRE_LOOKUP_POINTS - Fatal error!\n";
4410  std::cerr << " Illegal value of N = " << n << "\n";
4411  std::cerr << " Legal values are 1 through 33.\n";
4412  std::exit(1);
4413  }
4414  return;
4415 }
4416 
4417 //****************************************************************************
4418 template<class Scalar>
4420 //****************************************************************************
4421 //
4422 // Purpose:
4423 //
4424 // LEGENDRE_LOOKUP_WEIGHTS looks up weights for Gauss-Legendre quadrature.
4425 //
4426 // Licensing:
4427 //
4428 // This code is distributed under the GNU LGPL license.
4429 //
4430 // Modified:
4431 //
4432 // 27 April 2010
4433 //
4434 // Author:
4435 //
4436 // John Burkardt
4437 //
4438 // Reference:
4439 //
4440 // Milton Abramowitz, Irene Stegun,
4441 // Handbook of Mathematical Functions,
4442 // National Bureau of Standards, 1964,
4443 // ISBN: 0-486-61272-4,
4444 // LC: QA47.A34.
4445 //
4446 // Vladimir Krylov,
4447 // Approximate Calculation of Integrals,
4448 // Dover, 2006,
4449 // ISBN: 0486445798.
4450 // LC: QA311.K713.
4451 //
4452 // Arthur Stroud, Don Secrest,
4453 // Gaussian Quadrature Formulas,
4454 // Prentice Hall, 1966,
4455 // LC: QA299.4G3S7.
4456 //
4457 // Stephen Wolfram,
4458 // The Mathematica Book,
4459 // Fourth Edition,
4460 // Cambridge University Press, 1999,
4461 // ISBN: 0-521-64314-7,
4462 // LC: QA76.95.W65.
4463 //
4464 // Daniel Zwillinger, editor,
4465 // CRC Standard Mathematical Tables and Formulae,
4466 // 30th Edition,
4467 // CRC Press, 1996,
4468 // ISBN: 0-8493-2479-3,
4469 // LC: QA47.M315.
4470 //
4471 // Parameters:
4472 //
4473 // Input, int N, the order.
4474 // N must be between 1 and 33.
4475 //
4476 // Output, Scalar W[N], the weights.
4477 //
4478 {
4479  if (n==1) {
4480  w[ 0] = 2.000000000000000000000000000000;
4481  }
4482  else if (n==2) {
4483  w[ 0] = 1.000000000000000000000000000000;
4484  w[ 1] = 1.000000000000000000000000000000;
4485  }
4486  else if (n==3) {
4487  w[ 0] = 0.555555555555555555555555555556;
4488  w[ 1] = 0.888888888888888888888888888889;
4489  w[ 2] = 0.555555555555555555555555555556;
4490  }
4491  else if (n==4) {
4492  w[ 0] = 0.347854845137453857373063949222;
4493  w[ 1] = 0.652145154862546142626936050778;
4494  w[ 2] = 0.652145154862546142626936050778;
4495  w[ 3] = 0.347854845137453857373063949222;
4496  }
4497  else if (n==5) {
4498  w[ 0] = 0.236926885056189087514264040720;
4499  w[ 1] = 0.478628670499366468041291514836;
4500  w[ 2] = 0.568888888888888888888888888889;
4501  w[ 3] = 0.478628670499366468041291514836;
4502  w[ 4] = 0.236926885056189087514264040720;
4503  }
4504  else if (n==6) {
4505  w[ 0] = 0.171324492379170345040296142173;
4506  w[ 1] = 0.360761573048138607569833513838;
4507  w[ 2] = 0.467913934572691047389870343990;
4508  w[ 3] = 0.467913934572691047389870343990;
4509  w[ 4] = 0.360761573048138607569833513838;
4510  w[ 5] = 0.171324492379170345040296142173;
4511  }
4512  else if (n==7) {
4513  w[ 0] = 0.129484966168869693270611432679;
4514  w[ 1] = 0.279705391489276667901467771424;
4515  w[ 2] = 0.381830050505118944950369775489;
4516  w[ 3] = 0.417959183673469387755102040816;
4517  w[ 4] = 0.381830050505118944950369775489;
4518  w[ 5] = 0.279705391489276667901467771424;
4519  w[ 6] = 0.129484966168869693270611432679;
4520  }
4521  else if (n==8) {
4522  w[ 0] = 0.101228536290376259152531354310;
4523  w[ 1] = 0.222381034453374470544355994426;
4524  w[ 2] = 0.313706645877887287337962201987;
4525  w[ 3] = 0.362683783378361982965150449277;
4526  w[ 4] = 0.362683783378361982965150449277;
4527  w[ 5] = 0.313706645877887287337962201987;
4528  w[ 6] = 0.222381034453374470544355994426;
4529  w[ 7] = 0.101228536290376259152531354310;
4530  }
4531  else if (n==9) {
4532  w[ 0] = 0.081274388361574411971892158111;
4533  w[ 1] = 0.18064816069485740405847203124;
4534  w[ 2] = 0.26061069640293546231874286942;
4535  w[ 3] = 0.31234707704000284006863040658;
4536  w[ 4] = 0.33023935500125976316452506929;
4537  w[ 5] = 0.31234707704000284006863040658;
4538  w[ 6] = 0.26061069640293546231874286942;
4539  w[ 7] = 0.18064816069485740405847203124;
4540  w[ 8] = 0.081274388361574411971892158111;
4541  }
4542  else if (n==10) {
4543  w[ 0] = 0.066671344308688137593568809893;
4544  w[ 1] = 0.14945134915058059314577633966;
4545  w[ 2] = 0.21908636251598204399553493423;
4546  w[ 3] = 0.26926671930999635509122692157;
4547  w[ 4] = 0.29552422471475287017389299465;
4548  w[ 5] = 0.29552422471475287017389299465;
4549  w[ 6] = 0.26926671930999635509122692157;
4550  w[ 7] = 0.21908636251598204399553493423;
4551  w[ 8] = 0.14945134915058059314577633966;
4552  w[ 9] = 0.066671344308688137593568809893;
4553  }
4554  else if (n==11) {
4555  w[ 0] = 0.055668567116173666482753720443;
4556  w[ 1] = 0.12558036946490462463469429922;
4557  w[ 2] = 0.18629021092773425142609764143;
4558  w[ 3] = 0.23319376459199047991852370484;
4559  w[ 4] = 0.26280454451024666218068886989;
4560  w[ 5] = 0.27292508677790063071448352834;
4561  w[ 6] = 0.26280454451024666218068886989;
4562  w[ 7] = 0.23319376459199047991852370484;
4563  w[ 8] = 0.18629021092773425142609764143;
4564  w[ 9] = 0.12558036946490462463469429922;
4565  w[10] = 0.055668567116173666482753720443;
4566  }
4567  else if (n==12) {
4568  w[ 0] = 0.047175336386511827194615961485;
4569  w[ 1] = 0.10693932599531843096025471819;
4570  w[ 2] = 0.16007832854334622633465252954;
4571  w[ 3] = 0.20316742672306592174906445581;
4572  w[ 4] = 0.23349253653835480876084989892;
4573  w[ 5] = 0.24914704581340278500056243604;
4574  w[ 6] = 0.24914704581340278500056243604;
4575  w[ 7] = 0.23349253653835480876084989892;
4576  w[ 8] = 0.20316742672306592174906445581;
4577  w[ 9] = 0.16007832854334622633465252954;
4578  w[10] = 0.10693932599531843096025471819;
4579  w[11] = 0.047175336386511827194615961485;
4580  }
4581  else if (n==13) {
4582  w[ 0] = 0.040484004765315879520021592201;
4583  w[ 1] = 0.092121499837728447914421775954;
4584  w[ 2] = 0.13887351021978723846360177687;
4585  w[ 3] = 0.17814598076194573828004669200;
4586  w[ 4] = 0.20781604753688850231252321931;
4587  w[ 5] = 0.22628318026289723841209018604;
4588  w[ 6] = 0.23255155323087391019458951527;
4589  w[ 7] = 0.22628318026289723841209018604;
4590  w[ 8] = 0.20781604753688850231252321931;
4591  w[ 9] = 0.17814598076194573828004669200;
4592  w[10] = 0.13887351021978723846360177687;
4593  w[11] = 0.092121499837728447914421775954;
4594  w[12] = 0.040484004765315879520021592201;
4595  }
4596  else if (n==14) {
4597  w[ 0] = 0.035119460331751863031832876138;
4598  w[ 1] = 0.08015808715976020980563327706;
4599  w[ 2] = 0.12151857068790318468941480907;
4600  w[ 3] = 0.15720316715819353456960193862;
4601  w[ 4] = 0.18553839747793781374171659013;
4602  w[ 5] = 0.20519846372129560396592406566;
4603  w[ 6] = 0.21526385346315779019587644332;
4604  w[ 7] = 0.21526385346315779019587644332;
4605  w[ 8] = 0.20519846372129560396592406566;
4606  w[ 9] = 0.18553839747793781374171659013;
4607  w[10] = 0.15720316715819353456960193862;
4608  w[11] = 0.12151857068790318468941480907;
4609  w[12] = 0.08015808715976020980563327706;
4610  w[13] = 0.035119460331751863031832876138;
4611  }
4612  else if (n==15) {
4613  w[ 0] = 0.030753241996117268354628393577;
4614  w[ 1] = 0.070366047488108124709267416451;
4615  w[ 2] = 0.107159220467171935011869546686;
4616  w[ 3] = 0.13957067792615431444780479451;
4617  w[ 4] = 0.16626920581699393355320086048;
4618  w[ 5] = 0.18616100001556221102680056187;
4619  w[ 6] = 0.19843148532711157645611832644;
4620  w[ 7] = 0.20257824192556127288062019997;
4621  w[ 8] = 0.19843148532711157645611832644;
4622  w[ 9] = 0.18616100001556221102680056187;
4623  w[10] = 0.16626920581699393355320086048;
4624  w[11] = 0.13957067792615431444780479451;
4625  w[12] = 0.107159220467171935011869546686;
4626  w[13] = 0.070366047488108124709267416451;
4627  w[14] = 0.030753241996117268354628393577;
4628  }
4629  else if (n==16) {
4630  w[ 0] = 0.027152459411754094851780572456;
4631  w[ 1] = 0.062253523938647892862843836994;
4632  w[ 2] = 0.09515851168249278480992510760;
4633  w[ 3] = 0.12462897125553387205247628219;
4634  w[ 4] = 0.14959598881657673208150173055;
4635  w[ 5] = 0.16915651939500253818931207903;
4636  w[ 6] = 0.18260341504492358886676366797;
4637  w[ 7] = 0.18945061045506849628539672321;
4638  w[ 8] = 0.18945061045506849628539672321;
4639  w[ 9] = 0.18260341504492358886676366797;
4640  w[10] = 0.16915651939500253818931207903;
4641  w[11] = 0.14959598881657673208150173055;
4642  w[12] = 0.12462897125553387205247628219;
4643  w[13] = 0.09515851168249278480992510760;
4644  w[14] = 0.062253523938647892862843836994;
4645  w[15] = 0.027152459411754094851780572456;
4646  }
4647  else if (n==17) {
4648  w[ 0] = 0.024148302868547931960110026288;
4649  w[ 1] = 0.055459529373987201129440165359;
4650  w[ 2] = 0.085036148317179180883535370191;
4651  w[ 3] = 0.111883847193403971094788385626;
4652  w[ 4] = 0.13513636846852547328631998170;
4653  w[ 5] = 0.15404576107681028808143159480;
4654  w[ 6] = 0.16800410215645004450997066379;
4655  w[ 7] = 0.17656270536699264632527099011;
4656  w[ 8] = 0.17944647035620652545826564426;
4657  w[ 9] = 0.17656270536699264632527099011;
4658  w[10] = 0.16800410215645004450997066379;
4659  w[11] = 0.15404576107681028808143159480;
4660  w[12] = 0.13513636846852547328631998170;
4661  w[13] = 0.111883847193403971094788385626;
4662  w[14] = 0.085036148317179180883535370191;
4663  w[15] = 0.055459529373987201129440165359;
4664  w[16] = 0.024148302868547931960110026288;
4665  }
4666  else if (n==18) {
4667  w[ 0] = 0.021616013526483310313342710266;
4668  w[ 1] = 0.049714548894969796453334946203;
4669  w[ 2] = 0.07642573025488905652912967762;
4670  w[ 3] = 0.10094204410628716556281398492;
4671  w[ 4] = 0.12255520671147846018451912680;
4672  w[ 5] = 0.14064291467065065120473130375;
4673  w[ 6] = 0.15468467512626524492541800384;
4674  w[ 7] = 0.16427648374583272298605377647;
4675  w[ 8] = 0.16914238296314359184065647013;
4676  w[ 9] = 0.16914238296314359184065647013;
4677  w[10] = 0.16427648374583272298605377647;
4678  w[11] = 0.15468467512626524492541800384;
4679  w[12] = 0.14064291467065065120473130375;
4680  w[13] = 0.12255520671147846018451912680;
4681  w[14] = 0.10094204410628716556281398492;
4682  w[15] = 0.07642573025488905652912967762;
4683  w[16] = 0.049714548894969796453334946203;
4684  w[17] = 0.021616013526483310313342710266;
4685  }
4686  else if (n==19) {
4687  w[ 0] = 0.019461788229726477036312041464;
4688  w[ 1] = 0.044814226765699600332838157402;
4689  w[ 2] = 0.069044542737641226580708258006;
4690  w[ 3] = 0.091490021622449999464462094124;
4691  w[ 4] = 0.111566645547333994716023901682;
4692  w[ 5] = 0.12875396253933622767551578486;
4693  w[ 6] = 0.14260670217360661177574610944;
4694  w[ 7] = 0.15276604206585966677885540090;
4695  w[ 8] = 0.15896884339395434764995643946;
4696  w[ 9] = 0.16105444984878369597916362532;
4697  w[10] = 0.15896884339395434764995643946;
4698  w[11] = 0.15276604206585966677885540090;
4699  w[12] = 0.14260670217360661177574610944;
4700  w[13] = 0.12875396253933622767551578486;
4701  w[14] = 0.111566645547333994716023901682;
4702  w[15] = 0.091490021622449999464462094124;
4703  w[16] = 0.069044542737641226580708258006;
4704  w[17] = 0.044814226765699600332838157402;
4705  w[18] = 0.019461788229726477036312041464;
4706  }
4707  else if (n==20) {
4708  w[ 0] = 0.017614007139152118311861962352;
4709  w[ 1] = 0.040601429800386941331039952275;
4710  w[ 2] = 0.062672048334109063569506535187;
4711  w[ 3] = 0.08327674157670474872475814322;
4712  w[ 4] = 0.10193011981724043503675013548;
4713  w[ 5] = 0.11819453196151841731237737771;
4714  w[ 6] = 0.13168863844917662689849449975;
4715  w[ 7] = 0.14209610931838205132929832507;
4716  w[ 8] = 0.14917298647260374678782873700;
4717  w[ 9] = 0.15275338713072585069808433195;
4718  w[10] = 0.15275338713072585069808433195;
4719  w[11] = 0.14917298647260374678782873700;
4720  w[12] = 0.14209610931838205132929832507;
4721  w[13] = 0.13168863844917662689849449975;
4722  w[14] = 0.11819453196151841731237737771;
4723  w[15] = 0.10193011981724043503675013548;
4724  w[16] = 0.08327674157670474872475814322;
4725  w[17] = 0.062672048334109063569506535187;
4726  w[18] = 0.040601429800386941331039952275;
4727  w[19] = 0.017614007139152118311861962352;
4728  }
4729  else if (n==21) {
4730  w[ 0] = 0.016017228257774333324224616858;
4731  w[ 1] = 0.036953789770852493799950668299;
4732  w[ 2] = 0.057134425426857208283635826472;
4733  w[ 3] = 0.076100113628379302017051653300;
4734  w[ 4] = 0.093444423456033861553289741114;
4735  w[ 5] = 0.108797299167148377663474578070;
4736  w[ 6] = 0.12183141605372853419536717713;
4737  w[ 7] = 0.13226893863333746178105257450;
4738  w[ 8] = 0.13988739479107315472213342387;
4739  w[ 9] = 0.14452440398997005906382716655;
4740  w[10] = 0.14608113364969042719198514768;
4741  w[11] = 0.14452440398997005906382716655;
4742  w[12] = 0.13988739479107315472213342387;
4743  w[13] = 0.13226893863333746178105257450;
4744  w[14] = 0.12183141605372853419536717713;
4745  w[15] = 0.108797299167148377663474578070;
4746  w[16] = 0.093444423456033861553289741114;
4747  w[17] = 0.076100113628379302017051653300;
4748  w[18] = 0.057134425426857208283635826472;
4749  w[19] = 0.036953789770852493799950668299;
4750  w[20] = 0.016017228257774333324224616858;
4751  }
4752  else if (n==22) {
4753  w[ 0] = 0.014627995298272200684991098047;
4754  w[ 1] = 0.033774901584814154793302246866;
4755  w[ 2] = 0.052293335152683285940312051273;
4756  w[ 3] = 0.06979646842452048809496141893;
4757  w[ 4] = 0.08594160621706772741444368137;
4758  w[ 5] = 0.10041414444288096493207883783;
4759  w[ 6] = 0.11293229608053921839340060742;
4760  w[ 7] = 0.12325237681051242428556098615;
4761  w[ 8] = 0.13117350478706237073296499253;
4762  w[ 9] = 0.13654149834601517135257383123;
4763  w[10] = 0.13925187285563199337541024834;
4764  w[11] = 0.13925187285563199337541024834;
4765  w[12] = 0.13654149834601517135257383123;
4766  w[13] = 0.13117350478706237073296499253;
4767  w[14] = 0.12325237681051242428556098615;
4768  w[15] = 0.11293229608053921839340060742;
4769  w[16] = 0.10041414444288096493207883783;
4770  w[17] = 0.08594160621706772741444368137;
4771  w[18] = 0.06979646842452048809496141893;
4772  w[19] = 0.052293335152683285940312051273;
4773  w[20] = 0.033774901584814154793302246866;
4774  w[21] = 0.014627995298272200684991098047;
4775  }
4776  else if (n==23) {
4777  w[ 0] = 0.013411859487141772081309493459;
4778  w[ 1] = 0.030988005856979444310694219642;
4779  w[ 2] = 0.048037671731084668571641071632;
4780  w[ 3] = 0.064232421408525852127169615159;
4781  w[ 4] = 0.079281411776718954922892524742;
4782  w[ 5] = 0.092915766060035147477018617370;
4783  w[ 6] = 0.104892091464541410074086185015;
4784  w[ 7] = 0.11499664022241136494164351293;
4785  w[ 8] = 0.12304908430672953046757840067;
4786  w[ 9] = 0.12890572218808214997859533940;
4787  w[10] = 0.13246203940469661737164246470;
4788  w[11] = 0.13365457218610617535145711055;
4789  w[12] = 0.13246203940469661737164246470;
4790  w[13] = 0.12890572218808214997859533940;
4791  w[14] = 0.12304908430672953046757840067;
4792  w[15] = 0.11499664022241136494164351293;
4793  w[16] = 0.104892091464541410074086185015;
4794  w[17] = 0.092915766060035147477018617370;
4795  w[18] = 0.079281411776718954922892524742;
4796  w[19] = 0.064232421408525852127169615159;
4797  w[20] = 0.048037671731084668571641071632;
4798  w[21] = 0.030988005856979444310694219642;
4799  w[22] = 0.013411859487141772081309493459;
4800  }
4801  else if (n==24) {
4802  w[ 0] = 0.012341229799987199546805667070;
4803  w[ 1] = 0.028531388628933663181307815952;
4804  w[ 2] = 0.044277438817419806168602748211;
4805  w[ 3] = 0.059298584915436780746367758500;
4806  w[ 4] = 0.07334648141108030573403361525;
4807  w[ 5] = 0.08619016153195327591718520298;
4808  w[ 6] = 0.09761865210411388826988066446;
4809  w[ 7] = 0.10744427011596563478257734245;
4810  w[ 8] = 0.11550566805372560135334448391;
4811  w[ 9] = 0.12167047292780339120446315348;
4812  w[10] = 0.12583745634682829612137538251;
4813  w[11] = 0.12793819534675215697405616522;
4814  w[12] = 0.12793819534675215697405616522;
4815  w[13] = 0.12583745634682829612137538251;
4816  w[14] = 0.12167047292780339120446315348;
4817  w[15] = 0.11550566805372560135334448391;
4818  w[16] = 0.10744427011596563478257734245;
4819  w[17] = 0.09761865210411388826988066446;
4820  w[18] = 0.08619016153195327591718520298;
4821  w[19] = 0.07334648141108030573403361525;
4822  w[20] = 0.059298584915436780746367758500;
4823  w[21] = 0.044277438817419806168602748211;
4824  w[22] = 0.028531388628933663181307815952;
4825  w[23] = 0.012341229799987199546805667070;
4826  }
4827  else if (n==25) {
4828  w[ 0] = 0.0113937985010262879479029641132;
4829  w[ 1] = 0.026354986615032137261901815295;
4830  w[ 2] = 0.040939156701306312655623487712;
4831  w[ 3] = 0.054904695975835191925936891541;
4832  w[ 4] = 0.068038333812356917207187185657;
4833  w[ 5] = 0.080140700335001018013234959669;
4834  w[ 6] = 0.091028261982963649811497220703;
4835  w[ 7] = 0.100535949067050644202206890393;
4836  w[ 8] = 0.108519624474263653116093957050;
4837  w[ 9] = 0.11485825914571164833932554587;
4838  w[10] = 0.11945576353578477222817812651;
4839  w[11] = 0.12224244299031004168895951895;
4840  w[12] = 0.12317605372671545120390287308;
4841  w[13] = 0.12224244299031004168895951895;
4842  w[14] = 0.11945576353578477222817812651;
4843  w[15] = 0.11485825914571164833932554587;
4844  w[16] = 0.108519624474263653116093957050;
4845  w[17] = 0.100535949067050644202206890393;
4846  w[18] = 0.091028261982963649811497220703;
4847  w[19] = 0.080140700335001018013234959669;
4848  w[20] = 0.068038333812356917207187185657;
4849  w[21] = 0.054904695975835191925936891541;
4850  w[22] = 0.040939156701306312655623487712;
4851  w[23] = 0.026354986615032137261901815295;
4852  w[24] = 0.0113937985010262879479029641132;
4853  }
4854  else if (n==26) {
4855  w[ 0] = 0.010551372617343007155651187685;
4856  w[ 1] = 0.024417851092631908789615827520;
4857  w[ 2] = 0.037962383294362763950303141249;
4858  w[ 3] = 0.050975825297147811998319900724;
4859  w[ 4] = 0.063274046329574835539453689907;
4860  w[ 5] = 0.07468414976565974588707579610;
4861  w[ 6] = 0.08504589431348523921044776508;
4862  w[ 7] = 0.09421380035591414846366488307;
4863  w[ 8] = 0.10205916109442542323841407025;
4864  w[ 9] = 0.10847184052857659065657942673;
4865  w[10] = 0.11336181654631966654944071844;
4866  w[11] = 0.11666044348529658204466250754;
4867  w[12] = 0.11832141527926227651637108570;
4868  w[13] = 0.11832141527926227651637108570;
4869  w[14] = 0.11666044348529658204466250754;
4870  w[15] = 0.11336181654631966654944071844;
4871  w[16] = 0.10847184052857659065657942673;
4872  w[17] = 0.10205916109442542323841407025;
4873  w[18] = 0.09421380035591414846366488307;
4874  w[19] = 0.08504589431348523921044776508;
4875  w[20] = 0.07468414976565974588707579610;
4876  w[21] = 0.063274046329574835539453689907;
4877  w[22] = 0.050975825297147811998319900724;
4878  w[23] = 0.037962383294362763950303141249;
4879  w[24] = 0.024417851092631908789615827520;
4880  w[25] = 0.010551372617343007155651187685;
4881  }
4882  else if (n==27) {
4883  w[ 0] = 0.0097989960512943602611500550912;
4884  w[ 1] = 0.022686231596180623196034206447;
4885  w[ 2] = 0.035297053757419711022578289305;
4886  w[ 3] = 0.047449412520615062704096710114;
4887  w[ 4] = 0.058983536859833599110300833720;
4888  w[ 5] = 0.069748823766245592984322888357;
4889  w[ 6] = 0.079604867773057771263074959010;
4890  w[ 7] = 0.088423158543756950194322802854;
4891  w[ 8] = 0.096088727370028507565652646558;
4892  w[ 9] = 0.102501637817745798671247711533;
4893  w[10] = 0.107578285788533187212162984427;
4894  w[11] = 0.111252488356845192672163096043;
4895  w[12] = 0.113476346108965148620369948092;
4896  w[13] = 0.11422086737895698904504573690;
4897  w[14] = 0.113476346108965148620369948092;
4898  w[15] = 0.111252488356845192672163096043;
4899  w[16] = 0.107578285788533187212162984427;
4900  w[17] = 0.102501637817745798671247711533;
4901  w[18] = 0.096088727370028507565652646558;
4902  w[19] = 0.088423158543756950194322802854;
4903  w[20] = 0.079604867773057771263074959010;
4904  w[21] = 0.069748823766245592984322888357;
4905  w[22] = 0.058983536859833599110300833720;
4906  w[23] = 0.047449412520615062704096710114;
4907  w[24] = 0.035297053757419711022578289305;
4908  w[25] = 0.022686231596180623196034206447;
4909  w[26] = 0.0097989960512943602611500550912;
4910  }
4911  else if (n==28) {
4912  w[ 0] = 0.009124282593094517738816153923;
4913  w[ 1] = 0.021132112592771259751500380993;
4914  w[ 2] = 0.032901427782304379977630819171;
4915  w[ 3] = 0.044272934759004227839587877653;
4916  w[ 4] = 0.055107345675716745431482918227;
4917  w[ 5] = 0.06527292396699959579339756678;
4918  w[ 6] = 0.07464621423456877902393188717;
4919  w[ 7] = 0.08311341722890121839039649824;
4920  w[ 8] = 0.09057174439303284094218603134;
4921  w[ 9] = 0.09693065799792991585048900610;
4922  w[10] = 0.10211296757806076981421663851;
4923  w[11] = 0.10605576592284641791041643700;
4924  w[12] = 0.10871119225829413525357151930;
4925  w[13] = 0.11004701301647519628237626560;
4926  w[14] = 0.11004701301647519628237626560;
4927  w[15] = 0.10871119225829413525357151930;
4928  w[16] = 0.10605576592284641791041643700;
4929  w[17] = 0.10211296757806076981421663851;
4930  w[18] = 0.09693065799792991585048900610;
4931  w[19] = 0.09057174439303284094218603134;
4932  w[20] = 0.08311341722890121839039649824;
4933  w[21] = 0.07464621423456877902393188717;
4934  w[22] = 0.06527292396699959579339756678;
4935  w[23] = 0.055107345675716745431482918227;
4936  w[24] = 0.044272934759004227839587877653;
4937  w[25] = 0.032901427782304379977630819171;
4938  w[26] = 0.021132112592771259751500380993;
4939  w[27] = 0.009124282593094517738816153923;
4940  }
4941  else if (n==29) {
4942  w[ 0] = 0.0085169038787464096542638133022;
4943  w[ 1] = 0.019732085056122705983859801640;
4944  w[ 2] = 0.030740492202093622644408525375;
4945  w[ 3] = 0.041402062518682836104830010114;
4946  w[ 4] = 0.051594826902497923912594381180;
4947  w[ 5] = 0.061203090657079138542109848024;
4948  w[ 6] = 0.070117933255051278569581486949;
4949  w[ 7] = 0.078238327135763783828144888660;
4950  w[ 8] = 0.085472257366172527545344849297;
4951  w[ 9] = 0.091737757139258763347966411077;
4952  w[10] = 0.096963834094408606301900074883;
4953  w[11] = 0.101091273759914966121820546907;
4954  w[12] = 0.104073310077729373913328471285;
4955  w[13] = 0.105876155097320941406591327852;
4956  w[14] = 0.10647938171831424424651112691;
4957  w[15] = 0.105876155097320941406591327852;
4958  w[16] = 0.104073310077729373913328471285;
4959  w[17] = 0.101091273759914966121820546907;
4960  w[18] = 0.096963834094408606301900074883;
4961  w[19] = 0.091737757139258763347966411077;
4962  w[20] = 0.085472257366172527545344849297;
4963  w[21] = 0.078238327135763783828144888660;
4964  w[22] = 0.070117933255051278569581486949;
4965  w[23] = 0.061203090657079138542109848024;
4966  w[24] = 0.051594826902497923912594381180;
4967  w[25] = 0.041402062518682836104830010114;
4968  w[26] = 0.030740492202093622644408525375;
4969  w[27] = 0.019732085056122705983859801640;
4970  w[28] = 0.0085169038787464096542638133022;
4971  }
4972  else if (n==30) {
4973  w[ 0] = 0.007968192496166605615465883475;
4974  w[ 1] = 0.018466468311090959142302131912;
4975  w[ 2] = 0.028784707883323369349719179611;
4976  w[ 3] = 0.038799192569627049596801936446;
4977  w[ 4] = 0.048402672830594052902938140423;
4978  w[ 5] = 0.057493156217619066481721689402;
4979  w[ 6] = 0.06597422988218049512812851512;
4980  w[ 7] = 0.07375597473770520626824385002;
4981  w[ 8] = 0.08075589522942021535469493846;
4982  w[ 9] = 0.08689978720108297980238753072;
4983  w[10] = 0.09212252223778612871763270709;
4984  w[11] = 0.09636873717464425963946862635;
4985  w[12] = 0.09959342058679526706278028210;
4986  w[13] = 0.10176238974840550459642895217;
4987  w[14] = 0.10285265289355884034128563671;
4988  w[15] = 0.10285265289355884034128563671;
4989  w[16] = 0.10176238974840550459642895217;
4990  w[17] = 0.09959342058679526706278028210;
4991  w[18] = 0.09636873717464425963946862635;
4992  w[19] = 0.09212252223778612871763270709;
4993  w[20] = 0.08689978720108297980238753072;
4994  w[21] = 0.08075589522942021535469493846;
4995  w[22] = 0.07375597473770520626824385002;
4996  w[23] = 0.06597422988218049512812851512;
4997  w[24] = 0.057493156217619066481721689402;
4998  w[25] = 0.048402672830594052902938140423;
4999  w[26] = 0.038799192569627049596801936446;
5000  w[27] = 0.028784707883323369349719179611;
5001  w[28] = 0.018466468311090959142302131912;
5002  w[29] = 0.007968192496166605615465883475;
5003  }
5004  else if (n==31) {
5005  w[ 0] = 0.0074708315792487758586968750322;
5006  w[ 1] = 0.017318620790310582463157996087;
5007  w[ 2] = 0.027009019184979421800608708092;
5008  w[ 3] = 0.036432273912385464024392010468;
5009  w[ 4] = 0.045493707527201102902315857895;
5010  w[ 5] = 0.054103082424916853711666259087;
5011  w[ 6] = 0.062174786561028426910343543687;
5012  w[ 7] = 0.069628583235410366167756126255;
5013  w[ 8] = 0.076390386598776616426357674901;
5014  w[ 9] = 0.082392991761589263903823367432;
5015  w[10] = 0.087576740608477876126198069695;
5016  w[11] = 0.091890113893641478215362871607;
5017  w[12] = 0.095290242912319512807204197488;
5018  w[13] = 0.097743335386328725093474010979;
5019  w[14] = 0.099225011226672307874875514429;
5020  w[15] = 0.09972054479342645142753383373;
5021  w[16] = 0.099225011226672307874875514429;
5022  w[17] = 0.097743335386328725093474010979;
5023  w[18] = 0.095290242912319512807204197488;
5024  w[19] = 0.091890113893641478215362871607;
5025  w[20] = 0.087576740608477876126198069695;
5026  w[21] = 0.082392991761589263903823367432;
5027  w[22] = 0.076390386598776616426357674901;
5028  w[23] = 0.069628583235410366167756126255;
5029  w[24] = 0.062174786561028426910343543687;
5030  w[25] = 0.054103082424916853711666259087;
5031  w[26] = 0.045493707527201102902315857895;
5032  w[27] = 0.036432273912385464024392010468;
5033  w[28] = 0.027009019184979421800608708092;
5034  w[29] = 0.017318620790310582463157996087;
5035  w[30] = 0.0074708315792487758586968750322;
5036  }
5037  else if (n==32) {
5038  w[ 0] = 0.007018610009470096600407063739;
5039  w[ 1] = 0.016274394730905670605170562206;
5040  w[ 2] = 0.025392065309262059455752589789;
5041  w[ 3] = 0.034273862913021433102687732252;
5042  w[ 4] = 0.042835898022226680656878646606;
5043  w[ 5] = 0.050998059262376176196163244690;
5044  w[ 6] = 0.058684093478535547145283637300;
5045  w[ 7] = 0.06582222277636184683765006371;
5046  w[ 8] = 0.07234579410884850622539935648;
5047  w[ 9] = 0.07819389578707030647174091883;
5048  w[10] = 0.08331192422694675522219907460;
5049  w[11] = 0.08765209300440381114277146275;
5050  w[12] = 0.09117387869576388471286857711;
5051  w[13] = 0.09384439908080456563918023767;
5052  w[14] = 0.09563872007927485941908200220;
5053  w[15] = 0.09654008851472780056676483006;
5054  w[16] = 0.09654008851472780056676483006;
5055  w[17] = 0.09563872007927485941908200220;
5056  w[18] = 0.09384439908080456563918023767;
5057  w[19] = 0.09117387869576388471286857711;
5058  w[20] = 0.08765209300440381114277146275;
5059  w[21] = 0.08331192422694675522219907460;
5060  w[22] = 0.07819389578707030647174091883;
5061  w[23] = 0.07234579410884850622539935648;
5062  w[24] = 0.06582222277636184683765006371;
5063  w[25] = 0.058684093478535547145283637300;
5064  w[26] = 0.050998059262376176196163244690;
5065  w[27] = 0.042835898022226680656878646606;
5066  w[28] = 0.034273862913021433102687732252;
5067  w[29] = 0.025392065309262059455752589789;
5068  w[30] = 0.016274394730905670605170562206;
5069  w[31] = 0.007018610009470096600407063739;
5070  }
5071  else if (n==33) {
5072  w[ 0] = 0.0066062278475873780586492352085;
5073  w[ 1] = 0.015321701512934676127945768534;
5074  w[ 2] = 0.023915548101749480350533257529;
5075  w[ 3] = 0.032300358632328953281561447250;
5076  w[ 4] = 0.040401541331669591563409790527;
5077  w[ 5] = 0.048147742818711695670146880138;
5078  w[ 6] = 0.055470846631663561284944495439;
5079  w[ 7] = 0.062306482530317480031627725771;
5080  w[ 8] = 0.068594572818656712805955073015;
5081  w[ 9] = 0.074279854843954149342472175919;
5082  w[10] = 0.079312364794886738363908384942;
5083  w[11] = 0.083647876067038707613928014518;
5084  w[12] = 0.087248287618844337607281670945;
5085  w[13] = 0.090081958660638577239743705500;
5086  w[14] = 0.092123986643316846213240977717;
5087  w[15] = 0.093356426065596116160999126274;
5088  w[16] = 0.09376844616020999656730454155;
5089  w[17] = 0.093356426065596116160999126274;
5090  w[18] = 0.092123986643316846213240977717;
5091  w[19] = 0.090081958660638577239743705500;
5092  w[20] = 0.087248287618844337607281670945;
5093  w[21] = 0.083647876067038707613928014518;
5094  w[22] = 0.079312364794886738363908384942;
5095  w[23] = 0.074279854843954149342472175919;
5096  w[24] = 0.068594572818656712805955073015;
5097  w[25] = 0.062306482530317480031627725771;
5098  w[26] = 0.055470846631663561284944495439;
5099  w[27] = 0.048147742818711695670146880138;
5100  w[28] = 0.040401541331669591563409790527;
5101  w[29] = 0.032300358632328953281561447250;
5102  w[30] = 0.023915548101749480350533257529;
5103  w[31] = 0.015321701512934676127945768534;
5104  w[32] = 0.0066062278475873780586492352085;
5105  }
5106  else {
5107  std::cerr << "\n";
5108  std::cerr << "LEGENDRE_LOOKUP_WEIGHTS - Fatal error!\n";
5109  std::cerr << " Illegal value of N = " << n << "\n";
5110  std::cerr << " Legal values are 1 through 33.\n";
5111  std::exit(1);
5112  }
5113  return;
5114 }
5115 
5116 //****************************************************************************
5117 template<class Scalar>
5118 void IntrepidBurkardtRules::patterson_lookup ( int n, Scalar x[], Scalar w[] )
5119 //****************************************************************************
5120 //
5121 // Purpose:
5122 //
5123 // PATTERSON_LOOKUP looks up Patterson quadrature points and weights.
5124 //
5125 // Discussion:
5126 //
5127 // Our convention is that the abscissas are numbered from left to right.
5128 //
5129 // The rule is defined on [-1,1],
5130 //
5131 // Licensing:
5132 //
5133 // This code is distributed under the GNU LGPL license.
5134 //
5135 // Modified:
5136 //
5137 // 11 February 2010
5138 //
5139 // Author:
5140 //
5141 // John Burkardt
5142 //
5143 // Reference:
5144 //
5145 // Prem Kythe, Michael Schaeferkotter,
5146 // Handbook of Computational Methods for Integration,
5147 // Chapman and Hall, 2004,
5148 // ISBN: 1-58488-428-2,
5149 // LC: QA299.3.K98.
5150 //
5151 // Thomas Patterson,
5152 // The Optimal Addition of Points to Quadrature Formulae,
5153 // Mathematics of Computation,
5154 // Volume 22, Number 104, October 1968, pages 847-856.
5155 //
5156 // Parameters:
5157 //
5158 // Input, int N, the order.
5159 // Legal values are 1, 3, 7, 15, 31, 63, 127 and 255.
5160 //
5161 // Output, Scalar X[N], the abscissas.
5162 //
5163 // Output, Scalar W[N], the weights.
5164 //
5165 {
5168 
5169  return;
5170 }
5171 
5172 //****************************************************************************
5173 template<class Scalar>
5175 //****************************************************************************
5176 //
5177 // Purpose:
5178 //
5179 // PATTERSON_LOOKUP_POINTS looks up Patterson quadrature points.
5180 //
5181 // Discussion:
5182 //
5183 // Our convention is that the abscissas are numbered from left to right.
5184 //
5185 // The rule is defined on [-1,1],
5186 //
5187 // Licensing:
5188 //
5189 // This code is distributed under the GNU LGPL license.
5190 //
5191 // Modified:
5192 //
5193 // 17 December 2009
5194 //
5195 // Author:
5196 //
5197 // John Burkardt
5198 //
5199 // Reference:
5200 //
5201 // Prem Kythe, Michael Schaeferkotter,
5202 // Handbook of Computational Methods for Integration,
5203 // Chapman and Hall, 2004,
5204 // ISBN: 1-58488-428-2,
5205 // LC: QA299.3.K98.
5206 //
5207 // Thomas Patterson,
5208 // The Optimal Addition of Points to Quadrature Formulae,
5209 // Mathematics of Computation,
5210 // Volume 22, Number 104, October 1968, pages 847-856.
5211 //
5212 // Parameters:
5213 //
5214 // Input, int N, the order.
5215 // Legal values are 1, 3, 7, 15, 31, 63, 127 and 255.
5216 //
5217 // Output, Scalar X[N], the abscissas.
5218 //
5219 {
5220  if (n==1) {
5221  x[ 0] = 0.0;
5222  }
5223  else if (n==3) {
5224  x[ 0] = -0.77459666924148337704;
5225  x[ 1] = 0.0;
5226  x[ 2] = 0.77459666924148337704;
5227  }
5228  else if (n==7) {
5229  x[ 0] = -0.96049126870802028342;
5230  x[ 1] = -0.77459666924148337704;
5231  x[ 2] = -0.43424374934680255800;
5232  x[ 3] = 0.0;
5233  x[ 4] = 0.43424374934680255800;
5234  x[ 5] = 0.77459666924148337704;
5235  x[ 6] = 0.96049126870802028342;
5236  }
5237  else if (n==15) {
5238  x[ 0] = -0.99383196321275502221;
5239  x[ 1] = -0.96049126870802028342;
5240  x[ 2] = -0.88845923287225699889;
5241  x[ 3] = -0.77459666924148337704;
5242  x[ 4] = -0.62110294673722640294;
5243  x[ 5] = -0.43424374934680255800;
5244  x[ 6] = -0.22338668642896688163;
5245  x[ 7] = 0.0;
5246  x[ 8] = 0.22338668642896688163;
5247  x[ 9] = 0.43424374934680255800;
5248  x[10] = 0.62110294673722640294;
5249  x[11] = 0.77459666924148337704;
5250  x[12] = 0.88845923287225699889;
5251  x[13] = 0.96049126870802028342;
5252  x[14] = 0.99383196321275502221;
5253  }
5254  else if (n==31) {
5255  x[ 0] = -0.99909812496766759766;
5256  x[ 1] = -0.99383196321275502221;
5257  x[ 2] = -0.98153114955374010687;
5258  x[ 3] = -0.96049126870802028342;
5259  x[ 4] = -0.92965485742974005667;
5260  x[ 5] = -0.88845923287225699889;
5261  x[ 6] = -0.83672593816886873550;
5262  x[ 7] = -0.77459666924148337704;
5263  x[ 8] = -0.70249620649152707861;
5264  x[ 9] = -0.62110294673722640294;
5265  x[10] = -0.53131974364437562397;
5266  x[11] = -0.43424374934680255800;
5267  x[12] = -0.33113539325797683309;
5268  x[13] = -0.22338668642896688163;
5269  x[14] = -0.11248894313318662575;
5270  x[15] = 0.0;
5271  x[16] = 0.11248894313318662575;
5272  x[17] = 0.22338668642896688163;
5273  x[18] = 0.33113539325797683309;
5274  x[19] = 0.43424374934680255800;
5275  x[20] = 0.53131974364437562397;
5276  x[21] = 0.62110294673722640294;
5277  x[22] = 0.70249620649152707861;
5278  x[23] = 0.77459666924148337704;
5279  x[24] = 0.83672593816886873550;
5280  x[25] = 0.88845923287225699889;
5281  x[26] = 0.92965485742974005667;
5282  x[27] = 0.96049126870802028342;
5283  x[28] = 0.98153114955374010687;
5284  x[29] = 0.99383196321275502221;
5285  x[30] = 0.99909812496766759766;
5286  }
5287  else if (n==63) {
5288  x[ 0] = -0.99987288812035761194;
5289  x[ 1] = -0.99909812496766759766;
5290  x[ 2] = -0.99720625937222195908;
5291  x[ 3] = -0.99383196321275502221;
5292  x[ 4] = -0.98868475754742947994;
5293  x[ 5] = -0.98153114955374010687;
5294  x[ 6] = -0.97218287474858179658;
5295  x[ 7] = -0.96049126870802028342;
5296  x[ 8] = -0.94634285837340290515;
5297  x[ 9] = -0.92965485742974005667;
5298  x[10] = -0.91037115695700429250;
5299  x[11] = -0.88845923287225699889;
5300  x[12] = -0.86390793819369047715;
5301  x[13] = -0.83672593816886873550;
5302  x[14] = -0.80694053195021761186;
5303  x[15] = -0.77459666924148337704;
5304  x[16] = -0.73975604435269475868;
5305  x[17] = -0.70249620649152707861;
5306  x[18] = -0.66290966002478059546;
5307  x[19] = -0.62110294673722640294;
5308  x[20] = -0.57719571005204581484;
5309  x[21] = -0.53131974364437562397;
5310  x[22] = -0.48361802694584102756;
5311  x[23] = -0.43424374934680255800;
5312  x[24] = -0.38335932419873034692;
5313  x[25] = -0.33113539325797683309;
5314  x[26] = -0.27774982202182431507;
5315  x[27] = -0.22338668642896688163;
5316  x[28] = -0.16823525155220746498;
5317  x[29] = -0.11248894313318662575;
5318  x[30] = -0.056344313046592789972;
5319  x[31] = 0.0;
5320  x[32] = 0.056344313046592789972;
5321  x[33] = 0.11248894313318662575;
5322  x[34] = 0.16823525155220746498;
5323  x[35] = 0.22338668642896688163;
5324  x[36] = 0.27774982202182431507;
5325  x[37] = 0.33113539325797683309;
5326  x[38] = 0.38335932419873034692;
5327  x[39] = 0.43424374934680255800;
5328  x[40] = 0.48361802694584102756;
5329  x[41] = 0.53131974364437562397;
5330  x[42] = 0.57719571005204581484;
5331  x[43] = 0.62110294673722640294;
5332  x[44] = 0.66290966002478059546;
5333  x[45] = 0.70249620649152707861;
5334  x[46] = 0.73975604435269475868;
5335  x[47] = 0.77459666924148337704;
5336  x[48] = 0.80694053195021761186;
5337  x[49] = 0.83672593816886873550;
5338  x[50] = 0.86390793819369047715;
5339  x[51] = 0.88845923287225699889;
5340  x[52] = 0.91037115695700429250;
5341  x[53] = 0.92965485742974005667;
5342  x[54] = 0.94634285837340290515;
5343  x[55] = 0.96049126870802028342;
5344  x[56] = 0.97218287474858179658;
5345  x[57] = 0.98153114955374010687;
5346  x[58] = 0.98868475754742947994;
5347  x[59] = 0.99383196321275502221;
5348  x[60] = 0.99720625937222195908;
5349  x[61] = 0.99909812496766759766;
5350  x[62] = 0.99987288812035761194;
5351  }
5352  else if (n==127) {
5353  x[ 0] = -0.99998243035489159858;
5354  x[ 1] = -0.99987288812035761194;
5355  x[ 2] = -0.99959879967191068325;
5356  x[ 3] = -0.99909812496766759766;
5357  x[ 4] = -0.99831663531840739253;
5358  x[ 5] = -0.99720625937222195908;
5359  x[ 6] = -0.99572410469840718851;
5360  x[ 7] = -0.99383196321275502221;
5361  x[ 8] = -0.99149572117810613240;
5362  x[ 9] = -0.98868475754742947994;
5363  x[ 10] = -0.98537149959852037111;
5364  x[ 11] = -0.98153114955374010687;
5365  x[ 12] = -0.97714151463970571416;
5366  x[ 13] = -0.97218287474858179658;
5367  x[ 14] = -0.96663785155841656709;
5368  x[ 15] = -0.96049126870802028342;
5369  x[ 16] = -0.95373000642576113641;
5370  x[ 17] = -0.94634285837340290515;
5371  x[ 18] = -0.93832039777959288365;
5372  x[ 19] = -0.92965485742974005667;
5373  x[ 20] = -0.92034002547001242073;
5374  x[ 21] = -0.91037115695700429250;
5375  x[ 22] = -0.89974489977694003664;
5376  x[ 23] = -0.88845923287225699889;
5377  x[ 24] = -0.87651341448470526974;
5378  x[ 25] = -0.86390793819369047715;
5379  x[ 26] = -0.85064449476835027976;
5380  x[ 27] = -0.83672593816886873550;
5381  x[ 28] = -0.82215625436498040737;
5382  x[ 29] = -0.80694053195021761186;
5383  x[ 30] = -0.79108493379984836143;
5384  x[ 31] = -0.77459666924148337704;
5385  x[ 32] = -0.75748396638051363793;
5386  x[ 33] = -0.73975604435269475868;
5387  x[ 34] = -0.72142308537009891548;
5388  x[ 35] = -0.70249620649152707861;
5389  x[ 36] = -0.68298743109107922809;
5390  x[ 37] = -0.66290966002478059546;
5391  x[ 38] = -0.64227664250975951377;
5392  x[ 39] = -0.62110294673722640294;
5393  x[ 40] = -0.59940393024224289297;
5394  x[ 41] = -0.57719571005204581484;
5395  x[ 42] = -0.55449513263193254887;
5396  x[ 43] = -0.53131974364437562397;
5397  x[ 44] = -0.50768775753371660215;
5398  x[ 45] = -0.48361802694584102756;
5399  x[ 46] = -0.45913001198983233287;
5400  x[ 47] = -0.43424374934680255800;
5401  x[ 48] = -0.40897982122988867241;
5402  x[ 49] = -0.38335932419873034692;
5403  x[ 50] = -0.35740383783153215238;
5404  x[ 51] = -0.33113539325797683309;
5405  x[ 52] = -0.30457644155671404334;
5406  x[ 53] = -0.27774982202182431507;
5407  x[ 54] = -0.25067873030348317661;
5408  x[ 55] = -0.22338668642896688163;
5409  x[ 56] = -0.19589750271110015392;
5410  x[ 57] = -0.16823525155220746498;
5411  x[ 58] = -0.14042423315256017459;
5412  x[ 59] = -0.11248894313318662575;
5413  x[ 60] = -0.084454040083710883710;
5414  x[ 61] = -0.056344313046592789972;
5415  x[ 62] = -0.028184648949745694339;
5416  x[ 63] = 0.0;
5417  x[ 64] = 0.028184648949745694339;
5418  x[ 65] = 0.056344313046592789972;
5419  x[ 66] = 0.084454040083710883710;
5420  x[ 67] = 0.11248894313318662575;
5421  x[ 68] = 0.14042423315256017459;
5422  x[ 69] = 0.16823525155220746498;
5423  x[ 70] = 0.19589750271110015392;
5424  x[ 71] = 0.22338668642896688163;
5425  x[ 72] = 0.25067873030348317661;
5426  x[ 73] = 0.27774982202182431507;
5427  x[ 74] = 0.30457644155671404334;
5428  x[ 75] = 0.33113539325797683309;
5429  x[ 76] = 0.35740383783153215238;
5430  x[ 77] = 0.38335932419873034692;
5431  x[ 78] = 0.40897982122988867241;
5432  x[ 79] = 0.43424374934680255800;
5433  x[ 80] = 0.45913001198983233287;
5434  x[ 81] = 0.48361802694584102756;
5435  x[ 82] = 0.50768775753371660215;
5436  x[ 83] = 0.53131974364437562397;
5437  x[ 84] = 0.55449513263193254887;
5438  x[ 85] = 0.57719571005204581484;
5439  x[ 86] = 0.59940393024224289297;
5440  x[ 87] = 0.62110294673722640294;
5441  x[ 88] = 0.64227664250975951377;
5442  x[ 89] = 0.66290966002478059546;
5443  x[ 90] = 0.68298743109107922809;
5444  x[ 91] = 0.70249620649152707861;
5445  x[ 92] = 0.72142308537009891548;
5446  x[ 93] = 0.73975604435269475868;
5447  x[ 94] = 0.75748396638051363793;
5448  x[ 95] = 0.77459666924148337704;
5449  x[ 96] = 0.79108493379984836143;
5450  x[ 97] = 0.80694053195021761186;
5451  x[ 98] = 0.82215625436498040737;
5452  x[ 99] = 0.83672593816886873550;
5453  x[100] = 0.85064449476835027976;
5454  x[101] = 0.86390793819369047715;
5455  x[102] = 0.87651341448470526974;
5456  x[103] = 0.88845923287225699889;
5457  x[104] = 0.89974489977694003664;
5458  x[105] = 0.91037115695700429250;
5459  x[106] = 0.92034002547001242073;
5460  x[107] = 0.92965485742974005667;
5461  x[108] = 0.93832039777959288365;
5462  x[109] = 0.94634285837340290515;
5463  x[110] = 0.95373000642576113641;
5464  x[111] = 0.96049126870802028342;
5465  x[112] = 0.96663785155841656709;
5466  x[113] = 0.97218287474858179658;
5467  x[114] = 0.97714151463970571416;
5468  x[115] = 0.98153114955374010687;
5469  x[116] = 0.98537149959852037111;
5470  x[117] = 0.98868475754742947994;
5471  x[118] = 0.99149572117810613240;
5472  x[119] = 0.99383196321275502221;
5473  x[120] = 0.99572410469840718851;
5474  x[121] = 0.99720625937222195908;
5475  x[122] = 0.99831663531840739253;
5476  x[123] = 0.99909812496766759766;
5477  x[124] = 0.99959879967191068325;
5478  x[125] = 0.99987288812035761194;
5479  x[126] = 0.99998243035489159858;
5480  }
5481  else if (n==255) {
5482  x[ 0] = -0.99999759637974846462;
5483  x[ 1] = -0.99998243035489159858;
5484  x[ 2] = -0.99994399620705437576;
5485  x[ 3] = -0.99987288812035761194;
5486  x[ 4] = -0.99976049092443204733;
5487  x[ 5] = -0.99959879967191068325;
5488  x[ 6] = -0.99938033802502358193;
5489  x[ 7] = -0.99909812496766759766;
5490  x[ 8] = -0.99874561446809511470;
5491  x[ 9] = -0.99831663531840739253;
5492  x[ 10] = -0.99780535449595727456;
5493  x[ 11] = -0.99720625937222195908;
5494  x[ 12] = -0.99651414591489027385;
5495  x[ 13] = -0.99572410469840718851;
5496  x[ 14] = -0.99483150280062100052;
5497  x[ 15] = -0.99383196321275502221;
5498  x[ 16] = -0.99272134428278861533;
5499  x[ 17] = -0.99149572117810613240;
5500  x[ 18] = -0.99015137040077015918;
5501  x[ 19] = -0.98868475754742947994;
5502  x[ 20] = -0.98709252795403406719;
5503  x[ 21] = -0.98537149959852037111;
5504  x[ 22] = -0.98351865757863272876;
5505  x[ 23] = -0.98153114955374010687;
5506  x[ 24] = -0.97940628167086268381;
5507  x[ 25] = -0.97714151463970571416;
5508  x[ 26] = -0.97473445975240266776;
5509  x[ 27] = -0.97218287474858179658;
5510  x[ 28] = -0.96948465950245923177;
5511  x[ 29] = -0.96663785155841656709;
5512  x[ 30] = -0.96364062156981213252;
5513  x[ 31] = -0.96049126870802028342;
5514  x[ 32] = -0.95718821610986096274;
5515  x[ 33] = -0.95373000642576113641;
5516  x[ 34] = -0.95011529752129487656;
5517  x[ 35] = -0.94634285837340290515;
5518  x[ 36] = -0.94241156519108305981;
5519  x[ 37] = -0.93832039777959288365;
5520  x[ 38] = -0.93406843615772578800;
5521  x[ 39] = -0.92965485742974005667;
5522  x[ 40] = -0.92507893290707565236;
5523  x[ 41] = -0.92034002547001242073;
5524  x[ 42] = -0.91543758715576504064;
5525  x[ 43] = -0.91037115695700429250;
5526  x[ 44] = -0.90514035881326159519;
5527  x[ 45] = -0.89974489977694003664;
5528  x[ 46] = -0.89418456833555902286;
5529  x[ 47] = -0.88845923287225699889;
5530  x[ 48] = -0.88256884024734190684;
5531  x[ 49] = -0.87651341448470526974;
5532  x[ 50] = -0.87029305554811390585;
5533  x[ 51] = -0.86390793819369047715;
5534  x[ 52] = -0.85735831088623215653;
5535  x[ 53] = -0.85064449476835027976;
5536  x[ 54] = -0.84376688267270860104;
5537  x[ 55] = -0.83672593816886873550;
5538  x[ 56] = -0.82952219463740140018;
5539  x[ 57] = -0.82215625436498040737;
5540  x[ 58] = -0.81462878765513741344;
5541  x[ 59] = -0.80694053195021761186;
5542  x[ 60] = -0.79909229096084140180;
5543  x[ 61] = -0.79108493379984836143;
5544  x[ 62] = -0.78291939411828301639;
5545  x[ 63] = -0.77459666924148337704;
5546  x[ 64] = -0.76611781930376009072;
5547  x[ 65] = -0.75748396638051363793;
5548  x[ 66] = -0.74869629361693660282;
5549  x[ 67] = -0.73975604435269475868;
5550  x[ 68] = -0.73066452124218126133;
5551  x[ 69] = -0.72142308537009891548;
5552  x[ 70] = -0.71203315536225203459;
5553  x[ 71] = -0.70249620649152707861;
5554  x[ 72] = -0.69281376977911470289;
5555  x[ 73] = -0.68298743109107922809;
5556  x[ 74] = -0.67301883023041847920;
5557  x[ 75] = -0.66290966002478059546;
5558  x[ 76] = -0.65266166541001749610;
5559  x[ 77] = -0.64227664250975951377;
5560  x[ 78] = -0.63175643771119423041;
5561  x[ 79] = -0.62110294673722640294;
5562  x[ 80] = -0.61031811371518640016;
5563  x[ 81] = -0.59940393024224289297;
5564  x[ 82] = -0.58836243444766254143;
5565  x[ 83] = -0.57719571005204581484;
5566  x[ 84] = -0.56590588542365442262;
5567  x[ 85] = -0.55449513263193254887;
5568  x[ 86] = -0.54296566649831149049;
5569  x[ 87] = -0.53131974364437562397;
5570  x[ 88] = -0.51955966153745702199;
5571  x[ 89] = -0.50768775753371660215;
5572  x[ 90] = -0.49570640791876146017;
5573  x[ 91] = -0.48361802694584102756;
5574  x[ 92] = -0.47142506587165887693;
5575  x[ 93] = -0.45913001198983233287;
5576  x[ 94] = -0.44673538766202847374;
5577  x[ 95] = -0.43424374934680255800;
5578  x[ 96] = -0.42165768662616330006;
5579  x[ 97] = -0.40897982122988867241;
5580  x[ 98] = -0.39621280605761593918;
5581  x[ 99] = -0.38335932419873034692;
5582  x[100] = -0.37042208795007823014;
5583  x[101] = -0.35740383783153215238;
5584  x[102] = -0.34430734159943802278;
5585  x[103] = -0.33113539325797683309;
5586  x[104] = -0.31789081206847668318;
5587  x[105] = -0.30457644155671404334;
5588  x[106] = -0.29119514851824668196;
5589  x[107] = -0.27774982202182431507;
5590  x[108] = -0.26424337241092676194;
5591  x[109] = -0.25067873030348317661;
5592  x[110] = -0.23705884558982972721;
5593  x[111] = -0.22338668642896688163;
5594  x[112] = -0.20966523824318119477;
5595  x[113] = -0.19589750271110015392;
5596  x[114] = -0.18208649675925219825;
5597  x[115] = -0.16823525155220746498;
5598  x[116] = -0.15434681148137810869;
5599  x[117] = -0.14042423315256017459;
5600  x[118] = -0.12647058437230196685;
5601  x[119] = -0.11248894313318662575;
5602  x[120] = -0.098482396598119202090;
5603  x[121] = -0.084454040083710883710;
5604  x[122] = -0.070406976042855179063;
5605  x[123] = -0.056344313046592789972;
5606  x[124] = -0.042269164765363603212;
5607  x[125] = -0.028184648949745694339;
5608  x[126] = -0.014093886410782462614;
5609  x[127] = 0.0;
5610  x[128] = 0.014093886410782462614;
5611  x[129] = 0.028184648949745694339;
5612  x[130] = 0.042269164765363603212;
5613  x[131] = 0.056344313046592789972;
5614  x[132] = 0.070406976042855179063;
5615  x[133] = 0.084454040083710883710;
5616  x[134] = 0.098482396598119202090;
5617  x[135] = 0.11248894313318662575;
5618  x[136] = 0.12647058437230196685;
5619  x[137] = 0.14042423315256017459;
5620  x[138] = 0.15434681148137810869;
5621  x[139] = 0.16823525155220746498;
5622  x[140] = 0.18208649675925219825;
5623  x[141] = 0.19589750271110015392;
5624  x[142] = 0.20966523824318119477;
5625  x[143] = 0.22338668642896688163;
5626  x[144] = 0.23705884558982972721;
5627  x[145] = 0.25067873030348317661;
5628  x[146] = 0.26424337241092676194;
5629  x[147] = 0.27774982202182431507;
5630  x[148] = 0.29119514851824668196;
5631  x[149] = 0.30457644155671404334;
5632  x[150] = 0.31789081206847668318;
5633  x[151] = 0.33113539325797683309;
5634  x[152] = 0.34430734159943802278;
5635  x[153] = 0.35740383783153215238;
5636  x[154] = 0.37042208795007823014;
5637  x[155] = 0.38335932419873034692;
5638  x[156] = 0.39621280605761593918;
5639  x[157] = 0.40897982122988867241;
5640  x[158] = 0.42165768662616330006;
5641  x[159] = 0.43424374934680255800;
5642  x[160] = 0.44673538766202847374;
5643  x[161] = 0.45913001198983233287;
5644  x[162] = 0.47142506587165887693;
5645  x[163] = 0.48361802694584102756;
5646  x[164] = 0.49570640791876146017;
5647  x[165] = 0.50768775753371660215;
5648  x[166] = 0.51955966153745702199;
5649  x[167] = 0.53131974364437562397;
5650  x[168] = 0.54296566649831149049;
5651  x[169] = 0.55449513263193254887;
5652  x[170] = 0.56590588542365442262;
5653  x[171] = 0.57719571005204581484;
5654  x[172] = 0.58836243444766254143;
5655  x[173] = 0.59940393024224289297;
5656  x[174] = 0.61031811371518640016;
5657  x[175] = 0.62110294673722640294;
5658  x[176] = 0.63175643771119423041;
5659  x[177] = 0.64227664250975951377;
5660  x[178] = 0.65266166541001749610;
5661  x[179] = 0.66290966002478059546;
5662  x[180] = 0.67301883023041847920;
5663  x[181] = 0.68298743109107922809;
5664  x[182] = 0.69281376977911470289;
5665  x[183] = 0.70249620649152707861;
5666  x[184] = 0.71203315536225203459;
5667  x[185] = 0.72142308537009891548;
5668  x[186] = 0.73066452124218126133;
5669  x[187] = 0.73975604435269475868;
5670  x[188] = 0.74869629361693660282;
5671  x[189] = 0.75748396638051363793;
5672  x[190] = 0.76611781930376009072;
5673  x[191] = 0.77459666924148337704;
5674  x[192] = 0.78291939411828301639;
5675  x[193] = 0.79108493379984836143;
5676  x[194] = 0.79909229096084140180;
5677  x[195] = 0.80694053195021761186;
5678  x[196] = 0.81462878765513741344;
5679  x[197] = 0.82215625436498040737;
5680  x[198] = 0.82952219463740140018;
5681  x[199] = 0.83672593816886873550;
5682  x[200] = 0.84376688267270860104;
5683  x[201] = 0.85064449476835027976;
5684  x[202] = 0.85735831088623215653;
5685  x[203] = 0.86390793819369047715;
5686  x[204] = 0.87029305554811390585;
5687  x[205] = 0.87651341448470526974;
5688  x[206] = 0.88256884024734190684;
5689  x[207] = 0.88845923287225699889;
5690  x[208] = 0.89418456833555902286;
5691  x[209] = 0.89974489977694003664;
5692  x[210] = 0.90514035881326159519;
5693  x[211] = 0.91037115695700429250;
5694  x[212] = 0.91543758715576504064;
5695  x[213] = 0.92034002547001242073;
5696  x[214] = 0.92507893290707565236;
5697  x[215] = 0.92965485742974005667;
5698  x[216] = 0.93406843615772578800;
5699  x[217] = 0.93832039777959288365;
5700  x[218] = 0.94241156519108305981;
5701  x[219] = 0.94634285837340290515;
5702  x[220] = 0.95011529752129487656;
5703  x[221] = 0.95373000642576113641;
5704  x[222] = 0.95718821610986096274;
5705  x[223] = 0.96049126870802028342;
5706  x[224] = 0.96364062156981213252;
5707  x[225] = 0.96663785155841656709;
5708  x[226] = 0.96948465950245923177;
5709  x[227] = 0.97218287474858179658;
5710  x[228] = 0.97473445975240266776;
5711  x[229] = 0.97714151463970571416;
5712  x[230] = 0.97940628167086268381;
5713  x[231] = 0.98153114955374010687;
5714  x[232] = 0.98351865757863272876;
5715  x[233] = 0.98537149959852037111;
5716  x[234] = 0.98709252795403406719;
5717  x[235] = 0.98868475754742947994;
5718  x[236] = 0.99015137040077015918;
5719  x[237] = 0.99149572117810613240;
5720  x[238] = 0.99272134428278861533;
5721  x[239] = 0.99383196321275502221;
5722  x[240] = 0.99483150280062100052;
5723  x[241] = 0.99572410469840718851;
5724  x[242] = 0.99651414591489027385;
5725  x[243] = 0.99720625937222195908;
5726  x[244] = 0.99780535449595727456;
5727  x[245] = 0.99831663531840739253;
5728  x[246] = 0.99874561446809511470;
5729  x[247] = 0.99909812496766759766;
5730  x[248] = 0.99938033802502358193;
5731  x[249] = 0.99959879967191068325;
5732  x[250] = 0.99976049092443204733;
5733  x[251] = 0.99987288812035761194;
5734  x[252] = 0.99994399620705437576;
5735  x[253] = 0.99998243035489159858;
5736  x[254] = 0.99999759637974846462;
5737  }
5738  else {
5739  std::cerr << "\n";
5740  std::cerr << "PATTERSON_LOOKUP_POINTS - Fatal error!\n";
5741  std::cerr << " Unexpected value of N = " << n << "\n";
5742  std::exit(1);
5743  }
5744  return;
5745 }
5746 
5747 //****************************************************************************
5748 template<class Scalar>
5750 //****************************************************************************
5751 //
5752 // Purpose:
5753 //
5754 // PATTERSON_LOOKUP_WEIGHTS looks up Patterson quadrature weights.
5755 //
5756 // Discussion:
5757 //
5758 // The allowed orders are 1, 3, 7, 15, 31, 63, 127 and 255.
5759 //
5760 // The weights are positive, symmetric and should sum to 2.
5761 //
5762 // The user must preallocate space for the output array W.
5763 //
5764 // Licensing:
5765 //
5766 // This code is distributed under the GNU LGPL license.
5767 //
5768 // Modified:
5769 //
5770 // 17 December 2009
5771 //
5772 // Author:
5773 //
5774 // John Burkardt
5775 //
5776 // Reference:
5777 //
5778 // Milton Abramowitz, Irene Stegun,
5779 // Handbook of Mathematical Functions,
5780 // National Bureau of Standards, 1964,
5781 // ISBN: 0-486-61272-4,
5782 // LC: QA47.A34.
5783 //
5784 // Arthur Stroud, Don Secrest,
5785 // Gaussian Quadrature Formulas,
5786 // Prentice Hall, 1966,
5787 // LC: QA299.4G3S7.
5788 //
5789 // Parameters:
5790 //
5791 // Input, int N, the order.
5792 // Legal values are 1, 3, 7, 15, 31, 63, 127 or 255.
5793 //
5794 // Output, Scalar W[N], the weights.
5795 //
5796 {
5797  if (n==1) {
5798  w[ 0] = 2.0;
5799  }
5800  else if (n==3) {
5801  w[ 0] = 0.555555555555555555556;
5802  w[ 1] = 0.888888888888888888889;
5803  w[ 2] = 0.555555555555555555556;
5804  }
5805  else if (n==7) {
5806  w[ 0] = 0.104656226026467265194;
5807  w[ 1] = 0.268488089868333440729;
5808  w[ 2] = 0.401397414775962222905;
5809  w[ 3] = 0.450916538658474142345;
5810  w[ 4] = 0.401397414775962222905;
5811  w[ 5] = 0.268488089868333440729;
5812  w[ 6] = 0.104656226026467265194;
5813  }
5814  else if (n==15) {
5815  w[ 0] = 0.0170017196299402603390;
5816  w[ 1] = 0.0516032829970797396969;
5817  w[ 2] = 0.0929271953151245376859;
5818  w[ 3] = 0.134415255243784220360;
5819  w[ 4] = 0.171511909136391380787;
5820  w[ 5] = 0.200628529376989021034;
5821  w[ 6] = 0.219156858401587496404;
5822  w[ 7] = 0.225510499798206687386;
5823  w[ 8] = 0.219156858401587496404;
5824  w[ 9] = 0.200628529376989021034;
5825  w[ 10] = 0.171511909136391380787;
5826  w[ 11] = 0.134415255243784220360;
5827  w[ 12] = 0.0929271953151245376859;
5828  w[ 13] = 0.0516032829970797396969;
5829  w[ 14] = 0.0170017196299402603390;
5830  }
5831  else if (n==31) {
5832  w[ 0] = 0.00254478079156187441540;
5833  w[ 1] = 0.00843456573932110624631;
5834  w[ 2] = 0.0164460498543878109338;
5835  w[ 3] = 0.0258075980961766535646;
5836  w[ 4] = 0.0359571033071293220968;
5837  w[ 5] = 0.0464628932617579865414;
5838  w[ 6] = 0.0569795094941233574122;
5839  w[ 7] = 0.0672077542959907035404;
5840  w[ 8] = 0.0768796204990035310427;
5841  w[ 9] = 0.0857559200499903511542;
5842  w[ 10] = 0.0936271099812644736167;
5843  w[ 11] = 0.100314278611795578771;
5844  w[ 12] = 0.105669893580234809744;
5845  w[ 13] = 0.109578421055924638237;
5846  w[ 14] = 0.111956873020953456880;
5847  w[ 15] = 0.112755256720768691607;
5848  w[ 16] = 0.111956873020953456880;
5849  w[ 17] = 0.109578421055924638237;
5850  w[ 18] = 0.105669893580234809744;
5851  w[ 19] = 0.100314278611795578771;
5852  w[ 20] = 0.0936271099812644736167;
5853  w[ 21] = 0.0857559200499903511542;
5854  w[ 22] = 0.0768796204990035310427;
5855  w[ 23] = 0.0672077542959907035404;
5856  w[ 24] = 0.0569795094941233574122;
5857  w[ 25] = 0.0464628932617579865414;
5858  w[ 26] = 0.0359571033071293220968;
5859  w[ 27] = 0.0258075980961766535646;
5860  w[ 28] = 0.0164460498543878109338;
5861  w[ 29] = 0.00843456573932110624631;
5862  w[ 30] = 0.00254478079156187441540;
5863  }
5864  else if (n==63) {
5865  w[ 0] = 0.000363221481845530659694;
5866  w[ 1] = 0.00126515655623006801137;
5867  w[ 2] = 0.00257904979468568827243;
5868  w[ 3] = 0.00421763044155885483908;
5869  w[ 4] = 0.00611550682211724633968;
5870  w[ 5] = 0.00822300795723592966926;
5871  w[ 6] = 0.0104982469096213218983;
5872  w[ 7] = 0.0129038001003512656260;
5873  w[ 8] = 0.0154067504665594978021;
5874  w[ 9] = 0.0179785515681282703329;
5875  w[ 10] = 0.0205942339159127111492;
5876  w[ 11] = 0.0232314466399102694433;
5877  w[ 12] = 0.0258696793272147469108;
5878  w[ 13] = 0.0284897547458335486125;
5879  w[ 14] = 0.0310735511116879648799;
5880  w[ 15] = 0.0336038771482077305417;
5881  w[ 16] = 0.0360644327807825726401;
5882  w[ 17] = 0.0384398102494555320386;
5883  w[ 18] = 0.0407155101169443189339;
5884  w[ 19] = 0.0428779600250077344929;
5885  w[ 20] = 0.0449145316536321974143;
5886  w[ 21] = 0.0468135549906280124026;
5887  w[ 22] = 0.0485643304066731987159;
5888  w[ 23] = 0.0501571393058995374137;
5889  w[ 24] = 0.0515832539520484587768;
5890  w[ 25] = 0.0528349467901165198621;
5891  w[ 26] = 0.0539054993352660639269;
5892  w[ 27] = 0.0547892105279628650322;
5893  w[ 28] = 0.0554814043565593639878;
5894  w[ 29] = 0.0559784365104763194076;
5895  w[ 30] = 0.0562776998312543012726;
5896  w[ 31] = 0.0563776283603847173877;
5897  w[ 32] = 0.0562776998312543012726;
5898  w[ 33] = 0.0559784365104763194076;
5899  w[ 34] = 0.0554814043565593639878;
5900  w[ 35] = 0.0547892105279628650322;
5901  w[ 36] = 0.0539054993352660639269;
5902  w[ 37] = 0.0528349467901165198621;
5903  w[ 38] = 0.0515832539520484587768;
5904  w[ 39] = 0.0501571393058995374137;
5905  w[ 40] = 0.0485643304066731987159;
5906  w[ 41] = 0.0468135549906280124026;
5907  w[ 42] = 0.0449145316536321974143;
5908  w[ 43] = 0.0428779600250077344929;
5909  w[ 44] = 0.0407155101169443189339;
5910  w[ 45] = 0.0384398102494555320386;
5911  w[ 46] = 0.0360644327807825726401;
5912  w[ 47] = 0.0336038771482077305417;
5913  w[ 48] = 0.0310735511116879648799;
5914  w[ 49] = 0.0284897547458335486125;
5915  w[ 50] = 0.0258696793272147469108;
5916  w[ 51] = 0.0232314466399102694433;
5917  w[ 52] = 0.0205942339159127111492;
5918  w[ 53] = 0.0179785515681282703329;
5919  w[ 54] = 0.0154067504665594978021;
5920  w[ 55] = 0.0129038001003512656260;
5921  w[ 56] = 0.0104982469096213218983;
5922  w[ 57] = 0.00822300795723592966926;
5923  w[ 58] = 0.00611550682211724633968;
5924  w[ 59] = 0.00421763044155885483908;
5925  w[ 60] = 0.00257904979468568827243;
5926  w[ 61] = 0.00126515655623006801137;
5927  w[ 62] = 0.000363221481845530659694;
5928  }
5929  else if (n==127) {
5930  w[ 0] = 0.0000505360952078625176247;
5931  w[ 1] = 0.000180739564445388357820;
5932  w[ 2] = 0.000377746646326984660274;
5933  w[ 3] = 0.000632607319362633544219;
5934  w[ 4] = 0.000938369848542381500794;
5935  w[ 5] = 0.00128952408261041739210;
5936  w[ 6] = 0.00168114286542146990631;
5937  w[ 7] = 0.00210881524572663287933;
5938  w[ 8] = 0.00256876494379402037313;
5939  w[ 9] = 0.00305775341017553113613;
5940  w[ 10] = 0.00357289278351729964938;
5941  w[ 11] = 0.00411150397865469304717;
5942  w[ 12] = 0.00467105037211432174741;
5943  w[ 13] = 0.00524912345480885912513;
5944  w[ 14] = 0.00584344987583563950756;
5945  w[ 15] = 0.00645190005017573692280;
5946  w[ 16] = 0.00707248999543355546805;
5947  w[ 17] = 0.00770337523327974184817;
5948  w[ 18] = 0.00834283875396815770558;
5949  w[ 19] = 0.00898927578406413572328;
5950  w[ 20] = 0.00964117772970253669530;
5951  w[ 21] = 0.0102971169579563555237;
5952  w[ 22] = 0.0109557333878379016480;
5953  w[ 23] = 0.0116157233199551347270;
5954  w[ 24] = 0.0122758305600827700870;
5955  w[ 25] = 0.0129348396636073734547;
5956  w[ 26] = 0.0135915710097655467896;
5957  w[ 27] = 0.0142448773729167743063;
5958  w[ 28] = 0.0148936416648151820348;
5959  w[ 29] = 0.0155367755558439824399;
5960  w[ 30] = 0.0161732187295777199419;
5961  w[ 31] = 0.0168019385741038652709;
5962  w[ 32] = 0.0174219301594641737472;
5963  w[ 33] = 0.0180322163903912863201;
5964  w[ 34] = 0.0186318482561387901863;
5965  w[ 35] = 0.0192199051247277660193;
5966  w[ 36] = 0.0197954950480974994880;
5967  w[ 37] = 0.0203577550584721594669;
5968  w[ 38] = 0.0209058514458120238522;
5969  w[ 39] = 0.0214389800125038672465;
5970  w[ 40] = 0.0219563663053178249393;
5971  w[ 41] = 0.0224572658268160987071;
5972  w[ 42] = 0.0229409642293877487608;
5973  w[ 43] = 0.0234067774953140062013;
5974  w[ 44] = 0.0238540521060385400804;
5975  w[ 45] = 0.0242821652033365993580;
5976  w[ 46] = 0.0246905247444876769091;
5977  w[ 47] = 0.0250785696529497687068;
5978  w[ 48] = 0.0254457699654647658126;
5979  w[ 49] = 0.0257916269760242293884;
5980  w[ 50] = 0.0261156733767060976805;
5981  w[ 51] = 0.0264174733950582599310;
5982  w[ 52] = 0.0266966229274503599062;
5983  w[ 53] = 0.0269527496676330319634;
5984  w[ 54] = 0.0271855132296247918192;
5985  w[ 55] = 0.0273946052639814325161;
5986  w[ 56] = 0.0275797495664818730349;
5987  w[ 57] = 0.0277407021782796819939;
5988  w[ 58] = 0.0278772514766137016085;
5989  w[ 59] = 0.0279892182552381597038;
5990  w[ 60] = 0.0280764557938172466068;
5991  w[ 61] = 0.0281388499156271506363;
5992  w[ 62] = 0.0281763190330166021307;
5993  w[ 63] = 0.0281888141801923586938;
5994  w[ 64] = 0.0281763190330166021307;
5995  w[ 65] = 0.0281388499156271506363;
5996  w[ 66] = 0.0280764557938172466068;
5997  w[ 67] = 0.0279892182552381597038;
5998  w[ 68] = 0.0278772514766137016085;
5999  w[ 69] = 0.0277407021782796819939;
6000  w[ 70] = 0.0275797495664818730349;
6001  w[ 71] = 0.0273946052639814325161;
6002  w[ 72] = 0.0271855132296247918192;
6003  w[ 73] = 0.0269527496676330319634;
6004  w[ 74] = 0.0266966229274503599062;
6005  w[ 75] = 0.0264174733950582599310;
6006  w[ 76] = 0.0261156733767060976805;
6007  w[ 77] = 0.0257916269760242293884;
6008  w[ 78] = 0.0254457699654647658126;
6009  w[ 79] = 0.0250785696529497687068;
6010  w[ 80] = 0.0246905247444876769091;
6011  w[ 81] = 0.0242821652033365993580;
6012  w[ 82] = 0.0238540521060385400804;
6013  w[ 83] = 0.0234067774953140062013;
6014  w[ 84] = 0.0229409642293877487608;
6015  w[ 85] = 0.0224572658268160987071;
6016  w[ 86] = 0.0219563663053178249393;
6017  w[ 87] = 0.0214389800125038672465;
6018  w[ 88] = 0.0209058514458120238522;
6019  w[ 89] = 0.0203577550584721594669;
6020  w[ 90] = 0.0197954950480974994880;
6021  w[ 91] = 0.0192199051247277660193;
6022  w[ 92] = 0.0186318482561387901863;
6023  w[ 93] = 0.0180322163903912863201;
6024  w[ 94] = 0.0174219301594641737472;
6025  w[ 95] = 0.0168019385741038652709;
6026  w[ 96] = 0.0161732187295777199419;
6027  w[ 97] = 0.0155367755558439824399;
6028  w[ 98] = 0.0148936416648151820348;
6029  w[ 99] = 0.0142448773729167743063;
6030  w[100] = 0.0135915710097655467896;
6031  w[101] = 0.0129348396636073734547;
6032  w[102] = 0.0122758305600827700870;
6033  w[103] = 0.0116157233199551347270;
6034  w[104] = 0.0109557333878379016480;
6035  w[105] = 0.0102971169579563555237;
6036  w[106] = 0.00964117772970253669530;
6037  w[107] = 0.00898927578406413572328;
6038  w[108] = 0.00834283875396815770558;
6039  w[109] = 0.00770337523327974184817;
6040  w[110] = 0.00707248999543355546805;
6041  w[111] = 0.00645190005017573692280;
6042  w[112] = 0.00584344987583563950756;
6043  w[113] = 0.00524912345480885912513;
6044  w[114] = 0.00467105037211432174741;
6045  w[115] = 0.00411150397865469304717;
6046  w[116] = 0.00357289278351729964938;
6047  w[117] = 0.00305775341017553113613;
6048  w[118] = 0.00256876494379402037313;
6049  w[119] = 0.00210881524572663287933;
6050  w[120] = 0.00168114286542146990631;
6051  w[121] = 0.00128952408261041739210;
6052  w[122] = 0.000938369848542381500794;
6053  w[123] = 0.000632607319362633544219;
6054  w[124] = 0.000377746646326984660274;
6055  w[125] = 0.000180739564445388357820;
6056  w[126] = 0.0000505360952078625176247;
6057  }
6058  else if (n==255) {
6059  w[ 0] = 0.69379364324108267170E-05;
6060  w[ 1] = 0.25157870384280661489E-04;
6061  w[ 2] = 0.53275293669780613125E-04;
6062  w[ 3] = 0.90372734658751149261E-04;
6063  w[ 4] = 0.13575491094922871973E-03;
6064  w[ 5] = 0.18887326450650491366E-03;
6065  w[ 6] = 0.24921240048299729402E-03;
6066  w[ 7] = 0.31630366082226447689E-03;
6067  w[ 8] = 0.38974528447328229322E-03;
6068  w[ 9] = 0.46918492424785040975E-03;
6069  w[ 10] = 0.55429531493037471492E-03;
6070  w[ 11] = 0.64476204130572477933E-03;
6071  w[ 12] = 0.74028280424450333046E-03;
6072  w[ 13] = 0.84057143271072246365E-03;
6073  w[ 14] = 0.94536151685852538246E-03;
6074  w[ 15] = 0.10544076228633167722E-02;
6075  w[ 16] = 0.11674841174299594077E-02;
6076  w[ 17] = 0.12843824718970101768E-02;
6077  w[ 18] = 0.14049079956551446427E-02;
6078  w[ 19] = 0.15288767050877655684E-02;
6079  w[ 20] = 0.16561127281544526052E-02;
6080  w[ 21] = 0.17864463917586498247E-02;
6081  w[ 22] = 0.19197129710138724125E-02;
6082  w[ 23] = 0.20557519893273465236E-02;
6083  w[ 24] = 0.21944069253638388388E-02;
6084  w[ 25] = 0.23355251860571608737E-02;
6085  w[ 26] = 0.24789582266575679307E-02;
6086  w[ 27] = 0.26245617274044295626E-02;
6087  w[ 28] = 0.27721957645934509940E-02;
6088  w[ 29] = 0.29217249379178197538E-02;
6089  w[ 30] = 0.30730184347025783234E-02;
6090  w[ 31] = 0.32259500250878684614E-02;
6091  w[ 32] = 0.33803979910869203823E-02;
6092  w[ 33] = 0.35362449977167777340E-02;
6093  w[ 34] = 0.36933779170256508183E-02;
6094  w[ 35] = 0.38516876166398709241E-02;
6095  w[ 36] = 0.40110687240750233989E-02;
6096  w[ 37] = 0.41714193769840788528E-02;
6097  w[ 38] = 0.43326409680929828545E-02;
6098  w[ 39] = 0.44946378920320678616E-02;
6099  w[ 40] = 0.46573172997568547773E-02;
6100  w[ 41] = 0.48205888648512683476E-02;
6101  w[ 42] = 0.49843645647655386012E-02;
6102  w[ 43] = 0.51485584789781777618E-02;
6103  w[ 44] = 0.53130866051870565663E-02;
6104  w[ 45] = 0.54778666939189508240E-02;
6105  w[ 46] = 0.56428181013844441585E-02;
6106  w[ 47] = 0.58078616599775673635E-02;
6107  w[ 48] = 0.59729195655081658049E-02;
6108  w[ 49] = 0.61379152800413850435E-02;
6109  w[ 50] = 0.63027734490857587172E-02;
6110  w[ 51] = 0.64674198318036867274E-02;
6111  w[ 52] = 0.66317812429018878941E-02;
6112  w[ 53] = 0.67957855048827733948E-02;
6113  w[ 54] = 0.69593614093904229394E-02;
6114  w[ 55] = 0.71224386864583871532E-02;
6115  w[ 56] = 0.72849479805538070639E-02;
6116  w[ 57] = 0.74468208324075910174E-02;
6117  w[ 58] = 0.76079896657190565832E-02;
6118  w[ 59] = 0.77683877779219912200E-02;
6119  w[ 60] = 0.79279493342948491103E-02;
6120  w[ 61] = 0.80866093647888599710E-02;
6121  w[ 62] = 0.82443037630328680306E-02;
6122  w[ 63] = 0.84009692870519326354E-02;
6123  w[ 64] = 0.85565435613076896192E-02;
6124  w[ 65] = 0.87109650797320868736E-02;
6125  w[ 66] = 0.88641732094824942641E-02;
6126  w[ 67] = 0.90161081951956431600E-02;
6127  w[ 68] = 0.91667111635607884067E-02;
6128  w[ 69] = 0.93159241280693950932E-02;
6129  w[ 70] = 0.94636899938300652943E-02;
6130  w[ 71] = 0.96099525623638830097E-02;
6131  w[ 72] = 0.97546565363174114611E-02;
6132  w[ 73] = 0.98977475240487497440E-02;
6133  w[ 74] = 0.10039172044056840798E-01;
6134  w[ 75] = 0.10178877529236079733E-01;
6135  w[ 76] = 0.10316812330947621682E-01;
6136  w[ 77] = 0.10452925722906011926E-01;
6137  w[ 78] = 0.10587167904885197931E-01;
6138  w[ 79] = 0.10719490006251933623E-01;
6139  w[ 80] = 0.10849844089337314099E-01;
6140  w[ 81] = 0.10978183152658912470E-01;
6141  w[ 82] = 0.11104461134006926537E-01;
6142  w[ 83] = 0.11228632913408049354E-01;
6143  w[ 84] = 0.11350654315980596602E-01;
6144  w[ 85] = 0.11470482114693874380E-01;
6145  w[ 86] = 0.11588074033043952568E-01;
6146  w[ 87] = 0.11703388747657003101E-01;
6147  w[ 88] = 0.11816385890830235763E-01;
6148  w[ 89] = 0.11927026053019270040E-01;
6149  w[ 90] = 0.12035270785279562630E-01;
6150  w[ 91] = 0.12141082601668299679E-01;
6151  w[ 92] = 0.12244424981611985899E-01;
6152  w[ 93] = 0.12345262372243838455E-01;
6153  w[ 94] = 0.12443560190714035263E-01;
6154  w[ 95] = 0.12539284826474884353E-01;
6155  w[ 96] = 0.12632403643542078765E-01;
6156  w[ 97] = 0.12722884982732382906E-01;
6157  w[ 98] = 0.12810698163877361967E-01;
6158  w[ 99] = 0.12895813488012114694E-01;
6159  w[100] = 0.12978202239537399286E-01;
6160  w[101] = 0.13057836688353048840E-01;
6161  w[102] = 0.13134690091960152836E-01;
6162  w[103] = 0.13208736697529129966E-01;
6163  w[104] = 0.13279951743930530650E-01;
6164  w[105] = 0.13348311463725179953E-01;
6165  w[106] = 0.13413793085110098513E-01;
6166  w[107] = 0.13476374833816515982E-01;
6167  w[108] = 0.13536035934956213614E-01;
6168  w[109] = 0.13592756614812395910E-01;
6169  w[110] = 0.13646518102571291428E-01;
6170  w[111] = 0.13697302631990716258E-01;
6171  w[112] = 0.13745093443001896632E-01;
6172  w[113] = 0.13789874783240936517E-01;
6173  w[114] = 0.13831631909506428676E-01;
6174  w[115] = 0.13870351089139840997E-01;
6175  w[116] = 0.13906019601325461264E-01;
6176  w[117] = 0.13938625738306850804E-01;
6177  w[118] = 0.13968158806516938516E-01;
6178  w[119] = 0.13994609127619079852E-01;
6179  w[120] = 0.14017968039456608810E-01;
6180  w[121] = 0.14038227896908623303E-01;
6181  w[122] = 0.14055382072649964277E-01;
6182  w[123] = 0.14069424957813575318E-01;
6183  w[124] = 0.14080351962553661325E-01;
6184  w[125] = 0.14088159516508301065E-01;
6185  w[126] = 0.14092845069160408355E-01;
6186  w[127] = 0.14094407090096179347E-01;
6187  w[128] = 0.14092845069160408355E-01;
6188  w[129] = 0.14088159516508301065E-01;
6189  w[130] = 0.14080351962553661325E-01;
6190  w[131] = 0.14069424957813575318E-01;
6191  w[132] = 0.14055382072649964277E-01;
6192  w[133] = 0.14038227896908623303E-01;
6193  w[134] = 0.14017968039456608810E-01;
6194  w[135] = 0.13994609127619079852E-01;
6195  w[136] = 0.13968158806516938516E-01;
6196  w[137] = 0.13938625738306850804E-01;
6197  w[138] = 0.13906019601325461264E-01;
6198  w[139] = 0.13870351089139840997E-01;
6199  w[140] = 0.13831631909506428676E-01;
6200  w[141] = 0.13789874783240936517E-01;
6201  w[142] = 0.13745093443001896632E-01;
6202  w[143] = 0.13697302631990716258E-01;
6203  w[144] = 0.13646518102571291428E-01;
6204  w[145] = 0.13592756614812395910E-01;
6205  w[146] = 0.13536035934956213614E-01;
6206  w[147] = 0.13476374833816515982E-01;
6207  w[148] = 0.13413793085110098513E-01;
6208  w[149] = 0.13348311463725179953E-01;
6209  w[150] = 0.13279951743930530650E-01;
6210  w[151] = 0.13208736697529129966E-01;
6211  w[152] = 0.13134690091960152836E-01;
6212  w[153] = 0.13057836688353048840E-01;
6213  w[154] = 0.12978202239537399286E-01;
6214  w[155] = 0.12895813488012114694E-01;
6215  w[156] = 0.12810698163877361967E-01;
6216  w[157] = 0.12722884982732382906E-01;
6217  w[158] = 0.12632403643542078765E-01;
6218  w[159] = 0.12539284826474884353E-01;
6219  w[160] = 0.12443560190714035263E-01;
6220  w[161] = 0.12345262372243838455E-01;
6221  w[162] = 0.12244424981611985899E-01;
6222  w[163] = 0.12141082601668299679E-01;
6223  w[164] = 0.12035270785279562630E-01;
6224  w[165] = 0.11927026053019270040E-01;
6225  w[166] = 0.11816385890830235763E-01;
6226  w[167] = 0.11703388747657003101E-01;
6227  w[168] = 0.11588074033043952568E-01;
6228  w[169] = 0.11470482114693874380E-01;
6229  w[170] = 0.11350654315980596602E-01;
6230  w[171] = 0.11228632913408049354E-01;
6231  w[172] = 0.11104461134006926537E-01;
6232  w[173] = 0.10978183152658912470E-01;
6233  w[174] = 0.10849844089337314099E-01;
6234  w[175] = 0.10719490006251933623E-01;
6235  w[176] = 0.10587167904885197931E-01;
6236  w[177] = 0.10452925722906011926E-01;
6237  w[178] = 0.10316812330947621682E-01;
6238  w[179] = 0.10178877529236079733E-01;
6239  w[180] = 0.10039172044056840798E-01;
6240  w[181] = 0.98977475240487497440E-02;
6241  w[182] = 0.97546565363174114611E-02;
6242  w[183] = 0.96099525623638830097E-02;
6243  w[184] = 0.94636899938300652943E-02;
6244  w[185] = 0.93159241280693950932E-02;
6245  w[186] = 0.91667111635607884067E-02;
6246  w[187] = 0.90161081951956431600E-02;
6247  w[188] = 0.88641732094824942641E-02;
6248  w[189] = 0.87109650797320868736E-02;
6249  w[190] = 0.85565435613076896192E-02;
6250  w[191] = 0.84009692870519326354E-02;
6251  w[192] = 0.82443037630328680306E-02;
6252  w[193] = 0.80866093647888599710E-02;
6253  w[194] = 0.79279493342948491103E-02;
6254  w[195] = 0.77683877779219912200E-02;
6255  w[196] = 0.76079896657190565832E-02;
6256  w[197] = 0.74468208324075910174E-02;
6257  w[198] = 0.72849479805538070639E-02;
6258  w[199] = 0.71224386864583871532E-02;
6259  w[200] = 0.69593614093904229394E-02;
6260  w[201] = 0.67957855048827733948E-02;
6261  w[202] = 0.66317812429018878941E-02;
6262  w[203] = 0.64674198318036867274E-02;
6263  w[204] = 0.63027734490857587172E-02;
6264  w[205] = 0.61379152800413850435E-02;
6265  w[206] = 0.59729195655081658049E-02;
6266  w[207] = 0.58078616599775673635E-02;
6267  w[208] = 0.56428181013844441585E-02;
6268  w[209] = 0.54778666939189508240E-02;
6269  w[210] = 0.53130866051870565663E-02;
6270  w[211] = 0.51485584789781777618E-02;
6271  w[212] = 0.49843645647655386012E-02;
6272  w[213] = 0.48205888648512683476E-02;
6273  w[214] = 0.46573172997568547773E-02;
6274  w[215] = 0.44946378920320678616E-02;
6275  w[216] = 0.43326409680929828545E-02;
6276  w[217] = 0.41714193769840788528E-02;
6277  w[218] = 0.40110687240750233989E-02;
6278  w[219] = 0.38516876166398709241E-02;
6279  w[220] = 0.36933779170256508183E-02;
6280  w[221] = 0.35362449977167777340E-02;
6281  w[222] = 0.33803979910869203823E-02;
6282  w[223] = 0.32259500250878684614E-02;
6283  w[224] = 0.30730184347025783234E-02;
6284  w[225] = 0.29217249379178197538E-02;
6285  w[226] = 0.27721957645934509940E-02;
6286  w[227] = 0.26245617274044295626E-02;
6287  w[228] = 0.24789582266575679307E-02;
6288  w[229] = 0.23355251860571608737E-02;
6289  w[230] = 0.21944069253638388388E-02;
6290  w[231] = 0.20557519893273465236E-02;
6291  w[232] = 0.19197129710138724125E-02;
6292  w[233] = 0.17864463917586498247E-02;
6293  w[234] = 0.16561127281544526052E-02;
6294  w[235] = 0.15288767050877655684E-02;
6295  w[236] = 0.14049079956551446427E-02;
6296  w[237] = 0.12843824718970101768E-02;
6297  w[238] = 0.11674841174299594077E-02;
6298  w[239] = 0.10544076228633167722E-02;
6299  w[240] = 0.94536151685852538246E-03;
6300  w[241] = 0.84057143271072246365E-03;
6301  w[242] = 0.74028280424450333046E-03;
6302  w[243] = 0.64476204130572477933E-03;
6303  w[244] = 0.55429531493037471492E-03;
6304  w[245] = 0.46918492424785040975E-03;
6305  w[246] = 0.38974528447328229322E-03;
6306  w[247] = 0.31630366082226447689E-03;
6307  w[248] = 0.24921240048299729402E-03;
6308  w[249] = 0.18887326450650491366E-03;
6309  w[250] = 0.13575491094922871973E-03;
6310  w[251] = 0.90372734658751149261E-04;
6311  w[252] = 0.53275293669780613125E-04;
6312  w[253] = 0.25157870384280661489E-04;
6313  w[254] = 0.69379364324108267170E-05;
6314  }
6315  else {
6316  std::cerr << "\n";
6317  std::cerr << "PATTERSON_LOOKUP_WEIGHTS - Fatal error!\n";
6318  std::cerr << " Unexpected value of N = " << n << ".\n";
6319  std::exit(1);
6320  }
6321  return;
6322 }
6323 
6324 //***************************************************************************
6325 template<class Scalar>
6326 void IntrepidBurkardtRules::trapezoidal_compute ( int n, Scalar x[], Scalar w[] )
6327 //***************************************************************************
6328 {
6329  if (n==1) {
6330  x[0] = 0.0;
6331  w[0] = 2.0;
6332  }
6333  else {
6334  Scalar h = 1.0/((Scalar)n-1.0);
6335  for (int i=0; i<n; i++) {
6336  x[i] = -1.0 + (Scalar)i*h*2.0;
6337  if (i==0||i==n-1) {
6338  w[i] = h;
6339  }
6340  else {
6341  w[i] = 2.0*h;
6342  }
6343  }
6344  }
6345  return;
6346 }
6347 
6348 //***************************************************************************
6349 template<class Scalar>
6351 //***************************************************************************
6352 {
6353  if (n==1) {
6354  x[0] = 0.0;
6355  }
6356  else {
6357  Scalar h = 1.0/((Scalar)n-1.0);
6358  for (int i=0; i<n; i++) {
6359  x[i] = -1.0 + (Scalar)i*h*2.0;
6360  }
6361  }
6362  return;
6363 }
6364 
6365 //***************************************************************************
6366 template<class Scalar>
6368 //***************************************************************************
6369 {
6370  if (n==1) {
6371  w[0] = 2.0;
6372  }
6373  else {
6374  Scalar h = 1.0/((Scalar)n-1.0);
6375  for (int i=0; i<n; i++) {
6376  if (i==0||i==n-1) {
6377  w[i] = h;
6378  }
6379  else {
6380  w[i] = 2.0*h;
6381  }
6382  }
6383  }
6384  return;
6385 }
6386 
6387 //****************************************************************************
6388 template<class Scalar>
6389 Scalar IntrepidBurkardtRules::r8_epsilon(Scalar one)
6390 //****************************************************************************
6391 //
6392 // Purpose:
6393 //
6394 // R8_EPSILON returns the R8 roundoff unit.
6395 //
6396 // Discussion:
6397 //
6398 // The roundoff unit is a number R which is a power of 2 with the
6399 // property that, to the precision of the computer's arithmetic,
6400 // 1 < 1 + R
6401 // but
6402 // 1 = ( 1 + R / 2 )
6403 //
6404 // Licensing:
6405 //
6406 // This code is distributed under the GNU LGPL license.
6407 //
6408 // Modified:
6409 //
6410 // 18 February 2008
6411 //
6412 // Author:
6413 //
6414 // John Burkardt
6415 //
6416 // Parameters:
6417 //
6418 // Output, Scalar R8_EPSILON, the R8 round-off unit.
6419 //
6420 {
6421  Scalar value; value = 1.0;
6422 
6423  while (1.0<(Scalar)(1.0+value)) {
6424  value = value / 2.0;
6425  }
6426 
6427  value = 2.0 * value;
6428 
6429  return value;
6430 }
6431 
6432 //****************************************************************************
6433 template<class Scalar>
6434 Scalar IntrepidBurkardtRules::r8_sign ( Scalar x )
6435 //****************************************************************************
6436 //
6437 // Purpose:
6438 //
6439 // R8_SIGN returns the sign of an R8.
6440 //
6441 // Licensing:
6442 //
6443 // This code is distributed under the GNU LGPL license.
6444 //
6445 // Modified:
6446 //
6447 // 18 October 2004
6448 //
6449 // Author:
6450 //
6451 // John Burkardt
6452 //
6453 // Parameters:
6454 //
6455 // Input, Scalar X, the number whose sign is desired.
6456 //
6457 // Output, Scalar R8_SIGN, the sign of X.
6458 //
6459 {
6460  Scalar value;
6461 
6462  if (x<0.0) {
6463  value = -1.0;
6464  }
6465  else {
6466  value = 1.0;
6467  }
6468  return value;
6469 }
6470 
6471 } // end of namespace Intrepid
6472 
static void chebyshev2_compute_weights(int order, Scalar w[])
Gauss-Chebyshev of Type 2; returns weights.
static void patterson_lookup_weights(int n, Scalar w[])
Gauss-Patterson; returns weights.
static void chebyshev1_compute(int order, Scalar x[], Scalar w[])
Gauss-Chebyshev of Type 1; returns points and weights.
static void laguerre_lookup_weights(int n, Scalar w[])
Gauss-Laguerre; returns weights.
static void chebyshev1_compute_weights(int order, Scalar w[])
Gauss-Chebyshev of Type 1; returns weights.
static void hermite_genz_keister_lookup_weights(int n, Scalar w[])
Hermite-Genz-Keister; returns weights.
static void clenshaw_curtis_compute_points(int order, Scalar x[])
Clenshaw-Curtis; returns points.
static void legendre_compute_weights(int order, Scalar w[])
Gauss-Legendre; returns weights.
static void hermite_lookup(int n, Scalar x[], Scalar w[])
Gauss-Hermite; returns points and weights.
static void laguerre_compute_points(int order, Scalar x[])
Gauss-Laguerre; returns points.
static void chebyshev2_compute_points(int order, Scalar x[])
Gauss-Chebyshev of Type 2; returns points.
static void patterson_lookup_points(int n, Scalar x[])
Gauss-Patterson; returns points.
static void fejer2_compute_weights(int order, Scalar w[])
Fejer type 2; returns weights.
static void hermite_genz_keister_lookup_points(int n, Scalar x[])
Hermite-Genz-Keister; returns points.
static void trapezoidal_compute(int n, Scalar x[], Scalar w[])
Trapezoidal rule; returns points and weights.
static void trapezoidal_compute_points(int order, Scalar x[])
Trapezoidal rule; returns points.
static void trapezoidal_compute_weights(int order, Scalar w[])
Trapezoidal rule; returns weights.
static void legendre_lookup_weights(int n, Scalar w[])
Gauss-Legendre; returns weights.
static void legendre_lookup_points(int n, Scalar x[])
Gauss-Legendre; returns points.
static void patterson_lookup(int n, Scalar x[], Scalar w[])
Gauss-Patterson; returns points and weights.
static void hermite_compute_weights(int order, Scalar w[])
Gauss-Hermite; returns weights.
static void fejer2_compute_points(int order, Scalar x[])
Fejer type 2; returns points.
static void hermite_compute_points(int order, Scalar x[])
Gauss-Hermite; returns points.
static void laguerre_compute(int n, Scalar x[], Scalar w[])
Gauss-Laguerre; returns points and weights.
static void hermite_lookup_weights(int n, Scalar w[])
Gauss-Hermite; returns weights.
static void legendre_compute_points(int order, Scalar x[])
Gauss-Legendre; returns points.
static void clenshaw_curtis_compute_weights(int order, Scalar w[])
Clenshaw-Curtis; returns weights.
static void hermite_compute(int order, Scalar x[], Scalar w[])
Gauss-Hermite; returns points and weights.
static void chebyshev2_compute(int order, Scalar x[], Scalar w[])
Gauss-Chebyshev of Type 2; returns points and weights.
static void chebyshev1_compute_points(int order, Scalar x[])
Gauss-Chebyshev of Type 1; returns points.
static void legendre_lookup(int n, Scalar x[], Scalar w[])
Gauss-Legendre; returns points and weights.
static void hermite_lookup_points(int n, Scalar x[])
Gauss-Hermite; returns points.
static void laguerre_compute_weights(int order, Scalar w[])
Gauss-Laguerre; returns weights.
static void laguerre_lookup(int n, Scalar x[], Scalar w[])
Gauss-Laguerre; returns points and weights.
static void legendre_compute(int n, Scalar x[], Scalar w[])
Gauss-Legendre; returns points and weights.
static void hermite_genz_keister_lookup(int n, Scalar x[], Scalar w[])
Hermite-Genz-Keister; returns points and weights.
static void laguerre_lookup_points(int n, Scalar x[])
Gauss-Laguerre; returns points.
static void clenshaw_curtis_compute(int order, Scalar x[], Scalar w[])
Clenshaw-Curtis; returns points and weights.
static void fejer2_compute(int order, Scalar x[], Scalar w[])
Fejer type 2; returns points and weights.