-
Notifications
You must be signed in to change notification settings - Fork 8
/
Copy pathsparse_grid_gl_dataset.cpp
2915 lines (2728 loc) · 78.1 KB
/
sparse_grid_gl_dataset.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
# include <cstdlib>
# include <cmath>
# include <iostream>
# include <fstream>
# include <iomanip>
# include <ctime>
using namespace std;
int main ( int argc, char *argv[] );
int choose ( int n, int k );
void comp_next ( int n, int k, int a[], bool *more, int *h, int *t );
void gl_abscissa ( int dim_num, int point_num, int grid_index[],
int grid_base[], double grid_point[] );
void gl_weights ( int order, double weight[] );
int i4_log_2 ( int i );
int i4_max ( int i1, int i2 );
int i4_min ( int i1, int i2 );
int i4_modp ( int i, int j );
int i4_power ( int i, int j );
string i4_to_string ( int i4, string format );
int i4vec_product ( int n, int a[] );
int *index_level_gl ( int level, int level_max, int dim_num, int point_num,
int grid_index[], int grid_base[] );
void level_to_order_open ( int dim_num, int level[], int order[] );
int *multigrid_index_z ( int dim_num, int order_1d[], int order_nd );
double *product_weight_gl ( int dim_num, int order_1d[], int order_nd );
double r8_epsilon ( );
double r8_huge ( );
void r8mat_transpose_print_some ( int m, int n, double a[], int ilo, int jlo,
int ihi, int jhi, string title );
void r8mat_write ( string output_filename, int m, int n, double table[] );
void r8vec_direct_product2 ( int factor_index, int factor_order,
double factor_value[], int factor_num, int point_num, double w[] );
void r8vec_print_some ( int n, double a[], int i_lo, int i_hi, string title );
int s_len_trim ( string s );
void sparse_grid_gl ( int dim_num, int level_max, int point_num,
double grid_weight[], double grid_point[] );
void sparse_grid_gl_index ( int dim_num, int level_max, int point_num,
int grid_index [], int grid_base[] );
int sparse_grid_gl_size ( int dim_num, int level_max );
void timestamp ( );
void vec_colex_next2 ( int dim_num, int base[], int a[], bool *more );
//****************************************************************************80
int main ( int argc, char *argv[] )
//****************************************************************************80
//
// Purpose:
//
// MAIN is the main program for SPARSE_GRID_GL_DATASET.
//
// Discussion:
//
// This program computes a sparse grid quadrature rule based on a 1D
// Gauss-Legendre rule and writes it to a file..
//
// The user specifies:
// * the spatial dimension of the quadrature region,
// * the level that defines the Smolyak grid.
//
// Licensing:
//
// This code is distributed under the GNU LGPL license.
//
// Modified:
//
// 09 July 2009
//
// Author:
//
// John Burkardt
//
// Reference:
//
// Fabio Nobile, Raul Tempone, Clayton Webster,
// A Sparse Grid Stochastic Collocation Method for Partial Differential
// Equations with Random Input Data,
// SIAM Journal on Numerical Analysis,
// Volume 46, Number 5, 2008, pages 2309-2345.
//
{
int dim;
int dim_num;
int level_max;
int level_min;
int point;
int point_num;
double *r;
string r_filename;
double *w;
string w_filename;
double weight_sum;
double *x;
string x_filename;
timestamp ( );
cout << "\n";
cout << "SPARSE_GRID_GL_DATASET\n";
cout << " C++ version\n";
cout << "\n";
cout << " Compiled on " << __DATE__ << " at " << __TIME__ << ".\n";
cout << "\n";
cout << " Compute the abscissas and weights of a quadrature rule\n";
cout << " associated with a sparse grid derived from a Smolyak\n";
cout << " construction based on a 1D Gauss-Legendre rule.\n";
cout << "\n";
cout << " Inputs to the program include:\n";
cout << "\n";
cout << " DIM_NUM, the spatial dimension.\n";
cout << " (typically in the range of 2 to 10)\n";
cout << "\n";
cout << " LEVEL_MAX, the level of the sparse grid.\n";
cout << " (typically in the range of 0, 1, 2, 3, ...\n";
cout << "\n";
cout << " Output from the program includes:\n";
cout << "\n";
cout << " * A printed table of the abscissas and weights.\n";
cout << "\n";
cout << " * A set of 3 files that define the quadrature rule.\n";
cout << "\n";
cout << " (1) gl_d?_level?_r.txt, the ranges;\n";
cout << " (2) gl_d?_level?_w.txt, the weights;\n";
cout << " (3) gl_d?_level?_x.txt, the abscissas.\n";
//
// Get the spatial dimension.
//
if ( 1 < argc )
{
dim_num = atoi ( argv[1] );
}
else
{
cout << "\n";
cout << " Enter the value of DIM_NUM (1 or greater)\n";
cin >> dim_num;
}
cout << "\n";
cout << " Spatial dimension requested is = " << dim_num << "\n";
//
// Get the level.
//
if ( 2 < argc )
{
level_max = atoi ( argv[2] );
}
else
{
cout << "\n";
cout << " Enter the value of LEVEL_MAX (0 or greater).\n";
cin >> level_max;
}
level_min = i4_max ( 0, level_max + 1 - dim_num );
cout << "\n";
cout << " LEVEL_MIN is = " << level_min << "\n";
cout << " LEVEL_MAX is = " << level_max << "\n";
//
// How many distinct points will there be?
//
point_num = sparse_grid_gl_size ( dim_num, level_max );
cout << "\n";
cout << " The number of distinct abscissas in the\n";
cout << " quadrature rule is determined from the spatial\n";
cout << " dimension DIM_NUM and the level LEVEL_MAX.\n";
cout << " For the given input, this value will be = " << point_num << "\n";
//
// Allocate memory.
//
r = new double[dim_num*2];
w = new double[point_num];
x = new double[dim_num*point_num];
//
// Compute the weights and points.
//
for ( dim = 0; dim < dim_num; dim++ )
{
r[dim+0*dim_num] = -1.0;
r[dim+1*dim_num] = +1.0;
}
sparse_grid_gl ( dim_num, level_max, point_num, w, x );
r8mat_transpose_print_some ( dim_num, point_num, x, 1, 1, dim_num,
10, " First 10 grid points:" );
r8vec_print_some ( point_num, w, 1, 10, " First 10 grid weights:" );
weight_sum = 0.0;
for ( point = 0; point < point_num; point++ )
{
weight_sum = weight_sum + w[point];
}
cout << "\n";
cout << " Weights sum to " << weight_sum << "\n";
cout << " Correct value is " << pow ( 2.0, dim_num ) << "\n";
//
// Construct appropriate file names.
//
r_filename = "gl_d" + i4_to_string ( dim_num, "%d" )
+ "_level" + i4_to_string ( level_max, "%d" ) + "_r.txt";
w_filename = "gl_d" + i4_to_string ( dim_num, "%d" )
+ "_level" + i4_to_string ( level_max, "%d" ) + "_w.txt";
x_filename = "gl_d" + i4_to_string ( dim_num, "%d" )
+ "_level" + i4_to_string ( level_max, "%d" ) + "_x.txt";
//
// Write the rule to files.
//
cout << "\n";
cout << " Creating R file = \"" << r_filename << "\".\n";
r8mat_write ( r_filename, dim_num, 2, r );
cout << " Creating W file = \"" << w_filename << "\".\n";
r8mat_write ( w_filename, 1, point_num, w );
cout << " Creating X file = \"" << x_filename << "\".\n";
r8mat_write ( x_filename, dim_num, point_num, x );
delete [] r;
delete [] w;
delete [] x;
cout << "\n";
cout << "SPARSE_GRID_GL_DATASET\n";
cout << " Normal end of execution.\n";
cout << "\n";
timestamp ( );
return 0;
}
//****************************************************************************80
int choose ( int n, int k )
//****************************************************************************80
//
// Purpose:
//
// CHOOSE computes the binomial coefficient C(N,K).
//
// Discussion:
//
// The value is calculated in such a way as to avoid overflow and
// roundoff. The calculation is done in integer arithmetic.
//
// The formula used is:
//
// C(N,K) = N! / ( K! * (N-K)! )
//
// Licensing:
//
// This code is distributed under the GNU LGPL license.
//
// Modified:
//
// 22 May 2007
//
// Author:
//
// John Burkardt
//
// Reference:
//
// ML Wolfson, HV Wright,
// Algorithm 160:
// Combinatorial of M Things Taken N at a Time,
// Communications of the ACM,
// Volume 6, Number 4, April 1963, page 161.
//
// Parameters:
//
// Input, int N, K, are the values of N and K.
//
// Output, int CHOOSE, the number of combinations of N
// things taken K at a time.
//
{
int i;
int mn;
int mx;
int value;
mn = i4_min ( k, n - k );
if ( mn < 0 )
{
value = 0;
}
else if ( mn == 0 )
{
value = 1;
}
else
{
mx = i4_max ( k, n - k );
value = mx + 1;
for ( i = 2; i <= mn; i++ )
{
value = ( value * ( mx + i ) ) / i;
}
}
return value;
}
//****************************************************************************80
void comp_next ( int n, int k, int a[], bool *more, int *h, int *t )
//****************************************************************************80
//
// Purpose:
//
// COMP_NEXT computes the compositions of the integer N into K parts.
//
// Discussion:
//
// A composition of the integer N into K parts is an ordered sequence
// of K nonnegative integers which sum to N. The compositions (1,2,1)
// and (1,1,2) are considered to be distinct.
//
// The routine computes one composition on each call until there are no more.
// For instance, one composition of 6 into 3 parts is
// 3+2+1, another would be 6+0+0.
//
// On the first call to this routine, set MORE = FALSE. The routine
// will compute the first element in the sequence of compositions, and
// return it, as well as setting MORE = TRUE. If more compositions
// are desired, call again, and again. Each time, the routine will
// return with a new composition.
//
// However, when the LAST composition in the sequence is computed
// and returned, the routine will reset MORE to FALSE, signaling that
// the end of the sequence has been reached.
//
// This routine originally used a SAVE statement to maintain the
// variables H and T. I have decided that it is safer
// to pass these variables as arguments, even though the user should
// never alter them. This allows this routine to safely shuffle
// between several ongoing calculations.
//
//
// There are 28 compositions of 6 into three parts. This routine will
// produce those compositions in the following order:
//
// I A
// - ---------
// 1 6 0 0
// 2 5 1 0
// 3 4 2 0
// 4 3 3 0
// 5 2 4 0
// 6 1 5 0
// 7 0 6 0
// 8 5 0 1
// 9 4 1 1
// 10 3 2 1
// 11 2 3 1
// 12 1 4 1
// 13 0 5 1
// 14 4 0 2
// 15 3 1 2
// 16 2 2 2
// 17 1 3 2
// 18 0 4 2
// 19 3 0 3
// 20 2 1 3
// 21 1 2 3
// 22 0 3 3
// 23 2 0 4
// 24 1 1 4
// 25 0 2 4
// 26 1 0 5
// 27 0 1 5
// 28 0 0 6
//
// Licensing:
//
// This code is distributed under the GNU LGPL license.
//
// Modified:
//
// 02 July 2008
//
// Author:
//
// Original FORTRAN77 version by Albert Nijenhuis, Herbert Wilf.
// C++ version by John Burkardt.
//
// Reference:
//
// Albert Nijenhuis, Herbert Wilf,
// Combinatorial Algorithms for Computers and Calculators,
// Second Edition,
// Academic Press, 1978,
// ISBN: 0-12-519260-6,
// LC: QA164.N54.
//
// Parameters:
//
// Input, int N, the integer whose compositions are desired.
//
// Input, int K, the number of parts in the composition.
//
// Input/output, int A[K], the parts of the composition.
//
// Input/output, bool *MORE.
// Set MORE = FALSE on first call. It will be reset to TRUE on return
// with a new composition. Each new call returns another composition until
// MORE is set to FALSE when the last composition has been computed
// and returned.
//
// Input/output, int *H, *T, two internal parameters needed for the
// computation. The user should allocate space for these in the calling
// program, include them in the calling sequence, but never alter them!
//
{
int i;
if ( !( *more ) )
{
*t = n;
*h = 0;
a[0] = n;
for ( i = 1; i < k; i++ )
{
a[i] = 0;
}
}
else
{
if ( 1 < *t )
{
*h = 0;
}
*h = *h + 1;
*t = a[*h-1];
a[*h-1] = 0;
a[0] = *t - 1;
a[*h] = a[*h] + 1;
}
*more = ( a[k-1] != n );
return;
}
//****************************************************************************80
void gl_abscissa ( int dim_num, int point_num, int grid_index[],
int grid_base[], double grid_point[] )
//****************************************************************************80
//
// Purpose:
//
// GL_ABSCISSA sets abscissas for multidimensional Gauss-Legendre quadrature.
//
// Discussion:
//
// The "nesting" as it occurs for Gauss-Legendre sparse grids simply
// involves the use of a specified set of permissible orders for the
// rule.
//
// The X array lists the (complete) Gauss-Legendre abscissas for rules
// of order 1, 3, 7, 15, 31, 63 or 127, in order.
//
// Licensing:
//
// This code is distributed under the GNU LGPL license.
//
// Modified:
//
// 02 October 2007
//
// Author:
//
// John Burkardt
//
// Parameters:
//
// Input, int DIM_NUM, the spatial dimension.
//
// Input, int POINT_NUM, the number of points.
//
// Input, int GRID_INDEX[DIM_NUM*POINT_NUM], the index of the abscissa
// from the Gauss-Legendre rule, for each dimension and point.
//
// Input, int GRID_BASE[DIM_NUM], the number of points used in the
// Gauss-Legendre rule for a given dimension.
//
// Output, double GRID_POINT[DIM_NUM], the grid points of
// Gauss-Legendre abscissas.
//
{
int dim;
int level;
int point;
int pointer;
int skip[8] = { 0, 1, 4, 11, 26, 57, 120, 247 };
double x[247] = {
0.0E+00,
- 0.774596669241483377035853079956E+00,
0.0E+00,
0.774596669241483377035853079956E+00,
- 0.949107912342758524526189684048E+00,
- 0.741531185599394439863864773281E+00,
- 0.405845151377397166906606412077E+00,
0.0E+00,
0.405845151377397166906606412077E+00,
0.741531185599394439863864773281E+00,
0.949107912342758524526189684048E+00,
- 0.987992518020485428489565718587E+00,
- 0.937273392400705904307758947710E+00,
- 0.848206583410427216200648320774E+00,
- 0.724417731360170047416186054614E+00,
- 0.570972172608538847537226737254E+00,
- 0.394151347077563369897207370981E+00,
- 0.201194093997434522300628303395E+00,
0.0E+00,
0.201194093997434522300628303395E+00,
0.394151347077563369897207370981E+00,
0.570972172608538847537226737254E+00,
0.724417731360170047416186054614E+00,
0.848206583410427216200648320774E+00,
0.937273392400705904307758947710E+00,
0.987992518020485428489565718587E+00,
-0.99708748181947707454263838179654,
-0.98468590966515248400211329970113,
-0.96250392509294966178905249675943,
-0.93075699789664816495694576311725,
-0.88976002994827104337419200908023,
-0.83992032014626734008690453594388,
-0.78173314841662494040636002019484,
-0.71577678458685328390597086536649,
-0.64270672292426034618441820323250,
-0.56324916140714926272094492359516,
-0.47819378204490248044059403935649,
-0.38838590160823294306135146128752,
-0.29471806998170161661790389767170,
-0.19812119933557062877241299603283,
-0.99555312152341520325174790118941E-01,
0.00000000000000000000000000000000,
0.99555312152341520325174790118941E-01,
0.19812119933557062877241299603283,
0.29471806998170161661790389767170,
0.38838590160823294306135146128752,
0.47819378204490248044059403935649,
0.56324916140714926272094492359516,
0.64270672292426034618441820323250,
0.71577678458685328390597086536649,
0.78173314841662494040636002019484,
0.83992032014626734008690453594388,
0.88976002994827104337419200908023,
0.93075699789664816495694576311725,
0.96250392509294966178905249675943,
0.98468590966515248400211329970113,
0.99708748181947707454263838179654,
-0.99928298402912378050701628988630E+00,
-0.99622401277797010860209018267357E+00,
-0.99072854689218946681089469460884E+00,
-0.98280881059372723486251140727639E+00,
-0.97248403469757002280196067864927E+00,
-0.95977944975894192707035416626398E+00,
-0.94472613404100980296637531962798E+00,
-0.92736092062184320544703138132518E+00,
-0.90772630277853155803695313291596E+00,
-0.88587032850785342629029845731337E+00,
-0.86184648236412371953961183943106E+00,
-0.83571355431950284347180776961571E+00,
-0.80753549577345676005146598636324E+00,
-0.77738126299037233556333018991104E+00,
-0.74532464831784741782932166103759E+00,
-0.71144409958484580785143153770401E+00,
-0.67582252811498609013110331596954E+00,
-0.63854710582136538500030695387338E+00,
-0.59970905187762523573900892686880E+00,
-0.55940340948628501326769780007005E+00,
-0.51772881329003324812447758452632E+00,
-0.47478724799480439992221230985149E+00,
-0.43068379879511160066208893391863E+00,
-0.38552639421224789247761502227440E+00,
-0.33942554197458440246883443159432E+00,
-0.29249405858625144003615715555067E+00,
-0.24484679324595336274840459392483E+00,
-0.19660034679150668455762745706572E+00,
-0.14787278635787196856983909655297E+00,
-0.98783356446945279529703669453922E-01,
-0.49452187116159627234233818051808E-01,
0.00000000000000000000000000000000E+00,
0.49452187116159627234233818051808E-01,
0.98783356446945279529703669453922E-01,
0.14787278635787196856983909655297E+00,
0.19660034679150668455762745706572E+00,
0.24484679324595336274840459392483E+00,
0.29249405858625144003615715555067E+00,
0.33942554197458440246883443159432E+00,
0.38552639421224789247761502227440E+00,
0.43068379879511160066208893391863E+00,
0.47478724799480439992221230985149E+00,
0.51772881329003324812447758452632E+00,
0.55940340948628501326769780007005E+00,
0.59970905187762523573900892686880E+00,
0.63854710582136538500030695387338E+00,
0.67582252811498609013110331596954E+00,
0.71144409958484580785143153770401E+00,
0.74532464831784741782932166103759E+00,
0.77738126299037233556333018991104E+00,
0.80753549577345676005146598636324E+00,
0.83571355431950284347180776961571E+00,
0.86184648236412371953961183943106E+00,
0.88587032850785342629029845731337E+00,
0.90772630277853155803695313291596E+00,
0.92736092062184320544703138132518E+00,
0.94472613404100980296637531962798E+00,
0.95977944975894192707035416626398E+00,
0.97248403469757002280196067864927E+00,
0.98280881059372723486251140727639E+00,
0.99072854689218946681089469460884E+00,
0.99622401277797010860209018267357E+00,
0.99928298402912378050701628988630E+00
};
for ( dim = 0; dim < dim_num; dim++ )
{
if ( grid_base[dim] < 0 )
{
cout << "\n";
cout << "GL_ABSCISSA - Fatal error!\n";
cout << " Some base values are less than 0.\n";
exit ( 1 );
}
}
for ( dim = 0; dim < dim_num; dim++ )
{
if ( 63 < grid_base[dim] )
{
cout << "\n";
cout << "GL_ABSCISSA - Fatal error!\n";
cout << " Some base values are greater than 63.\n";
exit ( 1 );
}
}
for ( point = 0; point < point_num; point++ )
{
for ( dim = 0; dim < dim_num; dim++ )
{
level = i4_log_2 ( grid_base[dim] + 1 );
pointer = skip[level] + ( grid_index[dim+point*dim_num] + grid_base[dim] );
grid_point[dim+point*dim_num] = x[pointer];
}
}
return;
}
//****************************************************************************80
void gl_weights ( int order, double weight[] )
//****************************************************************************80
//
// Purpose:
//
// GL_WEIGHTS returns weights for certain Gauss-Legendre quadrature rules.
//
// Discussion:
//
// The allowed orders are 1, 3, 7, 15, 31, 63 and 127.
//
// Licensing:
//
// This code is distributed under the GNU LGPL license.
//
// Modified:
//
// 02 October 2007
//
// Author:
//
// John Burkardt
//
// Reference:
//
// Milton Abramowitz, Irene Stegun,
// Handbook of Mathematical Functions,
// National Bureau of Standards, 1964,
// ISBN: 0-486-61272-4,
// LC: QA47.A34.
//
// Arthur Stroud, Don Secrest,
// Gaussian Quadrature Formulas,
// Prentice Hall, 1966,
// LC: QA299.4G3S7.
//
// Parameters:
//
// Input, int ORDER, the order of the rule.
// ORDER must be 1, 3, 7, 15, 31, 63 or 127.
//
// Output, double WEIGHT[ORDER], the weights.
// The weights are positive, symmetric and should sum to 2.
//
{
if ( order == 1 )
{
weight[1-1] = 2.0E+00;
}
else if ( order == 3 )
{
weight[1-1] = 5.0E+00 / 9.0E+00;
weight[2-1] = 8.0E+00 / 9.0E+00;
weight[3-1] = 5.0E+00 / 9.0E+00;
}
else if ( order == 7 )
{
weight[1-1] = 0.129484966168869693270611432679E+00;
weight[2-1] = 0.279705391489276667901467771424E+00;
weight[3-1] = 0.381830050505118944950369775489E+00;
weight[4-1] = 0.417959183673469387755102040816E+00;
weight[5-1] = 0.381830050505118944950369775489E+00;
weight[6-1] = 0.279705391489276667901467771424E+00;
weight[7-1] = 0.129484966168869693270611432679E+00;
}
else if ( order == 15 )
{
weight[1-1] = 0.307532419961172683546283935772E-01;
weight[2-1] = 0.703660474881081247092674164507E-01;
weight[3-1] = 0.107159220467171935011869546686E+00;
weight[4-1] = 0.139570677926154314447804794511E+00;
weight[5-1] = 0.166269205816993933553200860481E+00;
weight[6-1] = 0.186161000015562211026800561866E+00;
weight[7-1] = 0.198431485327111576456118326444E+00;
weight[8-1] = 0.202578241925561272880620199968E+00;
weight[9-1] = 0.198431485327111576456118326444E+00;
weight[10-1] = 0.186161000015562211026800561866E+00;
weight[11-1] = 0.166269205816993933553200860481E+00;
weight[12-1] = 0.139570677926154314447804794511E+00;
weight[13-1] = 0.107159220467171935011869546686E+00;
weight[14-1] = 0.703660474881081247092674164507E-01;
weight[15-1] = 0.307532419961172683546283935772E-01;
}
else if ( order == 31 )
{
weight[ 1-1] = 0.74708315792487746093913218970494E-02;
weight[ 2-1] = 0.17318620790310582463552990782414E-01;
weight[ 3-1] = 0.27009019184979421800608642617676E-01;
weight[ 4-1] = 0.36432273912385464024392008749009E-01;
weight[ 5-1] = 0.45493707527201102902315857856518E-01;
weight[ 6-1] = 0.54103082424916853711666259085477E-01;
weight[ 7-1] = 0.62174786561028426910343543686657E-01;
weight[ 8-1] = 0.69628583235410366167756126255124E-01;
weight[ 9-1] = 0.76390386598776616426357674901331E-01;
weight[10-1] = 0.82392991761589263903823367431962E-01;
weight[11-1] = 0.87576740608477876126198069695333E-01;
weight[12-1] = 0.91890113893641478215362871607150E-01;
weight[13-1] = 0.95290242912319512807204197487597E-01;
weight[14-1] = 0.97743335386328725093474010978997E-01;
weight[15-1] = 0.99225011226672307874875514428615E-01;
weight[16-1] = 0.99720544793426451427533833734349E-01;
weight[17-1] = 0.99225011226672307874875514428615E-01;
weight[18-1] = 0.97743335386328725093474010978997E-01;
weight[19-1] = 0.95290242912319512807204197487597E-01;
weight[20-1] = 0.91890113893641478215362871607150E-01;
weight[21-1] = 0.87576740608477876126198069695333E-01;
weight[22-1] = 0.82392991761589263903823367431962E-01;
weight[23-1] = 0.76390386598776616426357674901331E-01;
weight[24-1] = 0.69628583235410366167756126255124E-01;
weight[25-1] = 0.62174786561028426910343543686657E-01;
weight[26-1] = 0.54103082424916853711666259085477E-01;
weight[27-1] = 0.45493707527201102902315857856518E-01;
weight[28-1] = 0.36432273912385464024392008749009E-01;
weight[29-1] = 0.27009019184979421800608642617676E-01;
weight[30-1] = 0.17318620790310582463552990782414E-01;
weight[31-1] = 0.74708315792487746093913218970494E-02;
}
else if ( order == 63 )
{
weight[ 1-1] = 0.18398745955770837880499331680577E-02;
weight[ 2-1] = 0.42785083468637618661951422543371E-02;
weight[ 3-1] = 0.67102917659601362519069109850892E-02;
weight[ 4-1] = 0.91259686763266563540586445877022E-02;
weight[ 5-1] = 0.11519376076880041750750606118707E-01;
weight[ 6-1] = 0.13884612616115610824866086365937E-01;
weight[ 7-1] = 0.16215878410338338882283672974995E-01;
weight[ 8-1] = 0.18507464460161270409260545805144E-01;
weight[ 9-1] = 0.20753761258039090775341953421471E-01;
weight[10-1] = 0.22949271004889933148942319561770E-01;
weight[11-1] = 0.25088620553344986618630138068443E-01;
weight[12-1] = 0.27166574359097933225189839439413E-01;
weight[13-1] = 0.29178047208280526945551502154029E-01;
weight[14-1] = 0.31118116622219817508215988557189E-01;
weight[15-1] = 0.32982034883779341765683179672459E-01;
weight[16-1] = 0.34765240645355877697180504642788E-01;
weight[17-1] = 0.36463370085457289630452409787542E-01;
weight[18-1] = 0.38072267584349556763638324927889E-01;
weight[19-1] = 0.39587995891544093984807928149202E-01;
weight[20-1] = 0.41006845759666398635110037009072E-01;
weight[21-1] = 0.42325345020815822982505485403028E-01;
weight[22-1] = 0.43540267083027590798964315704401E-01;
weight[23-1] = 0.44648638825941395370332669516813E-01;
weight[24-1] = 0.45647747876292608685885992608542E-01;
weight[25-1] = 0.46535149245383696510395418746953E-01;
weight[26-1] = 0.47308671312268919080604988338844E-01;
weight[27-1] = 0.47966421137995131411052756195132E-01;
weight[28-1] = 0.48506789097883847864090099145802E-01;
weight[29-1] = 0.48928452820511989944709361549215E-01;
weight[30-1] = 0.49230380423747560785043116988145E-01;
weight[31-1] = 0.49411833039918178967039646116705E-01;
weight[32-1] = 0.49472366623931020888669360420926E-01;
weight[33-1] = 0.49411833039918178967039646116705E-01;
weight[34-1] = 0.49230380423747560785043116988145E-01;
weight[35-1] = 0.48928452820511989944709361549215E-01;
weight[36-1] = 0.48506789097883847864090099145802E-01;
weight[37-1] = 0.47966421137995131411052756195132E-01;
weight[38-1] = 0.47308671312268919080604988338844E-01;
weight[39-1] = 0.46535149245383696510395418746953E-01;
weight[40-1] = 0.45647747876292608685885992608542E-01;
weight[41-1] = 0.44648638825941395370332669516813E-01;
weight[42-1] = 0.43540267083027590798964315704401E-01;
weight[43-1] = 0.42325345020815822982505485403028E-01;
weight[44-1] = 0.41006845759666398635110037009072E-01;
weight[45-1] = 0.39587995891544093984807928149202E-01;
weight[46-1] = 0.38072267584349556763638324927889E-01;
weight[47-1] = 0.36463370085457289630452409787542E-01;
weight[48-1] = 0.34765240645355877697180504642788E-01;
weight[49-1] = 0.32982034883779341765683179672459E-01;
weight[50-1] = 0.31118116622219817508215988557189E-01;
weight[51-1] = 0.29178047208280526945551502154029E-01;
weight[52-1] = 0.27166574359097933225189839439413E-01;
weight[53-1] = 0.25088620553344986618630138068443E-01;
weight[54-1] = 0.22949271004889933148942319561770E-01;
weight[55-1] = 0.20753761258039090775341953421471E-01;
weight[56-1] = 0.18507464460161270409260545805144E-01;
weight[57-1] = 0.16215878410338338882283672974995E-01;
weight[58-1] = 0.13884612616115610824866086365937E-01;
weight[59-1] = 0.11519376076880041750750606118707E-01;
weight[60-1] = 0.91259686763266563540586445877022E-02;
weight[61-1] = 0.67102917659601362519069109850892E-02;
weight[62-1] = 0.42785083468637618661951422543371E-02;
weight[63-1] = 0.18398745955770837880499331680577E-02;
}
else if ( order == 127 )
{
weight[ 1-1] = 0.45645726109586654495731936146574E-03;
weight[ 2-1] = 0.10622766869538486959954760554099E-02;
weight[ 3-1] = 0.16683488125171936761028811985672E-02;
weight[ 4-1] = 0.22734860707492547802810838362671E-02;
weight[ 5-1] = 0.28772587656289004082883197417581E-02;
weight[ 6-1] = 0.34792893810051465908910894094105E-02;
weight[ 7-1] = 0.40792095178254605327114733456293E-02;
weight[ 8-1] = 0.46766539777779034772638165662478E-02;
weight[ 9-1] = 0.52712596565634400891303815906251E-02;
weight[ 10-1] = 0.58626653903523901033648343751367E-02;
weight[ 11-1] = 0.64505120486899171845442463868748E-02;
weight[ 12-1] = 0.70344427036681608755685893032552E-02;
weight[ 13-1] = 0.76141028256526859356393930849227E-02;
weight[ 14-1] = 0.81891404887415730817235884718726E-02;
weight[ 15-1] = 0.87592065795403145773316804234385E-02;
weight[ 16-1] = 0.93239550065309714787536985834029E-02;
weight[ 17-1] = 0.98830429087554914716648010899606E-02;
weight[ 18-1] = 0.10436130863141005225673171997668E-01;
weight[ 19-1] = 0.10982883090068975788799657376065E-01;
weight[ 20-1] = 0.11522967656921087154811609734510E-01;
weight[ 21-1] = 0.12056056679400848183529562144697E-01;
weight[ 22-1] = 0.12581826520465013101514365424172E-01;
weight[ 23-1] = 0.13099957986718627426172681912499E-01;
weight[ 24-1] = 0.13610136522139249906034237533759E-01;
weight[ 25-1] = 0.14112052399003395774044161633613E-01;
weight[ 26-1] = 0.14605400905893418351737288078952E-01;
weight[ 27-1] = 0.15089882532666922992635733981431E-01;
weight[ 28-1] = 0.15565203152273955098532590262975E-01;
weight[ 29-1] = 0.16031074199309941802254151842763E-01;
weight[ 30-1] = 0.16487212845194879399346060358146E-01;
weight[ 31-1] = 0.16933342169871654545878815295200E-01;
weight[ 32-1] = 0.17369191329918731922164721250350E-01;
weight[ 33-1] = 0.17794495722974774231027912900351E-01;
weight[ 34-1] = 0.18208997148375106468721469154479E-01;
weight[ 35-1] = 0.18612443963902310429440419898958E-01;
weight[ 36-1] = 0.19004591238555646611148901044533E-01;
weight[ 37-1] = 0.19385200901246454628112623489471E-01;
weight[ 38-1] = 0.19754041885329183081815217323169E-01;
weight[ 39-1] = 0.20110890268880247225644623956287E-01;
weight[ 40-1] = 0.20455529410639508279497065713301E-01;
weight[ 41-1] = 0.20787750081531811812652137291250E-01;
weight[ 42-1] = 0.21107350591688713643523847921658E-01;
weight[ 43-1] = 0.21414136912893259295449693233545E-01;
weight[ 44-1] = 0.21707922796373466052301324695331E-01;
weight[ 45-1] = 0.21988529885872983756478409758807E-01;
weight[ 46-1] = 0.22255787825930280235631416460158E-01;
weight[ 47-1] = 0.22509534365300608085694429903050E-01;
weight[ 48-1] = 0.22749615455457959852242553240982E-01;
weight[ 49-1] = 0.22975885344117206754377437838947E-01;
weight[ 50-1] = 0.23188206663719640249922582981729E-01;
weight[ 51-1] = 0.23386450514828194170722043496950E-01;
weight[ 52-1] = 0.23570496544381716050033676844306E-01;
weight[ 53-1] = 0.23740233018760777777714726703424E-01;
weight[ 54-1] = 0.23895556891620665983864481754172E-01;
weight[ 55-1] = 0.24036373866450369675132086026456E-01;
weight[ 56-1] = 0.24162598453819584716522917710986E-01;
weight[ 57-1] = 0.24274154023278979833195063936748E-01;
weight[ 58-1] = 0.24370972849882214952813561907241E-01;
weight[ 59-1] = 0.24452996155301467956140198471529E-01;
weight[ 60-1] = 0.24520174143511508275183033290175E-01;
weight[ 61-1] = 0.24572466031020653286354137335186E-01;
weight[ 62-1] = 0.24609840071630254092545634003360E-01;
weight[ 63-1] = 0.24632273575707679066033370218017E-01;
weight[ 64-1] = 0.24639752923961094419579417477503E-01;
weight[ 65-1] = 0.24632273575707679066033370218017E-01;
weight[ 66-1] = 0.24609840071630254092545634003360E-01;
weight[ 67-1] = 0.24572466031020653286354137335186E-01;
weight[ 68-1] = 0.24520174143511508275183033290175E-01;
weight[ 69-1] = 0.24452996155301467956140198471529E-01;
weight[ 70-1] = 0.24370972849882214952813561907241E-01;
weight[ 71-1] = 0.24274154023278979833195063936748E-01;
weight[ 72-1] = 0.24162598453819584716522917710986E-01;
weight[ 73-1] = 0.24036373866450369675132086026456E-01;
weight[ 74-1] = 0.23895556891620665983864481754172E-01;
weight[ 75-1] = 0.23740233018760777777714726703424E-01;
weight[ 76-1] = 0.23570496544381716050033676844306E-01;
weight[ 77-1] = 0.23386450514828194170722043496950E-01;
weight[ 78-1] = 0.23188206663719640249922582981729E-01;
weight[ 79-1] = 0.22975885344117206754377437838947E-01;
weight[ 80-1] = 0.22749615455457959852242553240982E-01;
weight[ 81-1] = 0.22509534365300608085694429903050E-01;
weight[ 82-1] = 0.22255787825930280235631416460158E-01;
weight[ 83-1] = 0.21988529885872983756478409758807E-01;
weight[ 84-1] = 0.21707922796373466052301324695331E-01;
weight[ 85-1] = 0.21414136912893259295449693233545E-01;
weight[ 86-1] = 0.21107350591688713643523847921658E-01;
weight[ 87-1] = 0.20787750081531811812652137291250E-01;
weight[ 88-1] = 0.20455529410639508279497065713301E-01;
weight[ 89-1] = 0.20110890268880247225644623956287E-01;
weight[ 90-1] = 0.19754041885329183081815217323169E-01;
weight[ 91-1] = 0.19385200901246454628112623489471E-01;
weight[ 92-1] = 0.19004591238555646611148901044533E-01;
weight[ 93-1] = 0.18612443963902310429440419898958E-01;
weight[ 94-1] = 0.18208997148375106468721469154479E-01;
weight[ 95-1] = 0.17794495722974774231027912900351E-01;
weight[ 96-1] = 0.17369191329918731922164721250350E-01;
weight[ 97-1] = 0.16933342169871654545878815295200E-01;
weight[ 98-1] = 0.16487212845194879399346060358146E-01;
weight[ 99-1] = 0.16031074199309941802254151842763E-01;
weight[100-1] = 0.15565203152273955098532590262975E-01;
weight[101-1] = 0.15089882532666922992635733981431E-01;
weight[102-1] = 0.14605400905893418351737288078952E-01;
weight[103-1] = 0.14112052399003395774044161633613E-01;
weight[104-1] = 0.13610136522139249906034237533759E-01;
weight[105-1] = 0.13099957986718627426172681912499E-01;
weight[106-1] = 0.12581826520465013101514365424172E-01;
weight[107-1] = 0.12056056679400848183529562144697E-01;
weight[108-1] = 0.11522967656921087154811609734510E-01;
weight[109-1] = 0.10982883090068975788799657376065E-01;
weight[110-1] = 0.10436130863141005225673171997668E-01;
weight[111-1] = 0.98830429087554914716648010899606E-02;
weight[112-1] = 0.93239550065309714787536985834029E-02;
weight[113-1] = 0.87592065795403145773316804234385E-02;
weight[114-1] = 0.81891404887415730817235884718726E-02;
weight[115-1] = 0.76141028256526859356393930849227E-02;
weight[116-1] = 0.70344427036681608755685893032552E-02;
weight[117-1] = 0.64505120486899171845442463868748E-02;
weight[118-1] = 0.58626653903523901033648343751367E-02;
weight[119-1] = 0.52712596565634400891303815906251E-02;
weight[120-1] = 0.46766539777779034772638165662478E-02;
weight[121-1] = 0.40792095178254605327114733456293E-02;
weight[122-1] = 0.34792893810051465908910894094105E-02;
weight[123-1] = 0.28772587656289004082883197417581E-02;
weight[124-1] = 0.22734860707492547802810838362671E-02;
weight[125-1] = 0.16683488125171936761028811985672E-02;
weight[126-1] = 0.10622766869538486959954760554099E-02;
weight[127-1] = 0.45645726109586654495731936146574E-03;
}
else
{
cout << "\n";
cout << "GL_WEIGHTS - Fatal error!\n";
cout << " Illegal value of ORDER = " << order << "\n";
cout << " Legal values are 1, 3, 7, 15, 31, 63 and 127.\n";
exit ( 1 );
}
return;
}
//****************************************************************************80
int i4_log_2 ( int i )
//****************************************************************************80
//