50 double level_weight[] )
124 double *level_weight2;
125 double level_weight_min;
128 if ( alpha_max < 1.0 )
131 std::cerr <<
"SGMGA_ANISO_BALANCE - Fatal error!\n";
132 std::cerr <<
" ALPHA_MAX < 1.0\n";
138 level_weight_min =
webbur->r8_huge ( );
141 for ( dim = 0; dim < dim_num; dim++ )
143 if ( 0.0 < level_weight[dim] )
145 if ( level_weight[dim] < level_weight_min )
147 level_weight_min = level_weight[dim];
148 nonzero_num = nonzero_num + 1;
153 if ( nonzero_num == 0 )
156 std::cerr <<
"SGMGA_ANISO_BALANCE - Fatal error!\n";
157 std::cerr <<
" Could not find a positive entry in LEVEL_WEIGHT.\n";
163 level_weight2 =
new double[dim_num];
164 for ( dim = 0; dim < dim_num; dim++ )
166 level_weight2[dim] = level_weight[dim] / level_weight_min;
171 for ( dim = 0; dim < dim_num; dim++ )
173 level_weight2[dim] =
webbur->r8_min ( alpha_max, level_weight2[dim] );
175 return level_weight2;
182 double level_weight[] )
250 double level_weight_min;
251 double level_weight_sum;
261 else if ( option == 1 )
263 level_weight_min =
webbur->r8_huge ( );
265 for ( dim = 0; dim < dim_num; dim++ )
267 if ( 0.0 < level_weight[dim] )
269 if ( level_weight[dim] < level_weight_min )
271 level_weight_min = level_weight[dim];
280 std::cerr <<
"SGMGA_ANISO_NORMALIZE - Fatal error!\n";
281 std::cerr <<
" Could not find a positive entry in LEVEL_WEIGHT.\n";
285 for ( dim = 0; dim < dim_num; dim++ )
287 level_weight[dim] = level_weight[dim] / level_weight_min;
293 else if ( option == 2 )
295 level_weight_sum =
webbur->r8vec_sum ( dim_num, level_weight );
297 if ( level_weight_sum <= 0.0 )
300 std::cerr <<
"SGMGA_ANISO_NORMALIZE - Fatal error!\n";
301 std::cerr <<
" Sum of level weights is not positive.\n";
304 for ( dim = 0; dim < dim_num; dim++ )
306 level_weight[dim] = ( ( double ) ( dim_num ) * level_weight[dim] )
318 double level_weight[] )
397 for ( dim = 0; dim < dim_num; dim++ )
399 if ( importance[dim] < 0.0 )
402 std::cerr <<
"SGMGA_IMPORTANCE_TO_ANISO - Fatal error!\n";
403 std::cerr <<
" Some IMPORTANCE entries are not positive.\n";
410 for ( dim = 0; dim < dim_num; dim++ )
412 if ( 0.0 < importance[dim] )
414 level_weight[dim] = 1.0 / importance[dim];
419 level_weight[dim] = 0.0;
426 std::cerr <<
"SGMGA_IMPORTANCE_TO_ANISO - Fatal error!\n";
427 std::cerr <<
" No importance entry is positive.\n";
437 double level_weight[],
442 int sparse_unique_index[],
549 double level_weight_min_pos;
567 if ( level_max == 0 )
570 for ( dim = 0; dim < dim_num; dim++ )
572 sparse_order[dim+point*dim_num] = 1;
573 sparse_index[dim+point*dim_num] = 1;
580 for ( point = 0; point < point_num; point++ )
582 for ( dim = 0; dim < dim_num; dim++ )
584 sparse_order[dim+point*dim_num] = -1;
585 sparse_index[dim+point*dim_num] = -1;
591 level_1d =
new int[dim_num];
592 level_1d_max =
new int[dim_num];
593 order_1d =
new int[dim_num];
594 point_index =
new int[dim_num];
598 level_weight_min_pos =
webbur->r8vec_min_pos ( dim_num, level_weight );
599 q_min = ( double ) ( level_max ) * level_weight_min_pos
600 -
webbur->r8vec_sum ( dim_num, level_weight );
601 q_max = ( double ) ( level_max ) * level_weight_min_pos;
602 for ( dim = 0; dim < dim_num; dim++ )
604 if ( 0.0 < level_weight[dim] )
606 level_1d_max[dim] = (int)
webbur->r8_floor(q_max/level_weight[dim])+1;
607 if ( q_max <= ( level_1d_max[dim] - 1 ) * level_weight[dim] )
609 level_1d_max[dim] = level_1d_max[dim] - 1;
614 level_1d_max[dim] = 0;
628 level_1d, q_min, q_max, &more_grids );
646 webbur->level_growth_to_order ( dim_num, level_1d, rule, growth, order_1d );
654 webbur->vec_colex_next3 ( dim_num, order_1d, point_index, &more_points );
660 point_unique = sparse_unique_index[point_count];
661 for ( dim = 0; dim < dim_num; dim++ )
663 sparse_order[dim+point_unique*dim_num] = order_1d[dim];
665 for ( dim = 0; dim < dim_num; dim++ )
667 sparse_index[dim+point_unique*dim_num] = point_index[dim];
669 point_count = point_count + 1;
674 delete [] level_1d_max;
676 delete [] point_index;
684 double level_weight[],
689 void ( *gw_compute_points[] ) (
int order,
697 double sparse_point[] )
800 double level_weight_min_pos;
807 for ( point = 0; point < point_num; point++ )
809 for ( dim = 0; dim < dim_num; dim++ )
811 sparse_point[dim+point*dim_num] = -
webbur->r8_huge ( );
817 level_1d_max =
new int[dim_num];
818 level_weight_min_pos =
webbur->r8vec_min_pos ( dim_num, level_weight );
819 q_max = ( double ) ( level_max ) * level_weight_min_pos;
823 for ( dim = 0; dim < dim_num; dim++ )
825 if ( 0.0 < level_weight[dim] )
827 level_1d_max[dim] = (int)
webbur->r8_floor(q_max/level_weight[dim]) + 1;
828 if ( q_max <= ( level_1d_max[dim] - 1 ) * level_weight[dim] )
830 level_1d_max[dim] = level_1d_max[dim] - 1;
835 level_1d_max[dim] = 0;
838 for ( level = 0; level <= level_1d_max[dim]; level++ )
840 webbur->level_growth_to_order ( 1, &level, rule+dim, growth+dim, &order );
842 points =
new double[order];
844 if ( rule[dim] == 1 )
846 webbur->clenshaw_curtis_compute_points_np (
847 order, np[dim], p+p_index, points );
849 else if ( rule[dim] == 2 )
851 webbur->fejer2_compute_points_np (
852 order, np[dim], p+p_index, points );
854 else if ( rule[dim] == 3 )
856 webbur->patterson_lookup_points_np (
857 order, np[dim], p+p_index, points );
859 else if ( rule[dim] == 4 )
861 webbur->legendre_compute_points_np (
862 order, np[dim], p+p_index, points );
864 else if ( rule[dim] == 5 )
866 webbur->hermite_compute_points_np (
867 order, np[dim], p+p_index, points );
869 else if ( rule[dim] == 6 )
871 webbur->gen_hermite_compute_points_np (
872 order, np[dim], p+p_index, points );
874 else if ( rule[dim] == 7 )
876 webbur->laguerre_compute_points_np (
877 order, np[dim], p+p_index, points );
879 else if ( rule[dim] == 8 )
881 webbur->gen_laguerre_compute_points_np (
882 order, np[dim], p+p_index, points );
884 else if ( rule[dim] == 9 )
886 webbur->jacobi_compute_points_np (
887 order, np[dim], p+p_index, points );
889 else if ( rule[dim] == 10 )
891 webbur->hermite_genz_keister_lookup_points_np (
892 order, np[dim], p+p_index, points );
894 else if ( rule[dim] == 11 )
896 gw_compute_points[dim] (
897 order, np[dim], p+p_index, points );
899 else if ( rule[dim] == 12 )
901 gw_compute_points[dim] (
902 order, np[dim], p+p_index, points );
907 std::cerr <<
"SGMGA_POINT - Fatal error!\n";
908 std::cerr <<
" Unexpected value of RULE[" << dim <<
"] = "
909 << rule[dim] <<
".\n";
913 for ( point = 0; point < point_num; point++ )
915 if ( sparse_order[dim+point*dim_num] == order )
917 sparse_point[dim+point*dim_num] =
918 points[sparse_index[dim+point*dim_num]-1];
923 p_index = p_index + np[dim];
928 for ( point = 0; point < point_num; point++ )
930 for ( dim = 0; dim < dim_num; dim++ )
932 if ( sparse_point[dim+point*dim_num] == -
webbur->r8_huge ( ) )
935 std::cerr <<
"SGMGA_POINT - Fatal error!\n";
936 std::cerr <<
" At least one point component was not assigned.\n";
937 std::cerr <<
" POINT = " << point <<
"\n";
938 std::cerr <<
" DIM = " << dim <<
"\n";
939 std::cerr <<
" SPARSE_ORDER(DIM,POINT) = "
940 << sparse_order[dim+point*dim_num] <<
"\n";
941 std::cerr <<
" LEVEL_WEIGHT(DIM) = " << level_weight[dim] <<
"\n";
947 delete [] level_1d_max;
960 void ( *gw_compute_weights[] )(
int order,
1043 for ( i = 0; i < order_nd; i++ )
1050 for ( dim = 0; dim < dim_num; dim++ )
1052 weight_1d =
new double[order_1d[dim]];
1054 if ( rule[dim] == 1 )
1056 webbur->clenshaw_curtis_compute_weights_np (
1057 order_1d[dim], np[dim], p+p_index, weight_1d );
1059 else if ( rule[dim] == 2 )
1061 webbur->fejer2_compute_weights_np (
1062 order_1d[dim], np[dim], p+p_index, weight_1d );
1064 else if ( rule[dim] == 3 )
1066 webbur->patterson_lookup_weights_np (
1067 order_1d[dim], np[dim], p+p_index, weight_1d );
1069 else if ( rule[dim] == 4 )
1071 webbur->legendre_compute_weights_np (
1072 order_1d[dim], np[dim], p+p_index, weight_1d );
1074 else if ( rule[dim] == 5 )
1076 webbur->hermite_compute_weights_np (
1077 order_1d[dim], np[dim], p+p_index, weight_1d );
1079 else if ( rule[dim] == 6 )
1081 webbur->gen_hermite_compute_weights_np (
1082 order_1d[dim], np[dim], p+p_index, weight_1d );
1084 else if ( rule[dim] == 7 )
1086 webbur->laguerre_compute_weights_np (
1087 order_1d[dim], np[dim], p+p_index, weight_1d );
1089 else if ( rule[dim] == 8 )
1091 webbur->gen_laguerre_compute_weights_np (
1092 order_1d[dim], np[dim], p+p_index, weight_1d );
1094 else if ( rule[dim] == 9 )
1096 webbur->jacobi_compute_weights_np (
1097 order_1d[dim], np[dim], p+p_index, weight_1d );
1099 else if ( rule[dim] == 10 )
1101 webbur->hermite_genz_keister_lookup_weights_np (
1102 order_1d[dim], np[dim], p+p_index, weight_1d );
1104 else if ( rule[dim] == 11 )
1106 gw_compute_weights[dim] (
1107 order_1d[dim], np[dim], p+p_index, weight_1d );
1109 else if ( rule[dim] == 12 )
1111 gw_compute_weights[dim] (
1112 order_1d[dim], np[dim], p+p_index, weight_1d );
1117 std::cerr <<
"SGMGA_PRODUCT_WEIGHT - Fatal error!\n";
1118 std::cerr <<
" Unexpected value of RULE[" << dim <<
"] = "
1119 << rule[dim] <<
".\n";
1123 p_index = p_index + np[dim];
1125 webbur->r8vec_direct_product2 ( dim, order_1d[dim], weight_1d,
1126 dim_num, order_nd, weight_nd );
1128 delete [] weight_1d;
1136 double level_weight[],
1141 void ( *gw_compute_points[] ) (
int order,
1245 double level_weight_min_pos;
1254 int point_total_num;
1255 int point_total_num2;
1260 int *sparse_total_index;
1261 int *sparse_total_order;
1262 double *sparse_total_point;
1266 if ( level_max < 0 )
1272 if ( level_max == 0 )
1286 sparse_total_order =
new int[dim_num*point_total_num];
1287 sparse_total_index =
new int[dim_num*point_total_num];
1289 point_total_num2 = 0;
1291 level_1d =
new int[dim_num];
1292 level_1d_max =
new int[dim_num];
1293 order_1d =
new int[dim_num];
1294 point_index =
new int[dim_num];
1298 level_weight_min_pos =
webbur->r8vec_min_pos ( dim_num, level_weight );
1299 q_min = ( double ) ( level_max ) * level_weight_min_pos
1300 -
webbur->r8vec_sum ( dim_num, level_weight );
1301 q_max = ( double ) ( level_max ) * level_weight_min_pos;
1302 for ( dim = 0; dim < dim_num; dim++ )
1304 if ( 0.0 < level_weight[dim] )
1306 level_1d_max[dim] = (int)
webbur->r8_floor (q_max/level_weight[dim]) + 1;
1307 if ( q_max <= ( level_1d_max[dim] - 1 ) * level_weight[dim] )
1309 level_1d_max[dim] = level_1d_max[dim] - 1;
1314 level_1d_max[dim] = 0;
1328 level_1d, q_min, q_max, &more_grids );
1337 coef =
sgmga_vcn_coef ( dim_num, level_weight, level_1d, q_max );
1346 webbur->level_growth_to_order ( dim_num, level_1d, rule, growth, order_1d );
1350 more_points =
false;
1354 webbur->vec_colex_next3 ( dim_num, order_1d, point_index, &more_points );
1360 for ( dim = 0; dim < dim_num; dim++ )
1362 sparse_total_order[dim+point_total_num2*dim_num] = order_1d[dim];
1364 for ( dim = 0; dim < dim_num; dim++ )
1366 sparse_total_index[dim+point_total_num2*dim_num] = point_index[dim];
1368 point_total_num2 = point_total_num2 + 1;
1373 delete [] point_index;
1377 sparse_total_point =
new double[dim_num*point_total_num];
1379 for ( point = 0; point < point_total_num; point++ )
1381 for ( dim = 0; dim < dim_num; dim++ )
1383 sparse_total_point[dim+point*dim_num] =
webbur->r8_huge ( );
1389 level_1d_max =
new int[dim_num];
1390 level_weight_min_pos =
webbur->r8vec_min_pos ( dim_num, level_weight );
1391 q_max = ( double ) ( level_max ) * level_weight_min_pos;
1394 for ( dim = 0; dim < dim_num; dim++ )
1396 if ( 0.0 < level_weight[dim] )
1398 level_1d_max[dim] = (int)
webbur->r8_floor(q_max/level_weight[dim]) + 1;
1399 if ( q_max <= ( level_1d_max[dim] - 1 ) * level_weight[dim] )
1401 level_1d_max[dim] = level_1d_max[dim] - 1;
1406 level_1d_max[dim] = 0;
1409 for ( level = 0; level <= level_1d_max[dim]; level++ )
1411 webbur->level_growth_to_order ( 1, &level, rule+dim, growth+dim, &order );
1413 points =
new double[order];
1415 if ( rule[dim] == 1 )
1417 webbur->clenshaw_curtis_compute_points_np (
1418 order, np[dim], p+p_index, points );
1420 else if ( rule[dim] == 2 )
1422 webbur->fejer2_compute_points_np (
1423 order, np[dim], p+p_index, points );
1425 else if ( rule[dim] == 3 )
1427 webbur->patterson_lookup_points_np (
1428 order, np[dim], p+p_index, points );
1430 else if ( rule[dim] == 4 )
1432 webbur->legendre_compute_points_np (
1433 order, np[dim], p+p_index, points );
1435 else if ( rule[dim] == 5 )
1437 webbur->hermite_compute_points_np (
1438 order, np[dim], p+p_index, points );
1440 else if ( rule[dim] == 6 )
1442 webbur->gen_hermite_compute_points_np (
1443 order, np[dim], p+p_index, points );
1445 else if ( rule[dim] == 7 )
1447 webbur->laguerre_compute_points_np (
1448 order, np[dim], p+p_index, points );
1450 else if ( rule[dim] == 8 )
1452 webbur->gen_laguerre_compute_points_np (
1453 order, np[dim], p+p_index, points );
1455 else if ( rule[dim] == 9 )
1457 webbur->jacobi_compute_points_np (
1458 order, np[dim], p+p_index, points );
1460 else if ( rule[dim] == 10 )
1462 webbur->hermite_genz_keister_lookup_points_np (
1463 order, np[dim], p+p_index, points );
1465 else if ( rule[dim] == 11 )
1467 gw_compute_points[dim] (
1468 order, np[dim], p+p_index, points );
1470 else if ( rule[dim] == 12 )
1472 gw_compute_points[dim] (
1473 order, np[dim], p+p_index, points );
1478 std::cerr <<
"SGMGA_SIZE - Fatal error!\n";
1479 std::cerr <<
" Unexpected value of RULE[" << dim <<
"] = "
1480 << rule[dim] <<
".\n";
1484 for ( point = 0; point < point_total_num; point++ )
1486 if ( sparse_total_order[dim+point*dim_num] == order )
1488 sparse_total_point[dim+point*dim_num] =
1489 points[sparse_total_index[dim+point*dim_num]-1];
1494 p_index = p_index + np[dim];
1501 point_num =
webbur->point_radial_tol_unique_count ( dim_num, point_total_num,
1502 sparse_total_point, tol, &seed );
1504 delete [] level_1d_max;
1505 delete [] sparse_total_index;
1506 delete [] sparse_total_order;
1507 delete [] sparse_total_point;
1515 double level_weight[],
1624 double level_weight_min_pos;
1627 int point_total_num;
1633 if ( level_max == 0 )
1635 point_total_num = 1;
1636 return point_total_num;
1639 point_total_num = 0;
1641 level_1d =
new int[dim_num];
1642 level_1d_max =
new int[dim_num];
1643 order_1d =
new int[dim_num];
1647 level_weight_min_pos =
webbur->r8vec_min_pos ( dim_num, level_weight );
1648 q_min = ( double ) ( level_max ) * level_weight_min_pos
1649 -
webbur->r8vec_sum ( dim_num, level_weight );
1650 q_max = ( double ) ( level_max ) * level_weight_min_pos;
1651 for ( dim = 0; dim < dim_num; dim++ )
1653 if ( 0.0 < level_weight[dim] )
1655 level_1d_max[dim] = (int)
webbur->r8_floor(q_max/level_weight[dim]) + 1;
1656 if ( q_max <= ( level_1d_max[dim] - 1 ) * level_weight[dim] )
1658 level_1d_max[dim] = level_1d_max[dim] - 1;
1663 level_1d_max[dim] = 0;
1677 level_1d, q_min, q_max, &more_grids );
1686 coef =
sgmga_vcn_coef ( dim_num, level_weight, level_1d, q_max );
1695 webbur->level_growth_to_order ( dim_num, level_1d, rule, growth, order_1d );
1697 point_total_num = point_total_num +
webbur->i4vec_product ( dim_num,
1701 delete [] level_1d_max;
1704 return point_total_num;
1710 double level_weight[],
1715 void ( *gw_compute_points[] ) (
int order,
1721 int point_total_num,
1723 int sparse_unique_index[] )
1833 double level_weight_min_pos;
1842 int point_total_num2;
1848 int *sparse_total_index;
1849 int *sparse_total_order;
1850 double *sparse_total_point;
1855 if ( level_max < 0 )
1860 if ( level_max == 0 )
1862 sparse_unique_index[0] = 0;
1869 sparse_total_order =
new int[dim_num*point_total_num];
1870 sparse_total_index =
new int[dim_num*point_total_num];
1872 level_1d =
new int[dim_num];
1873 level_1d_max =
new int[dim_num];
1874 order_1d =
new int[dim_num];
1875 point_index =
new int[dim_num];
1877 point_total_num2 = 0;
1881 level_weight_min_pos =
webbur->r8vec_min_pos ( dim_num, level_weight );
1882 q_min = ( double ) ( level_max ) * level_weight_min_pos
1883 -
webbur->r8vec_sum ( dim_num, level_weight );
1884 q_max = ( double ) ( level_max ) * level_weight_min_pos;
1885 for ( dim = 0; dim < dim_num; dim++ )
1887 if ( 0.0 < level_weight[dim] )
1889 level_1d_max[dim] = (int)
webbur->r8_floor(q_max/level_weight[dim]) + 1;
1890 if ( q_max <= ( level_1d_max[dim] - 1 ) * level_weight[dim] )
1892 level_1d_max[dim] = level_1d_max[dim] - 1;
1897 level_1d_max[dim] = 0;
1911 level_1d, q_min, q_max, &more_grids );
1920 coef =
sgmga_vcn_coef ( dim_num, level_weight, level_1d, q_max );
1929 webbur->level_growth_to_order ( dim_num, level_1d, rule, growth, order_1d );
1933 more_points =
false;
1937 webbur->vec_colex_next3 ( dim_num, order_1d, point_index, &more_points );
1943 for ( dim = 0; dim < dim_num; dim++ )
1945 sparse_total_order[dim+point_total_num2*dim_num] = order_1d[dim];
1947 for ( dim = 0; dim < dim_num; dim++ )
1949 sparse_total_index[dim+point_total_num2*dim_num] = point_index[dim];
1951 point_total_num2 = point_total_num2 + 1;
1955 delete [] level_1d_max;
1957 delete [] point_index;
1961 sparse_total_point =
new double[dim_num*point_total_num];
1963 for ( point = 0; point < point_total_num; point++ )
1965 for ( dim = 0; dim < dim_num; dim++ )
1967 sparse_total_point[dim+point*dim_num] =
webbur->r8_huge ( );
1973 level_1d_max =
new int[dim_num];
1974 level_weight_min_pos =
webbur->r8vec_min_pos ( dim_num, level_weight );
1975 q_max = ( double ) ( level_max ) * level_weight_min_pos;
1978 for ( dim = 0; dim < dim_num; dim++ )
1980 if ( 0.0 < level_weight[dim] )
1982 level_1d_max[dim] = (int)
webbur->r8_floor(q_max/level_weight[dim]) + 1;
1983 if ( q_max <= ( level_1d_max[dim] - 1 ) * level_weight[dim] )
1985 level_1d_max[dim] = level_1d_max[dim] - 1;
1990 level_1d_max[dim] = 0;
1993 for ( level = 0; level <= level_1d_max[dim]; level++ )
1995 webbur->level_growth_to_order ( 1, &level, rule+dim, growth+dim, &order );
1997 points =
new double[order];
1999 if ( rule[dim] == 1 )
2001 webbur->clenshaw_curtis_compute_points_np (
2002 order, np[dim], p+p_index, points );
2004 else if ( rule[dim] == 2 )
2006 webbur->fejer2_compute_points_np (
2007 order, np[dim], p+p_index, points );
2009 else if ( rule[dim] == 3 )
2011 webbur->patterson_lookup_points_np (
2012 order, np[dim], p+p_index, points );
2014 else if ( rule[dim] == 4 )
2016 webbur->legendre_compute_points_np (
2017 order, np[dim], p+p_index, points );
2019 else if ( rule[dim] == 5 )
2021 webbur->hermite_compute_points_np (
2022 order, np[dim], p+p_index, points );
2024 else if ( rule[dim] == 6 )
2026 webbur->gen_hermite_compute_points_np (
2027 order, np[dim], p+p_index, points );
2029 else if ( rule[dim] == 7 )
2031 webbur->laguerre_compute_points_np (
2032 order, np[dim], p+p_index, points );
2034 else if ( rule[dim] == 8 )
2036 webbur->gen_laguerre_compute_points_np (
2037 order, np[dim], p+p_index, points );
2039 else if ( rule[dim] == 9 )
2041 webbur->jacobi_compute_points_np (
2042 order, np[dim], p+p_index, points );
2044 else if ( rule[dim] == 10 )
2046 webbur->hermite_genz_keister_lookup_points_np (
2047 order, np[dim], p+p_index, points );
2049 else if ( rule[dim] == 11 )
2051 gw_compute_points[dim] (
2052 order, np[dim], p+p_index, points );
2054 else if ( rule[dim] == 12 )
2056 gw_compute_points[dim] (
2057 order, np[dim], p+p_index, points );
2062 std::cerr <<
"SGMGA_UNIQUE_INDEX - Fatal error!\n";
2063 std::cerr <<
" Unexpected value of RULE[" << dim <<
"] = "
2064 << rule[dim] <<
".\n";
2068 for ( point = 0; point < point_total_num; point++ )
2070 if ( sparse_total_order[dim+point*dim_num] == order )
2072 sparse_total_point[dim+point*dim_num] =
2073 points[sparse_total_index[dim+point*dim_num]-1];
2078 p_index = p_index + np[dim];
2085 undx =
new int[point_num];
2089 webbur->point_radial_tol_unique_index ( dim_num, point_total_num,
2090 sparse_total_point, tol, &seed, undx, sparse_unique_index );
2092 for ( point = 0; point < point_total_num; point++ )
2094 rep = undx[sparse_unique_index[point]];
2097 for ( dim = 0; dim < dim_num; dim++ )
2099 sparse_total_point[dim+point*dim_num] = sparse_total_point[dim+rep*dim_num];
2106 webbur->point_unique_index ( dim_num, point_total_num, sparse_total_point,
2107 point_num, undx, sparse_unique_index );
2111 delete [] sparse_total_index;
2112 delete [] sparse_total_order;
2113 delete [] sparse_total_point;
2252 for ( i = 0; i < n; i++ )
2265 if ( nstart == - 1 )
2268 std::cerr <<
" SGMGA_VCN - Fatal error!\n";
2269 std::cerr <<
" No weight is positive.\n";
2275 for ( i = 0; i < n; i++ )
2303 else if ( dir == - 1 || dir == 0 )
2315 else if ( nstart < n2 )
2319 for ( i = n2 + 1; i < n; i++ )
2321 t = t - w[i] * ( double ) x[i];
2323 xmax[n2] = (int)
webbur->r8_floor ( t / w[n2] );
2325 else if ( n2 == nstart && dir == - 1 )
2328 for ( i = n2 + 1; i < n; i++ )
2330 t = t - w[i] * ( double ) x[i];
2332 xmin = (int)
webbur->r8_ceiling ( t / w[n2] );
2338 for ( i = 0; i < n2; i++ )
2340 t = t + w[i] * ( double ) x[i];
2342 t = t + w[n2] * xmin;
2343 for ( i = n2 + 1; i < n; i++ )
2345 t = t + w[i] * ( double ) x[i];
2353 for ( i = n2 + 1; i < n; i++ )
2355 t = t - w[i] * ( double ) x[i];
2357 xmax[n2] = (int)
webbur->r8_floor ( t / w[n2] );
2360 if ( xmax[n2] < xmin )
2373 else if ( dir == 0 )
2376 if ( x[n2] <= xmax[n2] )
2396 else if ( dir == + 1 )
2412 if ( x[n2] < xmax[n2] )
2428 double level_weight[],
2567 b =
new int[dim_num];
2569 for ( i = 0; i < dim_num; i++ )
2581 while ( i < dim_num - 1 )
2587 if ( level_weight[i] == 0.0 )
2593 else if ( b[i] == 1 )
2610 for ( j = 0; j < dim_num; j++ )
2612 q = q + level_weight[j] * ( double ) ( x[j] + b[j] );
2627 while ( i < dim_num - 1 )
2631 if ( level_weight[i] == 0.0 )
2634 else if ( b[i] == 1 )
2649 for ( j = 0; j < dim_num; j++ )
2651 b_sum = b_sum + b[j];
2656 c = c + 1 - 2 * ( b_sum % 2 );
2665 coef = ( double ) ( c );
2674 double level_weight[],
2761 b =
new int[dim_num];
2763 for ( i = 0; i < dim_num; i++ )
2774 webbur->binary_vector_next ( dim_num, b );
2775 b_sum =
webbur->i4vec_sum ( dim_num, b );
2788 for ( i = 0; i < dim_num; i++ )
2790 if ( x_max[i] < x[i] + b[i] )
2804 for ( i = 0; i < dim_num; i++ )
2806 q = q + level_weight[i] * ( double ) ( x[i] + b[i] );
2811 coef = coef +
webbur->r8_mop ( b_sum );
2823 double level_weight[],
2954 for ( i = 0; i < dim_num; i++ )
2960 for ( i = 0; i < dim_num; i++ )
2962 q = q + level_weight[i] * ( double ) ( x[i] );
2965 if ( q_min < q && q <= q_max )
2977 if ( x[j] < x_max[j] )
2982 if ( dim_num - 1 <= j )
2991 for ( i = 0; i < j; i++ )
2997 for ( i = 0; i < dim_num; i++ )
2999 q = q + level_weight[i] * ( double ) ( x[i] );
3002 if ( q_min < q && q <= q_max )
3014 double level_weight[],
3146 static double q_max2;
3147 static double q_min2;
3154 q_max2 =
webbur->r8_min ( q_min + 1.0, q_max );
3161 sgmga_vcn ( dim_num, level_weight, x, q_min2, q_max2, more );
3172 if ( q_max2 < q_max )
3175 q_max2 =
webbur->r8_min ( q_max2 + 1.0, q_max );
3192 double level_weight[],
3324 static double q_max2;
3325 static double q_min2;
3332 q_max2 =
webbur->r8_min ( q_min + 1.0, q_max );
3351 if ( q_max2 < q_max )
3354 q_max2 =
webbur->r8_min ( q_max2 + 1.0, q_max );
3371 double level_weight[],
3376 void ( *gw_compute_weights[] ) (
int order,
3381 int point_total_num,
3382 int sparse_unique_index[],
3384 double sparse_weight[] )
3477 double *grid_weight;
3481 double level_weight_min_pos;
3492 for ( point = 0; point < point_num; point++ )
3494 sparse_weight[point] = 0.0;
3499 level_1d =
new int[dim_num];
3500 order_1d =
new int[dim_num];
3501 level_1d_max =
new int[dim_num];
3505 level_weight_min_pos =
webbur->r8vec_min_pos ( dim_num, level_weight );
3506 q_min = ( double ) ( level_max ) * level_weight_min_pos
3507 -
webbur->r8vec_sum ( dim_num, level_weight );
3508 q_max = ( double ) ( level_max ) * level_weight_min_pos;
3509 for ( dim = 0; dim < dim_num; dim++ )
3511 if ( 0.0 < level_weight[dim] )
3513 level_1d_max[dim] = (int)
webbur->r8_floor(q_max/level_weight[dim]) + 1;
3514 if ( q_max <= ( level_1d_max[dim] - 1 ) * level_weight[dim] )
3516 level_1d_max[dim] = level_1d_max[dim] - 1;
3521 level_1d_max[dim] = 0;
3535 level_1d, q_min, q_max, &more_grids );
3544 coef =
sgmga_vcn_coef ( dim_num, level_weight, level_1d, q_max );
3553 webbur->level_growth_to_order ( dim_num, level_1d, rule, growth, order_1d );
3557 order_nd =
webbur->i4vec_product ( dim_num, order_1d );
3566 grid_weight =
new double[order_nd];
3569 np, p, gw_compute_weights, grid_weight );
3573 for ( order = 0; order < order_nd; order++ )
3575 point_unique = sparse_unique_index[point_total];
3577 point_total = point_total + 1;
3579 sparse_weight[point_unique] = sparse_weight[point_unique]
3580 + coef * grid_weight[order];
3583 delete [] grid_weight;
3587 delete [] level_1d_max;
3596 double level_weight[],
3601 double sparse_weight[],
3602 double sparse_point[],
3603 std::string file_name )
3685 std::string file_name_a;
3686 std::string file_name_n;
3687 std::string file_name_p;
3688 std::string file_name_r;
3689 std::string file_name_w;
3690 std::string file_name_x;
3693 double *sparse_region;
3697 sparse_region =
new double[dim_num*2];
3699 for ( dim = 0; dim < dim_num; dim++ )
3701 if ( rule[dim] == 1 )
3703 sparse_region[dim+0*dim_num] = -1.0;
3704 sparse_region[dim+1*dim_num] = +1.0;
3706 else if ( rule[dim] == 2 )
3708 sparse_region[dim+0*dim_num] = -1.0;
3709 sparse_region[dim+1*dim_num] = +1.0;
3711 else if ( rule[dim] == 3 )
3713 sparse_region[dim+0*dim_num] = -1.0;
3714 sparse_region[dim+1*dim_num] = +1.0;
3716 else if ( rule[dim] == 4 )
3718 sparse_region[dim+0*dim_num] = -1.0;
3719 sparse_region[dim+1*dim_num] = +1.0;
3721 else if ( rule[dim] == 5 )
3723 sparse_region[dim+0*dim_num] = -
webbur->r8_huge ( );
3724 sparse_region[dim+1*dim_num] = +
webbur->r8_huge ( );
3726 else if ( rule[dim] == 6 )
3728 sparse_region[dim+0*dim_num] = -
webbur->r8_huge ( );
3729 sparse_region[dim+1*dim_num] = +
webbur->r8_huge ( );
3731 else if ( rule[dim] == 7 )
3733 sparse_region[dim+0*dim_num] = 0.0;
3734 sparse_region[dim+1*dim_num] =
webbur->r8_huge ( );
3736 else if ( rule[dim] == 8 )
3738 sparse_region[dim+0*dim_num] = 0.0;
3739 sparse_region[dim+1*dim_num] =
webbur->r8_huge ( );
3741 else if ( rule[dim] == 9 )
3743 sparse_region[dim+0*dim_num] = -1.0;
3744 sparse_region[dim+1*dim_num] = +1.0;
3746 else if ( rule[dim] == 10 )
3748 sparse_region[dim+0*dim_num] = -
webbur->r8_huge ( );
3749 sparse_region[dim+1*dim_num] = +
webbur->r8_huge ( );
3754 else if ( rule[dim] == 11 )
3756 t1 =
webbur->r8_huge ( );
3757 t2 = -
webbur->r8_huge ( );
3758 for ( point = 0; point < point_num; point++ )
3760 t1 =
webbur->r8_min ( t1, sparse_point[dim+point*dim_num] );
3761 t2 =
webbur->r8_max ( t2, sparse_point[dim+point*dim_num] );
3763 sparse_region[dim+0*dim_num] = t1;
3764 sparse_region[dim+1*dim_num] = t2;
3766 else if ( rule[dim] == 12 )
3768 t1 =
webbur->r8_huge ( );
3769 t2 = -
webbur->r8_huge ( );
3770 for ( point = 0; point < point_num; point++ )
3772 t1 =
webbur->r8_min ( t1, sparse_point[dim+point*dim_num] );
3773 t2 =
webbur->r8_max ( t2, sparse_point[dim+point*dim_num] );
3775 sparse_region[dim+0*dim_num] = t1;
3776 sparse_region[dim+1*dim_num] = t2;
3781 std::cerr <<
"SGMGA_WRITE - Fatal error!\n";
3782 std::cerr <<
" Unexpected value of RULE[" << dim <<
"] = "
3783 << rule[dim] <<
".\n";
3788 std::cout <<
"SGMGA_WRITE:\n";
3790 file_name_a = file_name +
"_a.txt";
3791 webbur->r8mat_write ( file_name_a, 1, dim_num, level_weight );
3792 std::cout <<
" Wrote the A file = \"" << file_name_a <<
"\".\n";
3794 file_name_n = file_name +
"_n.txt";
3795 webbur->i4mat_write ( file_name_n, 1, dim_num, np );
3796 std::cout <<
" Wrote the N file = \"" << file_name_n <<
"\".\n";
3798 np_sum =
webbur->i4vec_sum ( dim_num, np );
3799 file_name_p = file_name +
"_p.txt";
3800 webbur->r8mat_write ( file_name_p, 1, np_sum, p );
3801 std::cout <<
" Wrote the P file = \"" << file_name_p <<
"\".\n";
3803 file_name_r = file_name +
"_r.txt";
3804 webbur->r8mat_write ( file_name_r, dim_num, 2, sparse_region );
3805 std::cout <<
" Wrote the R file = \"" << file_name_r <<
"\".\n";
3807 file_name_w = file_name +
"_w.txt";
3808 webbur->r8mat_write ( file_name_w, 1, point_num, sparse_weight );
3809 std::cout <<
" Wrote the W file = \"" << file_name_w <<
"\".\n";
3811 file_name_x = file_name +
"_x.txt";
3812 webbur->r8mat_write ( file_name_x, dim_num, point_num, sparse_point );
3813 std::cout <<
" Wrote the X file = \"" << file_name_x <<
"\".\n";
3815 delete [] sparse_region;
void sgmga_index(int dim_num, double level_weight[], int level_max, int rule[], int point_num, int point_total_num, int sparse_unique_index[], int growth[], int sparse_order[], int sparse_index[])
void sgmga_point(int dim_num, double level_weight[], int level_max, int rule[], int np[], double p[], void(*gw_compute_points[])(int order, int np, double p[], double x[]), int point_num, int sparse_order[], int sparse_index[], int growth[], double sparse_point[])
double sgmga_vcn_coef_naive(int n, double level_weight[], int x_max[], int x[], double q_min, double q_max)
void sgmga_importance_to_aniso(int dim_num, double importance[], double level_weight[])
int sgmga_size_total(int dim_num, double level_weight[], int level_max, int rule[], int growth[])
void sgmga_product_weight(int dim_num, int order_1d[], int order_nd, int rule[], int np[], double p[], void(*gw_compute_weights[])(int order, int np, double p[], double w[]), double weight_nd[])
double sgmga_vcn_coef(int n, double level_weight[], int x[], double q_max)
double * sgmga_aniso_balance(double alpha_max, int dim_num, double level_weight[])
void sgmga_weight(int dim_num, double level_weight[], int level_max, int rule[], int np[], double p[], void(*gw_compute_weights[])(int order, int np, double p[], double w[]), int point_num, int point_total_num, int sparse_unique_index[], int growth[], double sparse_weight[])
void sgmga_unique_index(int dim_num, double level_weight[], int level_max, int rule[], int np[], double p[], void(*gw_compute_points[])(int order, int np, double p[], double x[]), double tol, int point_num, int point_total_num, int growth[], int sparse_unique_index[])
Teuchos::RCP< SandiaRules > webbur
void sgmga_aniso_normalize(int option, int dim_num, double level_weight[])
void sgmga_write(int dim_num, double level_weight[], int rule[], int np[], double p[], int point_num, double sparse_weight[], double sparse_point[], std::string file_name)
int sgmga_size(int dim_num, double level_weight[], int level_max, int rule[], int np[], double p[], void(*gw_compute_points[])(int order, int np, double p[], double x[]), double tol, int growth[])
void sgmga_vcn_ordered(int dim_num, double level_weight[], int x_max[], int x[], double q_min, double q_max, bool *more)
void sgmga_vcn_ordered_naive(int dim_num, double level_weight[], int x_max[], int x[], double q_min, double q_max, bool *more)
void sgmga_vcn_naive(int n, double level_weight[], int x_max[], int x[], double q_min, double q_max, bool *more)
void sgmga_vcn(int n, double level_weight[], int x[], double q_min, double q_max, bool *more)