# include # include # include # include # include # include using namespace std; # include "sparse_grid_open.hpp" //****************************************************************************80 int *abscissa_level_open_nd ( int level_max, int dim_num, int test_num, int test_val[] ) //****************************************************************************80 // // Purpose: // // ABSCISSA_LEVEL_OPEN_ND: first level at which given abscissa is generated. // // Discussion: // // We assume an underlying product grid. In each dimension, this product // grid has order 2**(LEVEL_MAX+1) - 1. // // We will say a sparse grid has total level LEVEL if each point in the // grid has a total level of LEVEL or less. // // The "level" of a point is determined as the sum of the levels of the // point in each spatial dimension. // // The level of a point in a single spatial dimension I is determined as // the level, between 0 and LEVEL_MAX, at which the point's I'th index // would have been generated. // // // This description is terse and perhaps unenlightening. Keep in mind // that the product grid is the product of 1D grids, // that the 1D grids are built up by levels, having // orders (total number of points ) 1, 3, 7, 15, 31 and so on, // and that these 1D grids are nested, so that each point in a 1D grid // has a first level at which it appears. // // Our procedure for generating the points of a sparse grid, then, is // to choose a value LEVEL_MAX, to generate the full product grid, // but then only to keep those points on the full product grid whose // LEVEL is less than or equal to LEVEL_MAX. // // // Note that this routine is really just testing out the idea of // determining the level. Our true desire is to be able to start // with a value LEVEL, and determine, in a straightforward manner, // all the points that are generated exactly at that level, or // all the points that are generated up to and including that level. // // This allows us to generate the new points to be added to one sparse // grid to get the next, or to generate a particular sparse grid at once. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 19 April 2007 // // 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. // // Parameters: // // Input, int LEVEL_MAX, controls the size of the final sparse grid. // // Input, int DIM_NUM, the spatial dimension. // // Input, int TEST_NUM, the number of points to be tested. // // Input, int TEST_VAL[DIM_NUM*TEST_NUM], the indices of the points // to be tested. Normally, each index would be between 0 and 2**LEVEL_MAX. // // Output, int ABSCISSA_OPEN_LEVEL_ND[TEST_NUM], the value of LEVEL at which the // point would first be generated, assuming that a standard sequence of // nested grids is used. // { int dim; int j; int level; int order; int t; int *test_level; test_level = new int[test_num]; if ( level_max == 0 ) { for ( j = 0; j < test_num; j++ ) { test_level[j] = 0; } return test_level; } order = i4_power ( 2, level_max ) + 1; for ( j = 0; j < test_num; j++ ) { test_level[j] = index_to_level_open ( dim_num, test_val+j*dim_num, order, level_max ); } return test_level; } //****************************************************************************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 double f2_abscissa ( int order, int i ) //****************************************************************************80 // // Purpose: // // F2_ABSCISSA returns the I-th abscissa for the Fejer type 2 rule. // // Discussion: // // Our convention is that the abscissas are numbered from left to // right. // // This rule is defined on [-1,1]. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 11 April 2007 // // Author: // // John Burkardt // // Parameters: // // Input, int ORDER, the order of the Fejer type 2 rule. // 1 <= ORDER. // // Input, int I, the index of the desired abscissa. 1 <= I <= ORDER. // // Output, double F2_ABSCISSA, the value of the I-th // abscissa in the Fejer type 2 rule of order ORDER. // { double pi = 3.141592653589793; double value; if ( order < 1 ) { value = - r8_huge ( ); return value; } if ( i < 1 || order < i ) { cout << "\n"; cout << "F2_ABSCISSA - Fatal error!\n"; cout << " 1 <= I <= ORDER is required.\n"; exit ( 1 ); } if ( order == 1 ) { value = 0.0; return value; } value = cos ( ( double ) ( order + 1 - i ) * pi / ( double ) ( order + 1 ) ); return value; } //****************************************************************************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, -0.99982213041530614629963254927125E+00, -0.99906293435531189513828920479421E+00, -0.99769756618980462107441703193392E+00, -0.99572655135202722663543337085008E+00, -0.99315104925451714736113079489080E+00, -0.98997261459148415760778669967548E+00, -0.98619317401693166671043833175407E+00, -0.98181502080381411003346312451200E+00, -0.97684081234307032681744391886221E+00, -0.97127356816152919228894689830512E+00, -0.96511666794529212109082507703391E+00, -0.95837384942523877114910286998060E+00, -0.95104920607788031054790764659636E+00, -0.94314718462481482734544963026201E+00, -0.93467258232473796857363487794906E+00, -0.92563054405623384912746466814259E+00, -0.91602655919146580931308861741716E+00, -0.90586645826182138280246131760282E+00, -0.89515640941708370896904382642451E+00, -0.88390291468002656994525794802849E+00, -0.87211280599856071141963753428864E+00, -0.85979324109774080981203134414483E+00, -0.84695169913409759845333931085437E+00, -0.83359597615489951437955716480123E+00, -0.81973418036507867415511910167470E+00, -0.80537472720468021466656079404644E+00, -0.79052633423981379994544995252740E+00, -0.77519801587020238244496276354566E+00, -0.75939907785653667155666366659810E+00, -0.74313911167095451292056688997595E+00, -0.72642798867407268553569290153270E+00, -0.70927585412210456099944463906757E+00, -0.69169312100770067015644143286666E+00, -0.67369046373825048534668253831602E+00, -0.65527881165548263027676505156852E+00, -0.63646934240029724134760815684175E+00, -0.61727347512685828385763916340822E+00, -0.59770286357006522938441201887478E+00, -0.57776938897061258000325165713764E+00, -0.55748515286193223292186190687872E+00, -0.53686246972339756745816636353452E+00, -0.51591385950424935727727729906662E+00, -0.49465204002278211739494017368636E+00, -0.47308991924540524164509989939699E+00, -0.45124058745026622733189858020729E+00, -0.42911730928019337626254405355418E+00, -0.40673351568978256340867288124339E+00, -0.38410279579151693577907781452239E+00, -0.36123888860586970607092484346723E+00, -0.33815567472039850137600027657095E+00, -0.31486716786289498148601475374890E+00, -0.29138750639370562079451875284568E+00, -0.26773094472238862088834352027938E+00, -0.24391184465391785797071324453138E+00, -0.21994466666968754245452337866940E+00, -0.19584396114861085150428162519610E+00, -0.17162435953364216500834492248954E+00, -0.14730056544908566938932929319807E+00, -0.12288734577408297172603365288567E+00, -0.98399521677698970751091751509101E-01, -0.73851959621048545273440409360569E-01, -0.49259562331926630315379321821927E-01, -0.24637259757420944614897071846088E-01, 0.00000000000000000000000000000000E+00, 0.24637259757420944614897071846088E-01, 0.49259562331926630315379321821927E-01, 0.73851959621048545273440409360569E-01, 0.98399521677698970751091751509101E-01, 0.12288734577408297172603365288567E+00, 0.14730056544908566938932929319807E+00, 0.17162435953364216500834492248954E+00, 0.19584396114861085150428162519610E+00, 0.21994466666968754245452337866940E+00, 0.24391184465391785797071324453138E+00, 0.26773094472238862088834352027938E+00, 0.29138750639370562079451875284568E+00, 0.31486716786289498148601475374890E+00, 0.33815567472039850137600027657095E+00, 0.36123888860586970607092484346723E+00, 0.38410279579151693577907781452239E+00, 0.40673351568978256340867288124339E+00, 0.42911730928019337626254405355418E+00, 0.45124058745026622733189858020729E+00, 0.47308991924540524164509989939699E+00, 0.49465204002278211739494017368636E+00, 0.51591385950424935727727729906662E+00, 0.53686246972339756745816636353452E+00, 0.55748515286193223292186190687872E+00, 0.57776938897061258000325165713764E+00, 0.59770286357006522938441201887478E+00, 0.61727347512685828385763916340822E+00, 0.63646934240029724134760815684175E+00, 0.65527881165548263027676505156852E+00, 0.67369046373825048534668253831602E+00, 0.69169312100770067015644143286666E+00, 0.70927585412210456099944463906757E+00, 0.72642798867407268553569290153270E+00, 0.74313911167095451292056688997595E+00, 0.75939907785653667155666366659810E+00, 0.77519801587020238244496276354566E+00, 0.79052633423981379994544995252740E+00, 0.80537472720468021466656079404644E+00, 0.81973418036507867415511910167470E+00, 0.83359597615489951437955716480123E+00, 0.84695169913409759845333931085437E+00, 0.85979324109774080981203134414483E+00, 0.87211280599856071141963753428864E+00, 0.88390291468002656994525794802849E+00, 0.89515640941708370896904382642451E+00, 0.90586645826182138280246131760282E+00, 0.91602655919146580931308861741716E+00, 0.92563054405623384912746466814259E+00, 0.93467258232473796857363487794906E+00, 0.94314718462481482734544963026201E+00, 0.95104920607788031054790764659636E+00, 0.95837384942523877114910286998060E+00, 0.96511666794529212109082507703391E+00, 0.97127356816152919228894689830512E+00, 0.97684081234307032681744391886221E+00, 0.98181502080381411003346312451200E+00, 0.98619317401693166671043833175407E+00, 0.98997261459148415760778669967548E+00, 0.99315104925451714736113079489080E+00, 0.99572655135202722663543337085008E+00, 0.99769756618980462107441703193392E+00, 0.99906293435531189513828920479421E+00, 0.99982213041530614629963254927125E+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 double gp_abscissa ( int order, int index ) //****************************************************************************80 // // Purpose: // // GP_ABSCISSA returns the I-th abscissa for a Gauss-Patterson rule. // // Discussion: // // The rule is specified by its level. // // The number of points in the rule, known as the order, is // related to the level by the formula: // // ORDER = 2^(LEVEL+1)-1. // // Only rules of order 1, 3, 7, 15, 31, 63, 127 and 255 are allowed. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 26 December 2009 // // Author: // // John Burkardt // // Reference: // // Prem Kythe, Michael Schaeferkotter, // Handbook of Computational Methods for Integration, // Chapman and Hall, 2004, // ISBN: 1-58488-428-2, // LC: QA299.3.K98. // // Thomas Patterson, // The Optimal Addition of Points to Quadrature Formulae, // Mathematics of Computation, // Volume 22, Number 104, October 1968, pages 847-856. // // Parameters: // // Input, int ORDER, the order of the rule. // // Input, int INDEX, the index of the point in the rule. // // Output, double GP_ABSCISSA, the value of the INDEX-th // abscissa in the rule of order ORDER. // { static double x_001[1] = { 0.0 }; static double x_003[3] = { -0.77459666924148337704, 0.0, 0.77459666924148337704 }; static double x_007[7] = { -0.96049126870802028342, -0.77459666924148337704, -0.43424374934680255800, 0.0, 0.43424374934680255800, 0.77459666924148337704, 0.96049126870802028342 }; static double x_015[15] = { -0.99383196321275502221, -0.96049126870802028342, -0.88845923287225699889, -0.77459666924148337704, -0.62110294673722640294, -0.43424374934680255800, -0.22338668642896688163, 0.0, 0.22338668642896688163, 0.43424374934680255800, 0.62110294673722640294, 0.77459666924148337704, 0.88845923287225699889, 0.96049126870802028342, 0.99383196321275502221 }; static double x_031[31] = { -0.99909812496766759766, -0.99383196321275502221, -0.98153114955374010687, -0.96049126870802028342, -0.92965485742974005667, -0.88845923287225699889, -0.83672593816886873550, -0.77459666924148337704, -0.70249620649152707861, -0.62110294673722640294, -0.53131974364437562397, -0.43424374934680255800, -0.33113539325797683309, -0.22338668642896688163, -0.11248894313318662575, 0.0, 0.11248894313318662575, 0.22338668642896688163, 0.33113539325797683309, 0.43424374934680255800, 0.53131974364437562397, 0.62110294673722640294, 0.70249620649152707861, 0.77459666924148337704, 0.83672593816886873550, 0.88845923287225699889, 0.92965485742974005667, 0.96049126870802028342, 0.98153114955374010687, 0.99383196321275502221, 0.99909812496766759766 }; static double x_063[63] = { -0.99987288812035761194, -0.99909812496766759766, -0.99720625937222195908, -0.99383196321275502221, -0.98868475754742947994, -0.98153114955374010687, -0.97218287474858179658, -0.96049126870802028342, -0.94634285837340290515, -0.92965485742974005667, -0.91037115695700429250, -0.88845923287225699889, -0.86390793819369047715, -0.83672593816886873550, -0.80694053195021761186, -0.77459666924148337704, -0.73975604435269475868, -0.70249620649152707861, -0.66290966002478059546, -0.62110294673722640294, -0.57719571005204581484, -0.53131974364437562397, -0.48361802694584102756, -0.43424374934680255800, -0.38335932419873034692, -0.33113539325797683309, -0.27774982202182431507, -0.22338668642896688163, -0.16823525155220746498, -0.11248894313318662575, -0.056344313046592789972, 0.0, 0.056344313046592789972, 0.11248894313318662575, 0.16823525155220746498, 0.22338668642896688163, 0.27774982202182431507, 0.33113539325797683309, 0.38335932419873034692, 0.43424374934680255800, 0.48361802694584102756, 0.53131974364437562397, 0.57719571005204581484, 0.62110294673722640294, 0.66290966002478059546, 0.70249620649152707861, 0.73975604435269475868, 0.77459666924148337704, 0.80694053195021761186, 0.83672593816886873550, 0.86390793819369047715, 0.88845923287225699889, 0.91037115695700429250, 0.92965485742974005667, 0.94634285837340290515, 0.96049126870802028342, 0.97218287474858179658, 0.98153114955374010687, 0.98868475754742947994, 0.99383196321275502221, 0.99720625937222195908, 0.99909812496766759766, 0.99987288812035761194 }; static double x_127[127] = { -0.99998243035489159858, -0.99987288812035761194, -0.99959879967191068325, -0.99909812496766759766, -0.99831663531840739253, -0.99720625937222195908, -0.99572410469840718851, -0.99383196321275502221, -0.99149572117810613240, -0.98868475754742947994, -0.98537149959852037111, -0.98153114955374010687, -0.97714151463970571416, -0.97218287474858179658, -0.96663785155841656709, -0.96049126870802028342, -0.95373000642576113641, -0.94634285837340290515, -0.93832039777959288365, -0.92965485742974005667, -0.92034002547001242073, -0.91037115695700429250, -0.89974489977694003664, -0.88845923287225699889, -0.87651341448470526974, -0.86390793819369047715, -0.85064449476835027976, -0.83672593816886873550, -0.82215625436498040737, -0.80694053195021761186, -0.79108493379984836143, -0.77459666924148337704, -0.75748396638051363793, -0.73975604435269475868, -0.72142308537009891548, -0.70249620649152707861, -0.68298743109107922809, -0.66290966002478059546, -0.64227664250975951377, -0.62110294673722640294, -0.59940393024224289297, -0.57719571005204581484, -0.55449513263193254887, -0.53131974364437562397, -0.50768775753371660215, -0.48361802694584102756, -0.45913001198983233287, -0.43424374934680255800, -0.40897982122988867241, -0.38335932419873034692, -0.35740383783153215238, -0.33113539325797683309, -0.30457644155671404334, -0.27774982202182431507, -0.25067873030348317661, -0.22338668642896688163, -0.19589750271110015392, -0.16823525155220746498, -0.14042423315256017459, -0.11248894313318662575, -0.084454040083710883710, -0.056344313046592789972, -0.028184648949745694339, 0.0, 0.028184648949745694339, 0.056344313046592789972, 0.084454040083710883710, 0.11248894313318662575, 0.14042423315256017459, 0.16823525155220746498, 0.19589750271110015392, 0.22338668642896688163, 0.25067873030348317661, 0.27774982202182431507, 0.30457644155671404334, 0.33113539325797683309, 0.35740383783153215238, 0.38335932419873034692, 0.40897982122988867241, 0.43424374934680255800, 0.45913001198983233287, 0.48361802694584102756, 0.50768775753371660215, 0.53131974364437562397, 0.55449513263193254887, 0.57719571005204581484, 0.59940393024224289297, 0.62110294673722640294, 0.64227664250975951377, 0.66290966002478059546, 0.68298743109107922809, 0.70249620649152707861, 0.72142308537009891548, 0.73975604435269475868, 0.75748396638051363793, 0.77459666924148337704, 0.79108493379984836143, 0.80694053195021761186, 0.82215625436498040737, 0.83672593816886873550, 0.85064449476835027976, 0.86390793819369047715, 0.87651341448470526974, 0.88845923287225699889, 0.89974489977694003664, 0.91037115695700429250, 0.92034002547001242073, 0.92965485742974005667, 0.93832039777959288365, 0.94634285837340290515, 0.95373000642576113641, 0.96049126870802028342, 0.96663785155841656709, 0.97218287474858179658, 0.97714151463970571416, 0.98153114955374010687, 0.98537149959852037111, 0.98868475754742947994, 0.99149572117810613240, 0.99383196321275502221, 0.99572410469840718851, 0.99720625937222195908, 0.99831663531840739253, 0.99909812496766759766, 0.99959879967191068325, 0.99987288812035761194, 0.99998243035489159858 }; static double x_255[255] = { -0.99999759637974846462, -0.99998243035489159858, -0.99994399620705437576, -0.99987288812035761194, -0.99976049092443204733, -0.99959879967191068325, -0.99938033802502358193, -0.99909812496766759766, -0.99874561446809511470, -0.99831663531840739253, -0.99780535449595727456, -0.99720625937222195908, -0.99651414591489027385, -0.99572410469840718851, -0.99483150280062100052, -0.99383196321275502221, -0.99272134428278861533, -0.99149572117810613240, -0.99015137040077015918, -0.98868475754742947994, -0.98709252795403406719, -0.98537149959852037111, -0.98351865757863272876, -0.98153114955374010687, -0.97940628167086268381, -0.97714151463970571416, -0.97473445975240266776, -0.97218287474858179658, -0.96948465950245923177, -0.96663785155841656709, -0.96364062156981213252, -0.96049126870802028342, -0.95718821610986096274, -0.95373000642576113641, -0.95011529752129487656, -0.94634285837340290515, -0.94241156519108305981, -0.93832039777959288365, -0.93406843615772578800, -0.92965485742974005667, -0.92507893290707565236, -0.92034002547001242073, -0.91543758715576504064, -0.91037115695700429250, -0.90514035881326159519, -0.89974489977694003664, -0.89418456833555902286, -0.88845923287225699889, -0.88256884024734190684, -0.87651341448470526974, -0.87029305554811390585, -0.86390793819369047715, -0.85735831088623215653, -0.85064449476835027976, -0.84376688267270860104, -0.83672593816886873550, -0.82952219463740140018, -0.82215625436498040737, -0.81462878765513741344, -0.80694053195021761186, -0.79909229096084140180, -0.79108493379984836143, -0.78291939411828301639, -0.77459666924148337704, -0.76611781930376009072, -0.75748396638051363793, -0.74869629361693660282, -0.73975604435269475868, -0.73066452124218126133, -0.72142308537009891548, -0.71203315536225203459, -0.70249620649152707861, -0.69281376977911470289, -0.68298743109107922809, -0.67301883023041847920, -0.66290966002478059546, -0.65266166541001749610, -0.64227664250975951377, -0.63175643771119423041, -0.62110294673722640294, -0.61031811371518640016, -0.59940393024224289297, -0.58836243444766254143, -0.57719571005204581484, -0.56590588542365442262, -0.55449513263193254887, -0.54296566649831149049, -0.53131974364437562397, -0.51955966153745702199, -0.50768775753371660215, -0.49570640791876146017, -0.48361802694584102756, -0.47142506587165887693, -0.45913001198983233287, -0.44673538766202847374, -0.43424374934680255800, -0.42165768662616330006, -0.40897982122988867241, -0.39621280605761593918, -0.38335932419873034692, -0.37042208795007823014, -0.35740383783153215238, -0.34430734159943802278, -0.33113539325797683309, -0.31789081206847668318, -0.30457644155671404334, -0.29119514851824668196, -0.27774982202182431507, -0.26424337241092676194, -0.25067873030348317661, -0.23705884558982972721, -0.22338668642896688163, -0.20966523824318119477, -0.19589750271110015392, -0.18208649675925219825, -0.16823525155220746498, -0.15434681148137810869, -0.14042423315256017459, -0.12647058437230196685, -0.11248894313318662575, -0.098482396598119202090, -0.084454040083710883710, -0.070406976042855179063, -0.056344313046592789972, -0.042269164765363603212, -0.028184648949745694339, -0.014093886410782462614, 0.0, 0.014093886410782462614, 0.028184648949745694339, 0.042269164765363603212, 0.056344313046592789972, 0.070406976042855179063, 0.084454040083710883710, 0.098482396598119202090, 0.11248894313318662575, 0.12647058437230196685, 0.14042423315256017459, 0.15434681148137810869, 0.16823525155220746498, 0.18208649675925219825, 0.19589750271110015392, 0.20966523824318119477, 0.22338668642896688163, 0.23705884558982972721, 0.25067873030348317661, 0.26424337241092676194, 0.27774982202182431507, 0.29119514851824668196, 0.30457644155671404334, 0.31789081206847668318, 0.33113539325797683309, 0.34430734159943802278, 0.35740383783153215238, 0.37042208795007823014, 0.38335932419873034692, 0.39621280605761593918, 0.40897982122988867241, 0.42165768662616330006, 0.43424374934680255800, 0.44673538766202847374, 0.45913001198983233287, 0.47142506587165887693, 0.48361802694584102756, 0.49570640791876146017, 0.50768775753371660215, 0.51955966153745702199, 0.53131974364437562397, 0.54296566649831149049, 0.55449513263193254887, 0.56590588542365442262, 0.57719571005204581484, 0.58836243444766254143, 0.59940393024224289297, 0.61031811371518640016, 0.62110294673722640294, 0.63175643771119423041, 0.64227664250975951377, 0.65266166541001749610, 0.66290966002478059546, 0.67301883023041847920, 0.68298743109107922809, 0.69281376977911470289, 0.70249620649152707861, 0.71203315536225203459, 0.72142308537009891548, 0.73066452124218126133, 0.73975604435269475868, 0.74869629361693660282, 0.75748396638051363793, 0.76611781930376009072, 0.77459666924148337704, 0.78291939411828301639, 0.79108493379984836143, 0.79909229096084140180, 0.80694053195021761186, 0.81462878765513741344, 0.82215625436498040737, 0.82952219463740140018, 0.83672593816886873550, 0.84376688267270860104, 0.85064449476835027976, 0.85735831088623215653, 0.86390793819369047715, 0.87029305554811390585, 0.87651341448470526974, 0.88256884024734190684, 0.88845923287225699889, 0.89418456833555902286, 0.89974489977694003664, 0.90514035881326159519, 0.91037115695700429250, 0.91543758715576504064, 0.92034002547001242073, 0.92507893290707565236, 0.92965485742974005667, 0.93406843615772578800, 0.93832039777959288365, 0.94241156519108305981, 0.94634285837340290515, 0.95011529752129487656, 0.95373000642576113641, 0.95718821610986096274, 0.96049126870802028342, 0.96364062156981213252, 0.96663785155841656709, 0.96948465950245923177, 0.97218287474858179658, 0.97473445975240266776, 0.97714151463970571416, 0.97940628167086268381, 0.98153114955374010687, 0.98351865757863272876, 0.98537149959852037111, 0.98709252795403406719, 0.98868475754742947994, 0.99015137040077015918, 0.99149572117810613240, 0.99272134428278861533, 0.99383196321275502221, 0.99483150280062100052, 0.99572410469840718851, 0.99651414591489027385, 0.99720625937222195908, 0.99780535449595727456, 0.99831663531840739253, 0.99874561446809511470, 0.99909812496766759766, 0.99938033802502358193, 0.99959879967191068325, 0.99976049092443204733, 0.99987288812035761194, 0.99994399620705437576, 0.99998243035489159858, 0.99999759637974846462 }; double value; if ( order < 1 ) { value = - r8_huge ( ); } else if ( index < 1 || order < index ) { value = - r8_huge ( ); } else if ( order == 1 ) { value = x_001[index-1]; } else if ( order == 3 ) { value = x_003[index-1]; } else if ( order == 7 ) { value = x_007[index-1]; } else if ( order == 15 ) { value = x_015[index-1]; } else if ( order == 31 ) { value = x_031[index-1]; } else if ( order == 63 ) { value = x_063[index-1]; } else if ( order == 127 ) { value = x_127[index-1]; } else if ( order == 255 ) { value = x_255[index-1]; } else { value = - r8_huge ( ); } return value; } //****************************************************************************80 int i4_choose ( int n, int k ) //****************************************************************************80 // // Purpose: // // I4_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: // // 02 June 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 I4_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 int i4_log_2 ( int i ) //****************************************************************************80 // // Purpose: // // I4_LOG_2 returns the integer part of the logarithm base 2 of an I4. // // Example: // // I I4_LOG_10 // ----- -------- // 0 0 // 1 0 // 2 1 // 3 1 // 4 2 // 5 2 // 7 2 // 8 3 // 9 3 // 1000 9 // 1024 10 // // Discussion: // // I4_LOG_2 ( I ) + 1 is the number of binary digits in I. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 04 January 2004 // // Author: // // John Burkardt // // Parameters: // // Input, int I, the number whose logarithm base 2 is desired. // // Output, int I4_LOG_2, the integer part of the logarithm base 2 of // the absolute value of X. // { int i_abs; int two_pow; int value; if ( i == 0 ) { value = 0; } else { value = 0; two_pow = 2; i_abs = abs ( i ); while ( two_pow <= i_abs ) { value = value + 1; two_pow = two_pow * 2; } } return value; } //****************************************************************************80 int i4_max ( int i1, int i2 ) //****************************************************************************80 // // Purpose: // // I4_MAX returns the maximum of two I4's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 13 October 1998 // // Author: // // John Burkardt // // Parameters: // // Input, int I1, I2, are two integers to be compared. // // Output, int I4_MAX, the larger of I1 and I2. // { int value; if ( i2 < i1 ) { value = i1; } else { value = i2; } return value; } //****************************************************************************80 int i4_min ( int i1, int i2 ) //****************************************************************************80 // // Purpose: // // I4_MIN returns the minimum of two I4's. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 13 October 1998 // // Author: // // John Burkardt // // Parameters: // // Input, int I1, I2, two integers to be compared. // // Output, int I4_MIN, the smaller of I1 and I2. // { int value; if ( i1 < i2 ) { value = i1; } else { value = i2; } return value; } //****************************************************************************80 int i4_modp ( int i, int j ) //****************************************************************************80 // // Purpose: // // I4_MODP returns the nonnegative remainder of I4 division. // // Formula: // // If // NREM = I4_MODP ( I, J ) // NMULT = ( I - NREM ) / J // then // I = J * NMULT + NREM // where NREM is always nonnegative. // // The MOD function computes a result with the same sign as the // quantity being divided. Thus, suppose you had an angle A, // and you wanted to ensure that it was between 0 and 360. // Then mod(A,360) would do, if A was positive, but if A // was negative, your result would be between -360 and 0. // // On the other hand, I4_MODP(A,360) is between 0 and 360, always. // // Example: // // I J MOD I4_MODP I4_MODP Factorization // // 107 50 7 7 107 = 2 * 50 + 7 // 107 -50 7 7 107 = -2 * -50 + 7 // -107 50 -7 43 -107 = -3 * 50 + 43 // -107 -50 -7 43 -107 = 3 * -50 + 43 // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 26 May 1999 // // Author: // // John Burkardt // // Parameters: // // Input, int I, the number to be divided. // // Input, int J, the number that divides I. // // Output, int I4_MODP, the nonnegative remainder when I is // divided by J. // { int value; if ( j == 0 ) { cout << "\n"; cout << "I4_MODP - Fatal error!\n"; cout << " I4_MODP ( I, J ) called with J = " << j << "\n"; exit ( 1 ); } value = i % j; if ( value < 0 ) { value = value + abs ( j ); } return value; } //****************************************************************************80 int i4_power ( int i, int j ) //****************************************************************************80 // // Purpose: // // I4_POWER returns the value of I^J. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 01 April 2004 // // Author: // // John Burkardt // // Parameters: // // Input, int I, J, the base and the power. J should be nonnegative. // // Output, int I4_POWER, the value of I^J. // { int k; int value; if ( j < 0 ) { if ( i == 1 ) { value = 1; } else if ( i == 0 ) { cout << "\n"; cout << "I4_POWER - Fatal error!\n"; cout << " I^J requested, with I = 0 and J negative.\n"; exit ( 1 ); } else { value = 0; } } else if ( j == 0 ) { if ( i == 0 ) { cout << "\n"; cout << "I4_POWER - Fatal error!\n"; cout << " I^J requested, with I = 0 and J = 0.\n"; exit ( 1 ); } else { value = 1; } } else if ( j == 1 ) { value = i; } else { value = 1; for ( k = 1; k <= j; k++ ) { value = value * i; } } return value; } //****************************************************************************80 string i4_to_string ( int i4, string format ) //****************************************************************************80 // // Purpose: // // I4_TO_STRING converts an I4 to a C++ string. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 09 July 2009 // // Author: // // John Burkardt // // Parameters: // // Input, int I4, an integer. // // Input, string FORMAT, the format string. // // Output, string I4_TO_STRING, the string. // { char i4_char[80]; string i4_string; sprintf ( i4_char, format.c_str ( ), i4 ); i4_string = string ( i4_char ); return i4_string; } //****************************************************************************80 int i4vec_product ( int n, int a[] ) //****************************************************************************80 // // Purpose: // // I4VEC_PRODUCT multiplies the entries of an I4VEC. // // Example: // // Input: // // A = ( 1, 2, 3, 4 ) // // Output: // // I4VEC_PRODUCT = 24 // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 17 May 2003 // // Author: // // John Burkardt // // Parameters: // // Input, int N, the number of entries in the vector. // // Input, int A[N], the vector // // Output, int I4VEC_PRODUCT, the product of the entries of A. // { int i; int product; product = 1; for ( i = 0; i < n; i++ ) { product = product * a[i]; } return product; } //****************************************************************************80 int index_to_level_open ( int dim_num, int t[], int order, int level_max ) //****************************************************************************80 // // Purpose: // // INDEX_TO_LEVEL_OPEN determines the level of a point given its index. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 19 April 2007 // // 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. // // Parameters: // // Input, int DIM_NUM, the spatial dimension. // // Input, int T[DIM_NUM], the grid index of a point. // // Input, int ORDER, the order of the rule. // // Input, int LEVEL_MAX, the level with respect to which the // index applies. // // Output, int INDEX_TO_LEVEL_OPEN, the first level on which // the point associated with the given index will appear. // { int dim; int level; int s; int value; value = 0; for ( dim = 0; dim < dim_num; dim++ ) { s = t[dim]; s = i4_modp ( s, order ); if ( s == 0 ) { level = 0; } else { level = level_max; while ( ( s % 2 ) == 0 ) { s = s / 2; level = level - 1; } } if ( level == 0 ) { level = 1; } else if ( level == 1 ) { level = 0; } value = value + level; } return value; } //****************************************************************************80 void level_to_order_open ( int dim_num, int level[], int order[] ) //****************************************************************************80 // // Purpose: // // LEVEL_TO_ORDER_OPEN converts a level to an order for open rules. // // Discussion: // // Sparse grids can naturally be nested. A natural scheme is to use // a series of one-dimensional rules arranged in a series of "levels" // whose order roughly doubles with each step. // // The arrangement described here works naturally for the // Fejer Type 2, Newton Cotes Open, // and Gauss-Patterson rules. It also can be used, partially, to describe // the growth of Gauss-Legendre rules. // // The idea is that we start with LEVEL = 0, ORDER = 1 indicating the single // point at the center, and for all values afterwards, we use the relationship // // ORDER = 2**(LEVEL+1) - 1. // // The following table shows how the growth will occur: // // Level Order // // 0 1 // 1 3 = 4 - 1 // 2 7 = 8 - 1 // 3 15 = 16 - 1 // 4 31 = 32 - 1 // 5 63 = 64 - 1 // // For the Fejer Type 2, Newton Cotes Open, // and Gauss-Patterson rules, the point growth is // nested. If we have ORDER points on a particular LEVEL, the next level // includes all these old points, plus ORDER+1 new points, formed in the // gaps between successive pairs of old points plus an extra point at each // end. // // Level Order = New + Old // // 0 1 = 1 + 0 // 1 3 = 2 + 1 // 2 7 = 4 + 3 // 3 15 = 8 + 7 // 4 31 = 16 + 15 // 5 63 = 32 + 31 // // If we use a series of Gauss-Legendre rules, then there is almost no // nesting, except that the central point is shared. If we insist on // producing a comparable series of such points, then the "nesting" behavior // is as follows: // // Level Order = New + Old // // 0 1 = 1 + 0 // 1 3 = 2 + 1 // 2 7 = 6 + 1 // 3 15 = 14 + 1 // 4 31 = 30 + 1 // 5 63 = 62 + 1 // // Moreover, if we consider ALL the points used in such a set of "nested" // Gauss-Legendre rules, then we must sum the "NEW" column, and we see that // we get roughly twice as many points as for the truly nested rules. // // In this routine, we assume that a vector of levels is given, // and the corresponding orders are desired. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 19 April 2007 // // 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. // // Parameters: // // Input, int DIM_NUM, the spatial dimension. // // Input, int LEVEL[DIM_NUM], the nesting level. // // Output, int ORDER[DIM_NUM], the order (number of points) // of the rule. // { int dim; for ( dim = 0; dim < dim_num; dim++ ) { if ( level[dim] < 0 ) { order[dim] = -1; } else if ( level[dim] == 0 ) { order[dim] = 1; } else { order[dim] = i4_power ( 2, level[dim] + 1 ) - 1 ; } } return; } //****************************************************************************80 int *levels_open_index ( int dim_num, int level_max, int point_num ) //****************************************************************************80 // // Purpose: // // LEVELS_OPEN_INDEX computes grids with 0 <= LEVEL <= LEVEL_MAX. // // Discussion: // // The necessary dimensions of GRID_INDEX can be // determined by calling LEVELS_OPEN_INDEX_SIZE first. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 05 July 2008 // // 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. // // Parameters: // // Input, int DIM_NUM, the spatial dimension. // // Input, int LEVEL_MAX, the maximum value of LEVEL. // // Input, int POINT_NUM, the total number of points in the grids. // // Output, int LEVELS_MAX_INDEX[DIM_NUM*POINT_NUM], a list of point indices, // representing a subset of the product grid of level LEVEL_MAX, // representing (exactly once) each point that will show up in a // sparse grid of level LEVEL_MAX. // { int dim; int factor; int *grid_index; int *grid_index2; int h; int level; int *level_1d; bool more; int *order_1d; int order_nd; int point; int point_num2; int t; bool test; // // The outer loop generates LEVELs from 0 to LEVEL_MAX. // grid_index = new int[dim_num*point_num]; level_1d = new int[dim_num]; order_1d = new int[dim_num]; point_num2 = 0; for ( level = 0; level <= level_max; level++ ) { // // The middle loop generates the next partition that adds up to LEVEL. // more = false; h = 0; t = 0; for ( ; ; ) { comp_next ( level, dim_num, level_1d, &more, &h, &t ); // // Transform each 1D level to a corresponding 1D order. // level_to_order_open ( dim_num, level_1d, order_1d ); // // The product of the 1D orders gives us the number of points in this grid. // order_nd = i4vec_product ( dim_num, order_1d ); // // The inner (hidden) loop generates all points corresponding to given grid. // grid_index2 = multigrid_index1 ( dim_num, order_1d, order_nd ); // // Only keep those points which first appear on this level. // If you keep a point, it is necessary to rescale each of its components // so that we save the coordinates as they apply on the final grid. // for ( point = 0; point < order_nd; point++ ) { test = true; for ( dim = 0; dim < dim_num; dim++ ) { if ( grid_index2[dim+point*dim_num] % 2 == 0 ) { test = false; } } if ( test ) { for ( dim = 0; dim < dim_num; dim++ ) { grid_index[dim+point_num2*dim_num] = i4_power ( 2, level_max - level_1d[dim] ) * grid_index2[dim+point*dim_num]; } point_num2 = point_num2 + 1; } } delete [] grid_index2; if ( !more ) { break; } } } delete [] level_1d; delete [] order_1d; return grid_index; } //****************************************************************************80 int *multigrid_index1 ( int dim_num, int order_1d[], int order_nd ) //****************************************************************************80 // // Purpose: // // MULTIGRID_INDEX1 returns an indexed multidimensional grid. // // Discussion: // // For dimension DIM, the second index of INDX may vary from // 1 to ORDER_1D[DIM]. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 23 May 2007 // // 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. // // Parameters: // // Input, int DIM_NUM, the spatial dimension of the points. // // Input, int ORDER_1D[DIM_NUM], the order of the // rule in each dimension. // // Input, int ORDER_ND, the product of the entries of ORDER_1D. // // Output, int INDX[DIM_NUM*ORDER_ND], the indices of the points in // the grid. The second dimension of this array is equal to the // product of the entries of ORDER_1D. // { int *a; int dim; bool more; int p; int *indx; indx = new int[dim_num*order_nd]; a = new int[dim_num]; more = false; p = 0; for ( ; ; ) { vec_colex_next2 ( dim_num, order_1d, a, &more ); if ( !more ) { break; } for ( dim = 0; dim < dim_num; dim++ ) { indx[dim+p*dim_num] = a[dim] + 1; } p = p + 1; } delete [] a; return indx; } //****************************************************************************80 void multigrid_scale_open ( int dim_num, int order_nd, int level_max, int level_1d[], int grid_index[] ) //****************************************************************************80 // // Purpose: // // MULTIGRID_SCALE_OPEN renumbers a grid as a subgrid on a higher level. // // Discussion: // // This routine takes a grid associated with a given value of // LEVEL, and multiplies all the indices by a power of 2, so that // the indices reflect the position of the same points, but in // a grid of level LEVEL_MAX. // // For an open grid, going from one level to the next, a set of indices // will be rescaled by 2*INDEX-1. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 08 June 2007 // // 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. // // Parameters: // // Input, int DIM_NUM, the spatial dimension. // // Input, int ORDER_ND, the number of points in the grid. // // Input, int LEVEL_MAX, the maximum value of LEVEL. // // Input, int LEVEL_1D[DIM_NUM], the level in each dimension. // // Input/output, int GRID_INDEX[DIM_NUM*POINT_NUM], the index // values for each grid point. On input, these indices are based in // the level for which the grid was generated; on output, the // indices are appropriate for the grid as a subgrid of a grid // of level LEVEL_MAX. // { int dim; int factor; int order; for ( dim = 0; dim < dim_num; dim++ ) { factor = i4_power ( 2, level_max - level_1d[dim] ); for ( order = 0; order < order_nd; order++ ) { grid_index[dim+order*dim_num] = grid_index[dim+order*dim_num] * factor; } } return; } //****************************************************************************80 double nco_abscissa ( int order, int i ) //****************************************************************************80 // // Purpose: // // NCO_ABSCISSA returns the I-th abscissa for the Newton Cotes open rule. // // Discussion: // // Our convention is that the abscissas are numbered from left to // right. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 25 May 2007 // // Author: // // John Burkardt // // Parameters: // // Input, int ORDER, the order of the rule. // 1 <= ORDER. // // Input, int I, the index of the desired abscissa. // 1 <= I <= ORDER. // // Output, double NCO_ABSCISSA, the value of the I-th // abscissa in the Newton Cotes open rule of order ORDER. // { double value; double x_max = +1.0; double x_min = -1.0; if ( order < 1 ) { value = - r8_huge ( ); return value; } if ( i < 1 || order < i ) { cout << "\n"; cout << "NCO_ABSCISSA - Fatal error!\n"; cout << " 1 <= I <= ORDER is required.\n"; exit ( 1 ); } if ( order == 1 ) { value = ( x_min + x_max ) / 2.0; return value; } value = ( ( double ) ( order - i + 1 ) * x_min + ( double ) ( i ) * x_max ) / ( double ) ( order + 1 ); return value; } //****************************************************************************80 double r8_epsilon ( ) //****************************************************************************80 // // Purpose: // // R8_EPSILON returns the R8 roundoff unit. // // Discussion: // // The roundoff unit is a number R which is a power of 2 with the // property that, to the precision of the computer's arithmetic, // 1 < 1 + R // but // 1 = ( 1 + R / 2 ) // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 01 September 2012 // // Author: // // John Burkardt // // Parameters: // // Output, double R8_EPSILON, the R8 round-off unit. // { const double value = 2.220446049250313E-016; return value; } //****************************************************************************80 double r8_huge ( ) //****************************************************************************80 // // Purpose: // // R8_HUGE returns a "huge" R8. // // Discussion: // // The value returned by this function is NOT required to be the // maximum representable R8. This value varies from machine to machine, // from compiler to compiler, and may cause problems when being printed. // We simply want a "very large" but non-infinite number. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 06 October 2007 // // Author: // // John Burkardt // // Parameters: // // Output, double R8_HUGE, a "huge" R8 value. // { double value; value = 1.0E+30; return value; } //****************************************************************************80 void r8mat_write ( string output_filename, int m, int n, double table[] ) //****************************************************************************80 // // Purpose: // // R8MAT_WRITE writes an R8MAT file. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 29 June 2009 // // Author: // // John Burkardt // // Parameters: // // Input, string OUTPUT_FILENAME, the output filename. // // Input, int M, the spatial dimension. // // Input, int N, the number of points. // // Input, double TABLE[M*N], the table data. // { int i; int j; ofstream output; // // Open the file. // output.open ( output_filename.c_str ( ) ); if ( !output ) { cerr << "\n"; cerr << "R8MAT_WRITE - Fatal error!\n"; cerr << " Could not open the output file.\n"; return; } // // Write the data. // for ( j = 0; j < n; j++ ) { for ( i = 0; i < m; i++ ) { output << " " << setw(24) << setprecision(16) << table[i+j*m]; } output << "\n"; } // // Close the file. // output.close ( ); return; } //****************************************************************************80 int sparse_grid_f2s_size ( int dim_num, int level_max ) //****************************************************************************80 // // Purpose: // // SPARSE_GRID_F2S_SIZE sizes a sparse grid using Fejer Type 2 Slow rules. // // Discussion: // // The grid is defined as the sum of the product rules whose LEVEL // satisfies: // // 0 <= LEVEL <= LEVEL_MAX. // // This calculation is much faster than a previous method. It simply // computes the number of new points that are added at each level in the // 1D rule, and then counts the new points at a given DIM_NUM dimensional // level vector as the product of the new points added in each dimension. // // This approach will work for nested families, and may be extensible // to other families, and to mixed rules. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 26 December 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. // // Parameters: // // Input, int DIM_NUM, the spatial dimension. // // Input, int LEVEL_MAX, the maximum value of LEVEL. // // Output, int SPARSE_GRID_F2S_SIZE, the number of points in the grid. // { int dim; int h; int l; int level; int *level_1d; bool more; int *new_1d; int o; int p; int point_num; int t; int v; // // Special case. // if ( level_max < 0 ) { point_num = 0; return point_num; } if ( level_max == 0 ) { point_num = 1; return point_num; } // // Construct the vector that counts the new points in the 1D rule. // new_1d = new int[level_max+1]; new_1d[0] = 1; p = 1; o = 1; for ( l = 1; l <= level_max; l++ ) { p = 2 * l + 1; if ( o < p ) { new_1d[l] = o + 1; o = 2 * o + 1; } else { new_1d[l] = 0; } } // // Count the number of points by counting the number of new points // associated with each level vector. // level_1d = new int[dim_num]; point_num = 0; for ( level = 0; level <= level_max; level++ ) { more = false; h = 0; t = 0; for ( ; ;) { comp_next ( level, dim_num, level_1d, &more, &h, &t ); v = 1; for ( dim = 0; dim < dim_num; dim++ ) { v = v * new_1d[level_1d[dim]]; } point_num = point_num + v; if ( !more ) { break; } } } delete [] level_1d; delete [] new_1d; return point_num; } //****************************************************************************80 int sparse_grid_gps_size ( int dim_num, int level_max ) //****************************************************************************80 // // Purpose: // // SPARSE_GRID_GPS_SIZE sizes a sparse grid using Gauss-Patterson-Slow rules. // // Discussion: // // The Gauss-Patterson-Slow family assumes that, for the underlying 1D // rules, a precision of 2*L+1 is needed at level L. Therefore, the // lowest possible order Gauss-Patterson rule is chosen that will achieve // that precision. This retains a combination of the advantages of // nestedness and high accuracy. // // The grid is defined as the sum of the product rules whose LEVEL // satisfies: // // 0 <= LEVEL <= LEVEL_MAX. // // This calculation is much faster than a previous method. It simply // computes the number of new points that are added at each level in the // 1D rule, and then counts the new points at a given DIM_NUM dimensional // level vector as the product of the new points added in each dimension. // // This approach will work for nested families, and may be extensible // to other families, and to mixed rules. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 25 December 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. // // Parameters: // // Input, int DIM_NUM, the spatial dimension. // // Input, int LEVEL_MAX, the maximum value of LEVEL. // // Output, int SPARSE_GRID_CC_SIZE, the number of points in the grid. // { int dim; int h; int j; int level; int *level_1d; bool more; int *new_1d; int o; int *order_1d; int p; int point_num; int t; int v; // // Special case. // if ( level_max < 0 ) { point_num = 0; return point_num; } if ( level_max == 0 ) { point_num = 1; return point_num; } // // Count the points in the 1D rule. // order_1d = new int[level_max+1]; order_1d[0] = 1; for ( level = 1; level <= level_max; level++ ) { p = 5; o = 3; while ( p < 2 * level + 1 ) { p = 2 * p + 1; o = 2 * o + 1; } order_1d[level] = o; } // // Count the new points in the 1D rule. // new_1d = new int[level_max+1]; new_1d[0] = 1; for ( level = 1; level <= level_max; level++ ) { new_1d[level] = order_1d[level] - order_1d[level-1]; } // // Count the number of points by counting the number of new points // associated with each level vector. // level_1d = new int[dim_num]; point_num = 0; for ( level = 0; level <= level_max; level++ ) { more = false; h = 0; t = 0; for ( ; ;) { comp_next ( level, dim_num, level_1d, &more, &h, &t ); v = 1; for ( dim = 0; dim < dim_num; dim++ ) { v = v * new_1d[level_1d[dim]]; } point_num = point_num + v; if ( !more ) { break; } } } delete [] level_1d; delete [] new_1d; delete [] order_1d; return point_num; } //****************************************************************************80 int sparse_grid_ofn_size ( int dim_num, int level_max ) //****************************************************************************80 // // Purpose: // // SPARSE_GRID_OFN_SIZE sizes a sparse grid using Open Fully Nested rules. // // Discussion: // // The grid is defined as the sum of the product rules whose LEVEL // satisfies: // // 0 <= LEVEL <= LEVEL_MAX. // // This calculation is much faster than a previous method. It simply // computes the number of new points that are added at each level in the // 1D rule, and then counts the new points at a given DIM_NUM dimensional // level vector as the product of the new points added in each dimension. // // This approach will work for nested families, and may be extensible // to other families, and to mixed rules. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 23 December 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. // // Parameters: // // Input, int DIM_NUM, the spatial dimension. // // Input, int LEVEL_MAX, the maximum value of LEVEL. // // Output, int SPARSE_GRID_OFN_SIZE, the number of points in the grid. // { int dim; int h; int j; int l; int level; int *level_1d; bool more; int *new_1d; int point_num; int t; int v; // // Special case. // if ( level_max < 0 ) { point_num = 0; return point_num; } if ( level_max == 0 ) { point_num = 1; return point_num; } // // Construct the vector that counts the new points in the 1D rule. // new_1d = new int[level_max+1]; new_1d[0] = 1; for ( l = 1; l <= level_max; l++ ) { new_1d[l] = 2 * new_1d[l-1]; } // // Count the number of points by counting the number of new points // associated with each level vector. // level_1d = new int[dim_num]; point_num = 0; for ( level = 0; level <= level_max; level++ ) { more = false; h = 0; t = 0; for ( ; ;) { comp_next ( level, dim_num, level_1d, &more, &h, &t ); v = 1; for ( dim = 0; dim < dim_num; dim++ ) { v = v * new_1d[level_1d[dim]]; } point_num = point_num + v; if ( !more ) { break; } } } delete [] level_1d; delete [] new_1d; return point_num; } //****************************************************************************80 int sparse_grid_onn_size ( int dim_num, int level_max ) //****************************************************************************80 // // Purpose: // // SPARSE_GRID_ONN_SIZE sizes a sparse grid using Open Non-Nested rules. // // Discussion: // // The grid is defined as the sum of the product rules whose LEVEL // satisfies: // // 0 <= LEVEL <= LEVEL_MAX. // // This calculation is much faster than a previous method. It simply // computes the number of new points that are added at each level in the // 1D rule, and then counts the new points at a given DIM_NUM dimensional // level vector as the product of the new points added in each dimension. // // This approach will work for nested families, and may be extensible // to other families, and to mixed rules. // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 15 January 2010 // // 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. // // Parameters: // // Input, int DIM_NUM, the spatial dimension. // // Input, int LEVEL_MAX, the maximum value of LEVEL. // // Output, int SPARSE_GRID_ONN_SIZE, the number of points in the grid. // { int dim; int h; int j; int l; int level; int *level_1d; int level_min; bool more; int *order_1d; int point_num; int t; int v; // // Special case. // if ( level_max < 0 ) { point_num = 0; return point_num; } if ( level_max == 0 ) { point_num = 1; return point_num; } // // Construct the 1D order vector. // order_1d = new int[level_max+1]; for ( l = 0; l <= level_max; l++ ) { order_1d[l] = 2 * l + 1; } level_1d = new int[dim_num]; level_min = i4_max ( 0, level_max + 1 - dim_num ); point_num = 0; for ( level = level_min; level <= level_max; level++ ) { more = false; h = 0; t = 0; for ( ; ;) { comp_next ( level, dim_num, level_1d, &more, &h, &t ); v = 1; for ( dim = 0; dim < dim_num; dim++ ) { v = v * order_1d[level_1d[dim]]; } point_num = point_num + v; if ( !more ) { break; } } } delete [] level_1d; delete [] order_1d; return point_num; } //****************************************************************************80 int sparse_grid_own_size ( int dim_num, int level_max ) //****************************************************************************80 // // Purpose: // // SPARSE_GRID_OWN_SIZE sizes a sparse grid using Open Weakly Nested rules. // // Discussion: // // This calculation is much faster than a previous method. // // This calculation assumes that a linear growth rule is being used, // that is, that the 1D rules have orders 1, 3, 5, 7, 9, and so on. // // This calculation assumes that the 1D family of quadrature rules // contains only one repeated point, presumably the value 0.0. // This assumption holds for Gauss-Legendre, Gauss-Hermite and // Generalized Gauss-Hermite rules. // // The routine then counts the number of unique abscissas that will // be generated for a sparse grid of given dimension and level. // // The computation is complicated. It starts by counting just those // abscissas which have no 0.0 in them. This is relatively easy, since // it is like counting the points in a sparse grid that uses open // non-nested rules, but for which the order of each rule is reduced by 1. // // Then we have to count the abscissas with one 0.0, two 0.0's and so // on to DIM_NUM zeros. We are assuming this is an isotropic grid, // so for a particular number K of zeroes we only need to count the case // where the first K entries are zero, and multiply by C(DIM_NUM,K). // // To count the number of entries with K zeroes, (and assuming 0 < K), // then, we essentially count the number of abscissas in an open // non-nested rule as before, but modifed so that the minimum level is 0, // rather than LEVEL_MAX - DIM_NUM + 1. // // I will mention that this was a rather difficult computation to // figure out! // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 12 January 2010 // // 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. // // Parameters: // // Input, int DIM_NUM, the spatial dimension. // // Input, int LEVEL_MAX, the maximum value of LEVEL. // // Output, int SPARSE_GRID_OWN_SIZE, the number of points in the grid. // { int dim; int dim_num2; int h; int j; int l; int level; int level_min; int *level_1d; bool more; int *new_1d; int point_num; int point_num2; int t; int v; // // Special case. // if ( level_max < 0 ) { point_num = 0; return point_num; } if ( level_max == 0 ) { point_num = 1; return point_num; } // // Construct the vector that counts the new points in the 1D rule. // new_1d = new int[level_max+1]; new_1d[0] = 0; for ( l = 1; l <= level_max; l++ ) { new_1d[l] = 2 * l; } // // Count the nonzero points in the full dimensional table with the usual // LEVEL_MIN restriction. // // Then count the points with 1, 2, 3, ... DIM_NUM zeroes, by counting // the nonzero points in a DIM_NUM2 table, with LEVEL_MIN set to 0, and // multiplying by the appropriate combinatorial coefficient. // point_num = 0; for ( dim_num2 = dim_num; 0 <= dim_num2; dim_num2-- ) { if ( dim_num2 == dim_num ) { level_min = i4_max ( 0, level_max - dim_num + 1 ); } else { level_min = 0; } if ( dim_num2 == 0 ) { point_num2 = 1; } else { level_1d = new int[dim_num2]; point_num2 = 0; for ( level = level_min; level <= level_max; level++ ) { more = false; h = 0; t = 0; for ( ; ; ) { comp_next ( level, dim_num2, level_1d, &more, &h, &t ); v = 1; for ( dim = 0; dim < dim_num2; dim++ ) { v = v * new_1d[level_1d[dim]]; } point_num2 = point_num2 + v; if ( !more ) { break; } } } delete [] level_1d; } point_num = point_num + i4_choose ( dim_num, dim_num2 ) * point_num2; } delete [] new_1d; return point_num; } //****************************************************************************80 void timestamp ( ) //****************************************************************************80 // // Purpose: // // TIMESTAMP prints the current YMDHMS date as a time stamp. // // Example: // // 31 May 2001 09:45:54 AM // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 08 July 2009 // // Author: // // John Burkardt // // Parameters: // // None // { # define TIME_SIZE 40 static char time_buffer[TIME_SIZE]; const struct std::tm *tm_ptr; size_t len; std::time_t now; now = std::time ( NULL ); tm_ptr = std::localtime ( &now ); len = std::strftime ( time_buffer, TIME_SIZE, "%d %B %Y %I:%M:%S %p", tm_ptr ); std::cout << time_buffer << "\n"; return; # undef TIME_SIZE } //****************************************************************************80 void vec_colex_next2 ( int dim_num, int base[], int a[], bool *more ) //****************************************************************************80 // // Purpose: // // VEC_COLEX_NEXT2 generates vectors in colex order. // // Discussion: // // The vectors are produced in colexical order, starting with // // (0, 0, ...,0), // (1, 0, ...,0), // ... // (BASE(1)-1,0, ...,0) // // (0, 1, ...,0) // (1, 1, ...,0) // ... // (BASE(1)-1,1, ...,0) // // (0, 2, ...,0) // (1, 2, ...,0) // ... // (BASE(1)-1,BASE(2)-1,...,BASE(DIM_NUM)-1). // // Examples: // // DIM_NUM = 2, // BASE = { 3, 3 } // // 0 0 // 1 0 // 2 0 // 0 1 // 1 1 // 2 1 // 0 2 // 1 2 // 2 2 // // Licensing: // // This code is distributed under the GNU LGPL license. // // Modified: // // 25 May 2007 // // Author: // // John Burkardt // // Parameters: // // Input, int DIM_NUM, the spatial dimension. // // Input, int BASE[DIM_NUM], the bases to be used in each dimension. // In dimension I, entries will range from 0 to BASE[I]-1. // // Output, int A[DIM_NUM], the next vector. // // Input/output, bool *MORE. Set this variable false before // the first call. On return, MORE is TRUE if another vector has // been computed. If MORE is returned FALSE, ignore the output // vector and stop calling the routine. // { int i; if ( !( *more ) ) { for ( i = 0; i < dim_num; i++ ) { a[i] = 0; } *more = true; } else { for ( i = 0; i < dim_num; i++ ) { a[i] = a[i] + 1; if ( a[i] < base[i] ) { return; } a[i] = 0; } *more = false; } return; }