function sparse_grid_hermite_dataset ( dim_num, level_max ) %*****************************************************************************80 % %% SPARSE_GRID_HERMITE_DATASET is the main program for SPARSE_GRID_HERMITE_DATASET. % % Discussion: % % This program computes a sparse grid quadrature rule based on a 1D % Gauss-Hermite rule and writes it to a file.. % % The user specifies: % * the spatial dimension of the quadrature region, % * the level that defines the Smolyak grid. % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 10 August 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. % timestamp ( ); fprintf ( 1, '\n' ); fprintf ( 1, 'SPARSE_GRID_HERMITE_DATASET\n' ); fprintf ( 1, ' MATLAB version\n' ); fprintf ( 1, '\n' ); fprintf ( 1, ' Compute the abscissas and weights of a quadrature rule\n' ); fprintf ( 1, ' associated with a sparse grid derived from a Smolyak\n' ); fprintf ( 1, ' construction based on 1D Gauss-Hermite rules.\n' ); fprintf ( 1, '\n' ); fprintf ( 1, ' Inputs to the program include:\n' ); fprintf ( 1, '\n' ); fprintf ( 1, ' DIM_NUM, the spatial dimension.\n' ); fprintf ( 1, ' (typically in the range of 2 to 10)\n' ); fprintf ( 1, '\n' ); fprintf ( 1, ' LEVEL_MAX, the "level" of the sparse grid.\n' ); fprintf ( 1, ' (typically in the range of 0, 1, 2, 3, ...\n' ); fprintf ( 1, '\n' ); fprintf ( 1, ' Output from the program includes:\n' ); fprintf ( 1, '\n' ); fprintf ( 1, ' A printed table of the abscissas and weights.\n' ); fprintf ( 1, '\n' ); fprintf ( 1, ' * A set of 3 files that define the quadrature rule.\n' ); fprintf ( 1, '\n' ); fprintf ( 1, ' "herm_d?_level?_r.txt", a file of the ranges.\n' ); fprintf ( 1, ' "herm_d?_level?_w.txt", a file of the weights;\n' ); fprintf ( 1, ' "herm_d?_level?_x.txt", a file of the abscissas;\n' ); % % Get the spatial dimension. % if ( nargin < 1) fprintf ( 1, '\n' ); dim_num = input ( ' Enter the value of DIM_NUM (1 or greater)' ); end fprintf ( 1, '\n' ); fprintf ( 1, ' Spatial dimension requested is = %d\n', dim_num ); % % Get the level. % if ( nargin < 2 ) fprintf ( 1, '\n' ); level_max = input ( ' Enter the value of LEVEL_MAX (0 or greater)' ); end level_min = max ( 0, level_max + 1 - dim_num ); fprintf ( 1, '\n' ); fprintf ( 1, ' LEVEL_MIN = %d\n', level_min ); fprintf ( 1, ' LEVEL_MAX = %d\n', level_max ); % % How many distinct points will there be? % point_num = sparse_grid_herm_size ( dim_num, level_max ); fprintf ( 1, '\n' ); fprintf ( 1, ' The number of distinct abscissas in the\n' ); fprintf ( 1, ' quadrature rule is determined from the spatial\n' ); fprintf ( 1, ' dimension DIM_NUM and the level LEVEL_MAX.\n' ); fprintf ( 1, ' For the given input, this value will be = %d\n', point_num ); r = zeros ( dim_num, 2 ); % % Compute the weights and points. % r(1:dim_num,1) = - r8_huge ( ); r(1:dim_num,2) = + r8_huge ( ); [ w, x ] = sparse_grid_herm ( dim_num, level_max, point_num ); r8mat_transpose_print_some ( dim_num, point_num, x, 1, 1, ... dim_num, 10, ' First 10 grid points:' ); r8vec_print_some ( point_num, w, 1, 10, ' First 10 weights:' ); weight_sum = sum ( w(1:point_num) ); fprintf ( 1, '\n' ); fprintf ( 1, ' Weights sum to %24.16f\n', weight_sum ); fprintf ( 1, ' Correct value is %24.16f\n', sqrt ( pi^dim_num ) ); % % Write the rule to files. % r_filename = sprintf ( 'herm_d%d_level%d_r.txt', dim_num, level_max ); w_filename = sprintf ( 'herm_d%d_level%d_w.txt', dim_num, level_max ); x_filename = sprintf ( 'herm_d%d_level%d_x.txt', dim_num, level_max ); fprintf ( 1, '\n' ); fprintf ( 1, ' Creating R file = "%s".\n', r_filename ); r8mat_write ( r_filename, dim_num, 2, r ); fprintf ( 1, ' Creating W file = "%s".\n', w_filename ); r8mat_write ( w_filename, 1, point_num, w ); fprintf ( 1, ' Creating X file = "%s".\n', x_filename ); r8mat_write ( x_filename, dim_num, point_num, x ); fprintf ( 1, '\n' ); fprintf ( 1, 'SPARSE_GRID_HERM_DATASET:\n' ); fprintf ( 1, ' Normal end of execution.\n' ); fprintf ( 1, '\n' ); timestamp ( ); return end function [ a, more, h, t ] = comp_next ( n, k, a, more, h, t ) %*****************************************************************************80 % %% 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. % MATLAB 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, integer N, the integer whose compositions are desired. % % Input, integer K, the number of parts in the composition. % % Input, integer A(K), the previous composition. On the first call, % with MORE = FALSE, set A = []. Thereafter, A should be the % value of A output from the previous call. % % Input, logical MORE. The input value of MORE on the first % call should be FALSE, which tells the program to initialize. % On subsequent calls, MORE should be TRUE, or simply the % output value of MORE from the previous call. % % Input, integer H, T, two internal parameters needed for the % computation. The user may need to initialize these before the % very first call, but these initial values are not important. % The user should not alter these parameters once the computation % begins. % % Output, integer A(K), the next composition. % % Output, logical MORE, will be TRUE unless the composition % that is being returned is the final one in the sequence. % % Output, integer H, T, the updated values of the two internal % parameters. % if ( ~more ) t = n; h = 0; a(1) = n; a(2:k) = 0; else if ( 1 < t ) h = 0; end h = h + 1; t = a(h); a(h) = 0; a(1) = t - 1; a(h+1) = a(h+1) + 1; end more = ( a(k) ~= n ); return end function grid_point = hermite_abscissa ( dim_num, point_num, grid_index, ... grid_base ) %*****************************************************************************80 % %% HERMITE_ABSCISSA sets abscissas for multidimensional Gauss-Hermite quadrature. % % Discussion: % % The "nesting" as it occurs for Gauss-Hermite sparse grids simply % involves the use of a specified set of permissible orders for the % rule. % % The X array lists the (complete) Gauss-Hermite abscissas for rules % of order 1, 3, 7, 15, 31, 63 and 127, in order. % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 08 October 2007 % % Author: % % John Burkardt % % Parameters: % % Input, integer DIM_NUM, the spatial dimension. % % Input, integer POINT_NUM, the number of points. % % Input, integer GRID_INDEX(DIM_NUM,POINT_NUM), for each % point and dimension, the index of the abscissa. % % Input, integer GRID_BASE(DIM_NUM), the order of the % rule being used in each dimension. % % Output, real GRID_POINT(DIM_NUM,POINT_NUM), the grid points of abscissas. % skip = [ 0, 1, 4, 11, 26, 57, 120, 247 ]; x = [ ... 0.0E+00, ... -0.122474487139158904909864203735E+01, ... 0.0E+00, ... 0.122474487139158904909864203735E+01, ... -0.265196135683523349244708200652E+01, ... -0.167355162876747144503180139830E+01, ... -0.816287882858964663038710959027E+00, ... 0.0E+00, ... 0.816287882858964663038710959027E+00, ... 0.167355162876747144503180139830E+01, ... 0.265196135683523349244708200652E+01, ... -0.449999070730939155366438053053E+01, ... -0.366995037340445253472922383312E+01, ... -0.296716692790560324848896036355E+01, ... -0.232573248617385774545404479449E+01, ... -0.171999257518648893241583152515E+01, ... -0.113611558521092066631913490556E+01, ... -0.565069583255575748526020337198E+00, ... 0.0E+00, ... 0.565069583255575748526020337198E+00, ... 0.113611558521092066631913490556E+01, ... 0.171999257518648893241583152515E+01, ... 0.232573248617385774545404479449E+01, ... 0.296716692790560324848896036355E+01, ... 0.366995037340445253472922383312E+01, ... 0.449999070730939155366438053053E+01, ... -6.9956801237185402753248521473232E+00, ... -6.2750787049428601427036567812530E+00, ... -5.6739614446185883296332558789276E+00, ... -5.1335955771123807045862968913996E+00, ... -4.6315595063128599420667997654336E+00, ... -4.1562717558181451724831352315314E+00, ... -3.7007434032314694224497164589673E+00, ... -3.2603207323135408104645401509648E+00, ... -2.8316804533902054557015640151425E+00, ... -2.4123177054804201051740184582119E+00, ... -2.0002585489356389657975562598571E+00, ... -1.5938858604721398261388419455550E+00, ... -1.1918269983500464260821358649242E+00, ... -0.79287697691530893968593032998830E+00, ... -0.39594273647142311094670041663436E+00, ... 0.0000000000000000000000000000000E+00, ... 0.39594273647142311094670041663436E+00, ... 0.79287697691530893968593032998830E+00, ... 1.1918269983500464260821358649242E+00, ... 1.5938858604721398261388419455550E+00, ... 2.0002585489356389657975562598571E+00, ... 2.4123177054804201051740184582119E+00, ... 2.8316804533902054557015640151425E+00, ... 3.2603207323135408104645401509648E+00, ... 3.7007434032314694224497164589673E+00, ... 4.1562717558181451724831352315314E+00, ... 4.6315595063128599420667997654336E+00, ... 5.1335955771123807045862968913996E+00, ... 5.6739614446185883296332558789276E+00, ... 6.2750787049428601427036567812530E+00, ... 6.9956801237185402753248521473232E+00, ... -10.435499877854168053468115427285E+00, ... -9.8028759912974963635223935286507E+00, ... -9.2792019543050391319404745506496E+00, ... -8.8118581437284546442526628275570E+00, ... -8.3807683451863219343010651043788E+00, ... -7.9755950801420373181541806298501E+00, ... -7.5901395198641066762479783194468E+00, ... -7.2203167078889678461161324222529E+00, ... -6.8632544331795368527353285876066E+00, ... -6.5168348106821160605273395854042E+00, ... -6.1794379922705969862418461787263E+00, ... -5.8497884000810673462526582961482E+00, ... -5.5268572526403031425047575122840E+00, ... -5.2097979830408354861575136416263E+00, ... -4.8979018644975742350745099214868E+00, ... -4.5905665744435190229271294569091E+00, ... -4.2872733352824404031727616199454E+00, ... -3.9875699104197157485227052068068E+00, ... -3.6910577000963465117322810559754E+00, ... -3.3973817713303911852755941806287E+00, ... -3.1062230279282566329138616746036E+00, ... -2.8172919672837977750747135657355E+00, ... -2.5303236304712010926855221718499E+00, ... -2.2450734604812066298995918179330E+00, ... -1.9613138583081485293922008411321E+00, ... -1.6788312791720137520802800622638E+00, ... -1.3974237486049625107570752063702E+00, ... -1.1168987050996462690510970277840E+00, ... -0.83707109558947615977737795461293E+00, ... -0.55776166427908221668763665253822E+00, ... -0.27879538567115223986687628627202E+00, ... 0.00000000000000000000000000000000E+00, ... 0.27879538567115223986687628627202E+00, ... 0.55776166427908221668763665253822E+00, ... 0.83707109558947615977737795461293E+00, ... 1.1168987050996462690510970277840E+00, ... 1.3974237486049625107570752063702E+00, ... 1.6788312791720137520802800622638E+00, ... 1.9613138583081485293922008411321E+00, ... 2.2450734604812066298995918179330E+00, ... 2.5303236304712010926855221718499E+00, ... 2.8172919672837977750747135657355E+00, ... 3.1062230279282566329138616746036E+00, ... 3.3973817713303911852755941806287E+00, ... 3.6910577000963465117322810559754E+00, ... 3.9875699104197157485227052068068E+00, ... 4.2872733352824404031727616199454E+00, ... 4.5905665744435190229271294569091E+00, ... 4.8979018644975742350745099214868E+00, ... 5.2097979830408354861575136416263E+00, ... 5.5268572526403031425047575122840E+00, ... 5.8497884000810673462526582961482E+00, ... 6.1794379922705969862418461787263E+00, ... 6.5168348106821160605273395854042E+00, ... 6.8632544331795368527353285876066E+00, ... 7.2203167078889678461161324222529E+00, ... 7.5901395198641066762479783194468E+00, ... 7.9755950801420373181541806298501E+00, ... 8.3807683451863219343010651043788E+00, ... 8.8118581437284546442526628275570E+00, ... 9.2792019543050391319404745506496E+00, ... 9.8028759912974963635223935286507E+00, ... 10.435499877854168053468115427285E+00, ... -15.228338148167350978246954433464E+00, ... -14.669595158833972632746354112896E+00, ... -14.209085995284870755168244250887E+00, ... -13.799722290211676634645246746673E+00, ... -13.423518590070950062438258321855E+00, ... -13.071208660474601901583995439649E+00, ... -12.737235652415686338138003924072E+00, ... -12.417939378869715805445879624069E+00, ... -12.110749020947747600132123508132E+00, ... -11.813772198267727195134584136191E+00, ... -11.525565112572696599167888588564E+00, ... -11.244994583785543445194384194300E+00, ... -10.971150569840247423423040263881E+00, ... -10.703288201027481347670940744690E+00, ... -10.440787957772772867742591798027E+00, ... -10.183127473450343888624126450357E+00, ... -9.9298610495114250736847004273684E+00, ... -9.6806044412474728038150712732737E+00, ... -9.4350233389881650135019598506287E+00, ... -9.1928244988460305715774195052527E+00, ... -8.9537488108565404323807890169970E+00, ... -8.7175658087076307363833999548548E+00, ... -8.4840692689832473326097180339984E+00, ... -8.2530736454457156579694124243888E+00, ... -8.0244111514703375578594739796798E+00, ... -7.7979293513870105420829120455591E+00, ... -7.5734891556083454022834960763301E+00, ... -7.3509631392269052701961258043733E+00, ... -7.1302341220350710668064025713431E+00, ... -6.9111939615465713197465633109366E+00, ... -6.6937425208758294190074417381666E+00, ... -6.4777867811645365448144903821487E+00, ... -6.2632400742737354345609723857092E+00, ... -6.0500214161419845694465474482388E+00, ... -5.8380549248774187386601690807757E+00, ... -5.6272693105464816659423455794909E+00, ... -5.4175974259243240722848425872924E+00, ... -5.2089758693153983587570258372239E+00, ... -5.0013446320386360038520809107373E+00, ... -4.7946467843764925009748509930857E+00, ... -4.5888281947698372951606485031212E+00, ... -4.3838372778464736294253744407459E+00, ... -4.1796247675352031349421189892408E+00, ... -3.9761435120673355916035814195920E+00, ... -3.7733482881250526721004678400057E+00, ... -3.5711956317782180447199756485249E+00, ... -3.3696436841717397896643629240035E+00, ... -3.1686520501953630191857798261495E+00, ... -2.9681816685955910267761649521505E+00, ... -2.7681946921824058801226545958892E+00, ... -2.5686543769473501723144013022363E+00, ... -2.3695249790490401080012474645702E+00, ... -2.1707716587411506879498498083695E+00, ... -1.9723603904195020079324743227565E+00, ... -1.7742578780516791584676442103681E+00, ... -1.5764314753267801315519597621879E+00, ... -1.3788491099261778091441557053728E+00, ... -1.1814792113700685848678583598423E+00, ... -0.98429064194027277726568984213773E+00, ... -0.78725263021825034151596831878971E+00, ... -0.59033470680942102142230439346102E+00, ... -0.39350664185130136568037826200185E+00, ... -0.19673838392423251964272239737078E+00, ... 0.0000000000000000000000000000000E+00, ... 0.19673838392423251964272239737078E+00, ... 0.39350664185130136568037826200185E+00, ... 0.59033470680942102142230439346102E+00, ... 0.78725263021825034151596831878971E+00, ... 0.98429064194027277726568984213773E+00, ... 1.1814792113700685848678583598423E+00, ... 1.3788491099261778091441557053728E+00, ... 1.5764314753267801315519597621879E+00, ... 1.7742578780516791584676442103681E+00, ... 1.9723603904195020079324743227565E+00, ... 2.1707716587411506879498498083695E+00, ... 2.3695249790490401080012474645702E+00, ... 2.5686543769473501723144013022363E+00, ... 2.7681946921824058801226545958892E+00, ... 2.9681816685955910267761649521505E+00, ... 3.1686520501953630191857798261495E+00, ... 3.3696436841717397896643629240035E+00, ... 3.5711956317782180447199756485249E+00, ... 3.7733482881250526721004678400057E+00, ... 3.9761435120673355916035814195920E+00, ... 4.1796247675352031349421189892408E+00, ... 4.3838372778464736294253744407459E+00, ... 4.5888281947698372951606485031212E+00, ... 4.7946467843764925009748509930857E+00, ... 5.0013446320386360038520809107373E+00, ... 5.2089758693153983587570258372239E+00, ... 5.4175974259243240722848425872924E+00, ... 5.6272693105464816659423455794909E+00, ... 5.8380549248774187386601690807757E+00, ... 6.0500214161419845694465474482388E+00, ... 6.2632400742737354345609723857092E+00, ... 6.4777867811645365448144903821487E+00, ... 6.6937425208758294190074417381666E+00, ... 6.9111939615465713197465633109366E+00, ... 7.1302341220350710668064025713431E+00, ... 7.3509631392269052701961258043733E+00, ... 7.5734891556083454022834960763301E+00, ... 7.7979293513870105420829120455591E+00, ... 8.0244111514703375578594739796798E+00, ... 8.2530736454457156579694124243888E+00, ... 8.4840692689832473326097180339984E+00, ... 8.7175658087076307363833999548548E+00, ... 8.9537488108565404323807890169970E+00, ... 9.1928244988460305715774195052527E+00, ... 9.4350233389881650135019598506287E+00, ... 9.6806044412474728038150712732737E+00, ... 9.9298610495114250736847004273684E+00, ... 10.183127473450343888624126450357E+00, ... 10.440787957772772867742591798027E+00, ... 10.703288201027481347670940744690E+00, ... 10.971150569840247423423040263881E+00, ... 11.244994583785543445194384194300E+00, ... 11.525565112572696599167888588564E+00, ... 11.813772198267727195134584136191E+00, ... 12.110749020947747600132123508132E+00, ... 12.417939378869715805445879624069E+00, ... 12.737235652415686338138003924072E+00, ... 13.071208660474601901583995439649E+00, ... 13.423518590070950062438258321855E+00, ... 13.799722290211676634645246746673E+00, ... 14.209085995284870755168244250887E+00, ... 14.669595158833972632746354112896E+00, ... 15.228338148167350978246954433464E+00 ... ]; if ( any ( grid_base(1:dim_num) < 0 ) ) fprintf ( 1, '\n' ); fprintf ( 1, 'HERMITE_ABSCISSA - Fatal error!\n' ); fprintf ( 1, ' Some base values are less than 0.\n' ); error ( 'HERMITE_ABSCISSA - Fatal error!' ); end if ( any ( 63 < grid_base(1:dim_num) ) ) fprintf ( 1, '\n' ); fprintf ( 1, 'HERMITE_ABSCISSA - Fatal error!\n' ); fprintf ( 1, ' Some base values are greater than 63.\n' ); error ( 'HERMITE_ABSCISSA - Fatal error!' ); end grid_point = zeros ( dim_num, point_num ); for point = 1 : point_num for dim = 1 : dim_num level = i4_log_2 ( grid_base(dim) + 1 ); pointer = skip(level+1) + ( grid_index(dim,point) + grid_base(dim) + 1 ); grid_point(dim,point) = x(pointer); end end return end function weight = hermite_weights ( order ) %*****************************************************************************80 % %% HERMITE_WEIGHTS returns weights for certain Gauss-Hermite quadrature rules. % % Discussion: % % The allowed orders are 1, 3, 7, 15, 31, 63 or 127. % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 08 October 2007 % % Author: % % John Burkardt % % Reference: % % Milton Abramowitz, Irene Stegun, % Handbook of Mathematical Functions, % National Bureau of Standards, 1964, % ISBN: 0-486-61272-4, % LC: QA47.A34. % % Arthur Stroud, Don Secrest, % Gaussian Quadrature Formulas, % Prentice Hall, 1966, % LC: QA299.4G3S7. % % Parameters: % % Input, integer ORDER, the order of the rule. % ORDER must be 1, 3, 7, 15, 31, 63, or 127. % % Output, real WEIGHT(1,ORDER), the weights. % The weights are positive, symmetric and should sum to sqrt(PI). % weight = zeros ( 1, order ); if ( order == 1 ) weight(1) = 1.77245385090551602729816748334E+00; elseif ( order == 3 ) weight(1) = 0.295408975150919337883027913890E+00; weight(2) = 0.118163590060367735153211165556E+01; weight(3) = 0.295408975150919337883027913890E+00; elseif ( order == 7 ) weight(1) = 0.971781245099519154149424255939E-03; weight(2) = 0.545155828191270305921785688417E-01; weight(3) = 0.425607252610127800520317466666E+00; weight(4) = 0.810264617556807326764876563813E+00; weight(5) = 0.425607252610127800520317466666E+00; weight(6) = 0.545155828191270305921785688417E-01; weight(7) = 0.971781245099519154149424255939E-03; elseif ( order == 15 ) weight(1) = 0.152247580425351702016062666965E-08; weight(2) = 0.105911554771106663577520791055E-05; weight(3) = 0.100004441232499868127296736177E-03; weight(4) = 0.277806884291277589607887049229E-02; weight(5) = 0.307800338725460822286814158758E-01; weight(6) = 0.158488915795935746883839384960E+00; weight(7) = 0.412028687498898627025891079568E+00; weight(8) = 0.564100308726417532852625797340E+00; weight(9) = 0.412028687498898627025891079568E+00; weight(10) = 0.158488915795935746883839384960E+00; weight(11) = 0.307800338725460822286814158758E-01; weight(12) = 0.277806884291277589607887049229E-02; weight(13) = 0.100004441232499868127296736177E-03; weight(14) = 0.105911554771106663577520791055E-05; weight(15) = 0.152247580425351702016062666965E-08; elseif ( order == 31 ) weight( 1) = 0.46189683944498305857470556847735E-21; weight( 2) = 0.51106090079112519643027197715274E-17; weight( 3) = 0.58995564987355133075257722133966E-14; weight( 4) = 0.18603735214463569590294465062239E-11; weight( 5) = 0.23524920032013205739850619940094E-09; weight( 6) = 0.14611988344865057576066495091513E-07; weight( 7) = 0.50437125589241034841778074689627E-06; weight( 8) = 0.10498602757642934202945441341697E-04; weight( 9) = 0.13952090395003623854995664958146E-03; weight( 10) = 0.12336833073030489880608311394968E-02; weight( 11) = 0.74827999140119116765002499116934E-02; weight( 12) = 0.31847230731201222775249585776902E-01; weight( 13) = 0.96717948160569462991143316029341E-01; weight( 14) = 0.21213278866810461318136114862419E+00; weight( 15) = 0.33877265789305344906000174083214E+00; weight( 16) = 0.39577855609737786462923720809676E+00; weight( 17) = 0.33877265789305344906000174083214E+00; weight( 18) = 0.21213278866810461318136114862419E+00; weight( 19) = 0.96717948160569462991143316029341E-01; weight( 20) = 0.31847230731201222775249585776902E-01; weight( 21) = 0.74827999140119116765002499116934E-02; weight( 22) = 0.12336833073030489880608311394968E-02; weight( 23) = 0.13952090395003623854995664958146E-03; weight( 24) = 0.10498602757642934202945441341697E-04; weight( 25) = 0.50437125589241034841778074689627E-06; weight( 26) = 0.14611988344865057576066495091513E-07; weight( 27) = 0.23524920032013205739850619940094E-09; weight( 28) = 0.18603735214463569590294465062239E-11; weight( 29) = 0.58995564987355133075257722133966E-14; weight( 30) = 0.51106090079112519643027197715274E-17; weight( 31) = 0.46189683944498305857470556847735E-21; elseif ( order == 63 ) weight( 1) = 0.37099206434787551197827130470031E-47; weight( 2) = 0.10400778615192299534481914814892E-41; weight( 3) = 0.19796804708258311251124226474396E-37; weight( 4) = 0.84687478191640015120141181138947E-34; weight( 5) = 0.13071305930779945903630127634063E-30; weight( 6) = 0.93437837175367456929765381518998E-28; weight( 7) = 0.36027426635173044862245783257252E-25; weight( 8) = 0.82963863115951789374753323156164E-23; weight( 9) = 0.12266629909105281472971700203949E-20; weight( 10) = 0.12288435628797061539461585325494E-18; weight( 11) = 0.86925536958188009075932426691516E-17; weight( 12) = 0.44857058689176221240330804981619E-15; weight( 13) = 0.17335817955735154599902643794700E-13; weight( 14) = 0.51265062385038307838565047455223E-12; weight( 15) = 0.11808921844532942490513037158404E-10; weight( 16) = 0.21508698297808025739828859845140E-09; weight( 17) = 0.31371929535285447801497640621672E-08; weight( 18) = 0.37041625984781705796752840204084E-07; weight( 19) = 0.35734732949879669663960738150956E-06; weight( 20) = 0.28393114498380927832990899215541E-05; weight( 21) = 0.18709113003730498008961134765721E-04; weight( 22) = 0.10284880800653635546698378640623E-03; weight( 23) = 0.47411702610173128107201781718693E-03; weight( 24) = 0.18409222622384813438539657470055E-02; weight( 25) = 0.60436044551187631655712178246467E-02; weight( 26) = 0.16829299199599730926458559757600E-01; weight( 27) = 0.39858264027692992170237391875317E-01; weight( 28) = 0.80467087993950415219587554532823E-01; weight( 29) = 0.13871950817615293377792092082674E+00; weight( 30) = 0.20448695346833761570957197160475E+00; weight( 31) = 0.25799889943058042204920467417642E+00; weight( 32) = 0.27876694884838411919175686949858E+00; weight( 33) = 0.25799889943058042204920467417642E+00; weight( 34) = 0.20448695346833761570957197160475E+00; weight( 35) = 0.13871950817615293377792092082674E+00 ; weight( 36) = 0.80467087993950415219587554532823E-01; weight( 37) = 0.39858264027692992170237391875317E-01; weight( 38) = 0.16829299199599730926458559757600E-01; weight( 39) = 0.60436044551187631655712178246467E-02; weight( 40) = 0.18409222622384813438539657470055E-02; weight( 41) = 0.47411702610173128107201781718693E-03; weight( 42) = 0.10284880800653635546698378640623E-03; weight( 43) = 0.18709113003730498008961134765721E-04; weight( 44) = 0.28393114498380927832990899215541E-05; weight( 45) = 0.35734732949879669663960738150956E-06; weight( 46) = 0.37041625984781705796752840204084E-07; weight( 47) = 0.31371929535285447801497640621672E-08; weight( 48) = 0.21508698297808025739828859845140E-09; weight( 49) = 0.11808921844532942490513037158404E-10; weight( 50) = 0.51265062385038307838565047455223E-12; weight( 51) = 0.17335817955735154599902643794700E-13; weight( 52) = 0.44857058689176221240330804981619E-15; weight( 53) = 0.86925536958188009075932426691516E-17; weight( 54) = 0.12288435628797061539461585325494E-18; weight( 55) = 0.12266629909105281472971700203949E-20; weight( 56) = 0.82963863115951789374753323156164E-23; weight( 57) = 0.36027426635173044862245783257252E-25; weight( 58) = 0.93437837175367456929765381518998E-28; weight( 59) = 0.13071305930779945903630127634063E-30; weight( 60) = 0.84687478191640015120141181138947E-34; weight( 61) = 0.19796804708258311251124226474396E-37; weight( 62) = 0.10400778615192299534481914814892E-41; weight( 63) = 0.37099206434787551197827130470031E-47; elseif ( order == 127 ) weight( 1) = 0.12504497577050595552677230002883E-100; weight( 2) = 0.17272798059419131415318615789672E-93; weight( 3) = 0.89321681571986548608031150791499E-88; weight( 4) = 0.77306185240893578449625186483810E-83; weight( 5) = 0.20143957652648255497735460506196E-78; weight( 6) = 0.21503714733610239701351039429345E-74; weight( 7) = 0.11341924208594594813715533569504E-70; weight( 8) = 0.33489139011795051950683388483136E-67; weight( 9) = 0.60486548964016681064424451668405E-64; weight( 10) = 0.71375092946352177824971347343892E-61; weight( 11) = 0.57884563374885556636801095624030E-58; weight( 12) = 0.33581166223858230300409326551248E-55; weight( 13) = 0.14394641949253923568603163698953E-52; weight( 14) = 0.46821808383216117724080263903889E-50; weight( 15) = 0.11817054440684264071348471955361E-47; weight( 16) = 0.23581659156008927203181682045005E-45; weight( 17) = 0.37814427940797540210712758405540E-43; weight( 18) = 0.49411031115771638145610738414006E-41; weight( 19) = 0.53255303775425059266087298458297E-39; weight( 20) = 0.47854390680131484999315199332765E-37; weight( 21) = 0.36191883445952356128627543209554E-35; weight( 22) = 0.23232083386343554805352497446119E-33; weight( 23) = 0.12753331411008716683688974281454E-31; weight( 24) = 0.60277753850758742112436095241270E-30; weight( 25) = 0.24679773241777200207460855084439E-28; weight( 26) = 0.88019567691698482573264198727415E-27; weight( 27) = 0.27482489212040561315005725890593E-25; weight( 28) = 0.75468218903085486125222816438456E-24; weight( 29) = 0.18303134636280466270545996891835E-22; weight( 30) = 0.39355990860860813085582448449811E-21; weight( 31) = 0.75293161638581191068419292570042E-20; weight( 32) = 0.12857997786722855037584105682618E-18; weight( 33) = 0.19659326888445857792541925311450E-17; weight( 34) = 0.26986511907214101894995783364250E-16; weight( 35) = 0.33344414303198856330118301113874E-15; weight( 36) = 0.37173303125150639885726463109574E-14; weight( 37) = 0.37473954472839737091885387788983E-13; weight( 38) = 0.34230094493397259538669512076007E-12; weight( 39) = 0.28385303724993373166810860630552E-11; weight( 40) = 0.21406920290454669208938772802828E-10; weight( 41) = 0.14706331273431716244229273183839E-09; weight( 42) = 0.92173940967434659264335883218167E-09; weight( 43) = 0.52781663936972714041837056042506E-08; weight( 44) = 0.27650497044951117835905283127679E-07; weight( 45) = 0.13267855842539464770913063113371E-06; weight( 46) = 0.58380944276113062188573331195042E-06; weight( 47) = 0.23581561724775629112332165335800E-05; weight( 48) = 0.87524468034280444703919485644809E-05; weight( 49) = 0.29876790535909012274846532159647E-04; weight( 50) = 0.93874435720072545206729594267039E-04; weight( 51) = 0.27170762627931172053444716883938E-03; weight( 52) = 0.72493929742498358979684249380921E-03; weight( 53) = 0.17841208326763432884316727108264E-02; weight( 54) = 0.40524855186046131499765636276283E-02; weight( 55) = 0.85000263041544110385806705526917E-02; weight( 56) = 0.16471142241609687824005585301760E-01; weight( 57) = 0.29499296248213632269675010319119E-01; weight( 58) = 0.48847387114300011006959603975676E-01; weight( 59) = 0.74807989768583731416517226905270E-01; weight( 60) = 0.10598520508090929403834368934301E+00; weight( 61) = 0.13893945309051540832066283010510E+00; weight( 62) = 0.16856236074207929740526975049765E+00; weight( 63) = 0.18927849580120432177170145550076E+00; weight( 64) = 0.19673340688823289786163676995151E+00; weight( 65) = 0.18927849580120432177170145550076E+00; weight( 66) = 0.16856236074207929740526975049765E+00; weight( 67) = 0.13893945309051540832066283010510E+00; weight( 68) = 0.10598520508090929403834368934301E+00; weight( 69) = 0.74807989768583731416517226905270E-01; weight( 70) = 0.48847387114300011006959603975676E-01; weight( 71) = 0.29499296248213632269675010319119E-01; weight( 72) = 0.16471142241609687824005585301760E-01; weight( 73) = 0.85000263041544110385806705526917E-02; weight( 74) = 0.40524855186046131499765636276283E-02; weight( 75) = 0.17841208326763432884316727108264E-02; weight( 76) = 0.72493929742498358979684249380921E-03; weight( 77) = 0.27170762627931172053444716883938E-03; weight( 78) = 0.93874435720072545206729594267039E-04; weight( 79) = 0.29876790535909012274846532159647E-04; weight( 80) = 0.87524468034280444703919485644809E-05; weight( 81) = 0.23581561724775629112332165335800E-05; weight( 82) = 0.58380944276113062188573331195042E-06; weight( 83) = 0.13267855842539464770913063113371E-06; weight( 84) = 0.27650497044951117835905283127679E-07; weight( 85) = 0.52781663936972714041837056042506E-08; weight( 86) = 0.92173940967434659264335883218167E-09; weight( 87) = 0.14706331273431716244229273183839E-09; weight( 88) = 0.21406920290454669208938772802828E-10; weight( 89) = 0.28385303724993373166810860630552E-11; weight( 90) = 0.34230094493397259538669512076007E-12; weight( 91) = 0.37473954472839737091885387788983E-13; weight( 92) = 0.37173303125150639885726463109574E-14; weight( 93) = 0.33344414303198856330118301113874E-15; weight( 94) = 0.26986511907214101894995783364250E-16; weight( 95) = 0.19659326888445857792541925311450E-17; weight( 96) = 0.12857997786722855037584105682618E-18; weight( 97) = 0.75293161638581191068419292570042E-20; weight( 98) = 0.39355990860860813085582448449811E-21; weight( 99) = 0.18303134636280466270545996891835E-22; weight(100) = 0.75468218903085486125222816438456E-24; weight(101) = 0.27482489212040561315005725890593E-25; weight(102) = 0.88019567691698482573264198727415E-27; weight(103) = 0.24679773241777200207460855084439E-28; weight(104) = 0.60277753850758742112436095241270E-30; weight(105) = 0.12753331411008716683688974281454E-31; weight(106) = 0.23232083386343554805352497446119E-33; weight(107) = 0.36191883445952356128627543209554E-35; weight(108) = 0.47854390680131484999315199332765E-37; weight(109) = 0.53255303775425059266087298458297E-39; weight(110) = 0.49411031115771638145610738414006E-41; weight(111) = 0.37814427940797540210712758405540E-43; weight(112) = 0.23581659156008927203181682045005E-45; weight(113) = 0.11817054440684264071348471955361E-47; weight(114) = 0.46821808383216117724080263903889E-50; weight(115) = 0.14394641949253923568603163698953E-52; weight(116) = 0.33581166223858230300409326551248E-55; weight(117) = 0.57884563374885556636801095624030E-58; weight(118) = 0.71375092946352177824971347343892E-61; weight(119) = 0.60486548964016681064424451668405E-64; weight(120) = 0.33489139011795051950683388483136E-67; weight(121) = 0.11341924208594594813715533569504E-70; weight(122) = 0.21503714733610239701351039429345E-74; weight(123) = 0.20143957652648255497735460506196E-78; weight(124) = 0.77306185240893578449625186483810E-83; weight(125) = 0.89321681571986548608031150791499E-88; weight(126) = 0.17272798059419131415318615789672E-93; weight(127) = 0.12504497577050595552677230002883E-100; else fprintf ( 1, '\n' ); fprintf ( 1, 'HERMITE_WEIGHTS - Fatal error!\n' ); fprintf ( 1, ' Illegal value of ORDER = %d\n', order ); fprintf ( 1, ' Legal values are 1, 3, 7, 15, 31, 63 or 127.\n' ); error ( 'HERMITE_WEIGHTS - Fatal error!' ); end return end function value = i4_choose ( n, k ) %*****************************************************************************80 % %% 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, integer N, K, are the values of N and K. % % Output, integer VALUE, the number of combinations of N % things taken K at a time. % mn = min ( k, n - k ); if ( mn < 0 ) value = 0; elseif ( mn == 0 ) value = 1; else mx = max ( k, n - k ); value = mx + 1; for i = 2 : mn value = ( value * ( mx + i ) ) / i; end end return end function value = i4_log_2 ( i ) %*****************************************************************************80 % %% I4_LOG_2 returns the integer part of the logarithm base 2 of |I|. % % Discussion: % % For positive I4_LOG_2(I), it should be true that % 2**I4_LOG_2(X) <= |I| < 2**(I4_LOG_2(I)+1). % The special case of I4_LOG_2(0) returns -HUGE(). % % Example: % % I Value % % 0 -1 % 1, 0 % 2, 1 % 3, 1 % 4, 2 % 5, 2 % 6, 2 % 7, 2 % 8, 3 % 9, 3 % 10, 3 % 127, 6 % 128, 7 % 129, 7 % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 09 April 2004 % % Author: % % John Burkardt % % Parameters: % % Input, integer I, the number whose logarithm base 2 is desired. % % Output, integer VALUE, the integer part of the logarithm base 2 of % the absolute value of I. % i = floor ( i ); if ( i == 0 ) value = -inf; else value = 0; i = abs ( i ); while ( 2 <= i ) i = i / 2; value = value + 1; end end return end function grid_level = index_level_herm ( level, level_max, dim_num, ... point_num, grid_index, grid_base ) %*****************************************************************************80 % %% INDEX_LEVEL_HERM: determine first level at which given index is generated. % % Discussion: % % We are constructing a sparse grid of Gauss-Hermite points. The grid % is built up of product grids, with a characteristic LEVEL. % % We are concerned with identifying points in this product grid which % have actually been generated previously, on a lower value of LEVEL. % % This routine determines the lowest value of LEVEL at which each of % the input points would be generated. % % In 1D, given LEVEL, the number of points is ORDER = 2**(LEVEL+1) + 1, % (except that LEVEL = 0 implies ORDER = 1), the BASE is (ORDER-1)/2, % and the point INDEX values range from -BASE to +BASE. % % The values of INDEX and BASE allow us to determine the abstract % properties of the point. In particular, if INDEX is 0, the corresponding % Gauss-Hermite abscissa is 0, the special "nested" value we need % to take care of. % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 08 October 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, integer LEVEL, the level at which these points were % generated. LEVEL_MIN <= LEVEL <= LEVEL_MAX. % % Input, integer LEVEL_MAX, the maximum level. % % Input, integer DIM_NUM, the spatial dimension. % % Input, integer POINT_NUM, the number of points to be tested. % % Input, integer GRID_INDEX(DIM_NUM,POINT_NUM), the indices of the % points to be tested. % % Input, integer GRID_BASE(DIM_NUM), the "base", which is essentially % the denominator of the index. % % Output, integer GRID_LEVEL(POINT_NUM,1), the value of LEVEL at % which the point would first be generated. This will be the same as % the input value of LEVEL, unless the point has an INDEX of 0 and % a corresponding BASE that is NOT zero. % level_min = max ( 0, level_max + 1 - dim_num ); % % If a point has a DIM-th component whose INDEX is 0, then the % value of LEVEL at which this point would first be generated is % less than LEVEL, unless the DIM-th component of GRID_BASE is 0. % grid_level = zeros ( point_num, 1 ); for point = 1 : point_num grid_level(point) = max ( level, level_min ); for dim = 1 : dim_num if ( grid_index(dim,point) == 0 ) grid_level(point) = max ( grid_level(point) - grid_base(dim), level_min ); end end end return end function order = level_to_order_open ( dim_num, level ) %*****************************************************************************80 % %% LEVEL_TO_ORDER 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 1, % Fejer Type 2, Newton Cotes Open, Newton Cotes Half 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 1, Fejer Type 2, Newton Cotes Open, % Newton Cotes Open Half, 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: % % 18 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, integer DIM_NUM, the spatial dimension. % % Input, integer LEVEL(DIM_NUM), the nesting level. % % Output, integer ORDER(DIM_NUM,1), the order (number of points) of the rule. % order = zeros ( dim_num, 1 ); for dim = 1 : dim_num if ( level(dim) < 0 ) order(dim) = -1; elseif ( level(dim) == 0 ) order(dim) = 1; else order(dim) = 2^( level(dim) + 1 ) - 1; end end return end function indx = multigrid_index_z ( dim_num, order_1d, order_nd ) %*****************************************************************************80 % %% MULTIGRID_INDEX_Z returns an indexed multidimensional grid. % % Discussion: % % For dimension DIM, the number of points is ORDER_1D(DIM). % % We assume that ORDER_1D(DIM) is an odd number, % ORDER_1D(DIM) = N = 2 * M + 1 % so that the points have coordinates % -M/M, -(M-1)/M, ..., -1/M, 0/M, 1/M, 2/M, 3/M, ..., (M-1)/M, M/M. % and we index them as % -M, -(M-1), -1, 0, 1, 2, 3, ..., M-1, M. % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 05 February 2011 % % 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, integer DIM_NUM, the spatial dimension of the points. % % Input, integer ORDER_1D(DIM_NUM), the order of the % rule in each dimension. % % Input, integer ORDER_ND, the product of the entries of ORDER_1D. % % Output, integer 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. % indx = zeros ( dim_num, order_nd ); a = zeros ( dim_num, 1 ); more = 0; p = 0; while ( 1 ) [ a, more ] = vec_colex_next2 ( dim_num, order_1d, a, more ); if ( ~more ) break end p = p + 1; % % The values of A(DIM) are between 0 and ORDER_1D(DIM)-1 = N - 1 = 2 * M. % Subtracting M sets the range to -M to +M, as we wish. % indx(1:dim_num,p) = a(1:dim_num) - round ( ( order_1d(1:dim_num) - 1 ) / 2 ); end return end function w_nd = product_weight_herm ( dim_num, order_1d, order_nd ) %*****************************************************************************80 % %% PRODUCT_WEIGHT_HERM: weights for a product Gauss-Hermite rule. % % Discussion: % % This routine computes the weights for a quadrature rule which is % a product of 1D Gauss-Hermite rules of varying order. % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 08 October 2007 % % Author: % % John Burkardt % % Parameters: % % Input, integer DIM_NUM, the spatial dimension. % % Input, integer ORDER_1D(DIM_NUM), the order of the 1D rules. % % Input, integer ORDER_ND, the order of the product rule. % % Output, real W_ND(1,ORDER_ND), the product rule weights. % w_nd = ones ( 1, order_nd ); for dim = 1 : dim_num w_1d = hermite_weights ( order_1d(dim) ); w_nd = r8vec_direct_product2 ( dim, order_1d(dim), w_1d, dim_num, ... order_nd, w_nd ); end return end function value = r8_huge ( ) %*****************************************************************************80 % %% R8_HUGE returns a "huge" real number. % % 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. % % MATLAB provides a built-in symbolic constant "inf" that can be used % if a huge number is really what you want! % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 27 January 2008 % % Author: % % John Burkardt % % Parameters: % % Output, real VALUE, a huge number. % value = 1.0E+30; return end function r8mat_transpose_print_some ( m, n, a, ilo, jlo, ihi, jhi, title ) %*****************************************************************************80 % %% R8MAT_TRANSPOSE_PRINT_SOME prints some of an R8MAT, transposed. % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 23 May 2005 % % Author: % % John Burkardt % % Parameters: % % Input, integer M, N, the number of rows and columns. % % Input, real A(M,N), an M by N matrix to be printed. % % Input, integer ILO, JLO, the first row and column to print. % % Input, integer IHI, JHI, the last row and column to print. % % Input, string TITLE, a title. % incx = 5; fprintf ( 1, '\n' ); fprintf ( 1, '%s\n', title ); for i2lo = max ( ilo, 1 ) : incx : min ( ihi, m ) i2hi = i2lo + incx - 1; i2hi = min ( i2hi, m ); i2hi = min ( i2hi, ihi ); inc = i2hi + 1 - i2lo; fprintf ( 1, '\n' ); fprintf ( 1, ' Row: ' ); for i = i2lo : i2hi fprintf ( 1, '%7d ', i ); end fprintf ( 1, '\n' ); fprintf ( 1, ' Col\n' ); j2lo = max ( jlo, 1 ); j2hi = min ( jhi, n ); for j = j2lo : j2hi fprintf ( 1, '%5d ', j ); for i2 = 1 : inc i = i2lo - 1 + i2; fprintf ( 1, '%12f', a(i,j) ); end fprintf ( 1, '\n' ); end end return end function r8mat_write ( output_filename, m, n, table ) %*****************************************************************************80 % %% R8MAT_WRITE writes an R8MAT file. % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 09 August 2009 % % Author: % % John Burkardt % % Parameters: % % Input, string OUTPUT_FILENAME, the output filename. % % Input, integer M, the spatial dimension. % % Input, integer N, the number of points. % % Input, real TABLE(M,N), the points. % % % Open the file. % output_unit = fopen ( output_filename, 'wt' ); if ( output_unit < 0 ) fprintf ( 1, '\n' ); fprintf ( 1, 'R8MAT_WRITE - Error!\n' ); fprintf ( 1, ' Could not open the output file.\n' ); error ( 'R8MAT_WRITE - Error!' ); end % % Write the data. % for j = 1 : n for i = 1 : m fprintf ( output_unit, ' %24.16f', table(i,j) ); end fprintf ( output_unit, '\n' ); end % % Close the file. % fclose ( output_unit ); return end function w = r8vec_direct_product2 ( factor_index, factor_order, ... factor_value, factor_num, point_num, w ) %*****************************************************************************80 % %% R8VEC_DIRECT_PRODUCT2 creates a direct product of R8VEC's. % % Discussion: % % To explain what is going on here, suppose we had to construct % a multidimensional quadrature rule as the product of K rules % for 1D quadrature. % % The product rule will be represented as a list of points and weights. % % The J-th item in the product rule will be associated with % item J1 of 1D rule 1, % item J2 of 1D rule 2, % ..., % item JK of 1D rule K. % % In particular, % X(J) = ( X(1,J1), X(2,J2), ..., X(K,JK)) % and % W(J) = W(1,J1) * W(2,J2) * ... * W(K,JK) % % So we can construct the quadrature rule if we can properly % distribute the information in the 1D quadrature rules. % % This routine carries out that task for the weights W. % % Another way to do this would be to compute, one by one, the % set of all possible indices (J1,J2,...,JK), and then index % the appropriate information. An advantage of the method shown % here is that you can process the K-th set of information and % then discard it. % % Example: % % Rule 1: % Order = 4 % W(1:4) = ( 2, 3, 5, 7 ) % % Rule 2: % Order = 3 % W(1:3) = ( 11, 13, 17 ) % % Rule 3: % Order = 2 % W(1:2) = ( 19, 23 ) % % Product Rule: % Order = 24 % W(1:24) = % ( 2 * 11 * 19 ) % ( 3 * 11 * 19 ) % ( 4 * 11 * 19 ) % ( 7 * 11 * 19 ) % ( 2 * 13 * 19 ) % ( 3 * 13 * 19 ) % ( 5 * 13 * 19 ) % ( 7 * 13 * 19 ) % ( 2 * 17 * 19 ) % ( 3 * 17 * 19 ) % ( 5 * 17 * 19 ) % ( 7 * 17 * 19 ) % ( 2 * 11 * 23 ) % ( 3 * 11 * 23 ) % ( 5 * 11 * 23 ) % ( 7 * 11 * 23 ) % ( 2 * 13 * 23 ) % ( 3 * 13 * 23 ) % ( 5 * 13 * 23 ) % ( 7 * 13 * 23 ) % ( 2 * 17 * 23 ) % ( 3 * 17 * 23 ) % ( 5 * 17 * 23 ) % ( 7 * 17 * 23 ) % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 18 April 2009 % % Author: % % John Burkardt % % Parameters: % % Input, integer FACTOR_INDEX, the index of the factor being processed. % The first factor processed must be factor 1. % % Input, integer FACTOR_ORDER, the order of the factor. % % Input, real FACTOR_VALUE(FACTOR_ORDER), the factor values for % factor FACTOR_INDEX. % % Input, integer FACTOR_NUM, the number of factors. % % Input, integer POINT_NUM, the number of elements in the direct product. % % Input/output, real W(POINT_NUM), the elements of the % direct product, updated by the latest factor. % % Local Parameters: % % Local, integer START, the first location of a block of values to set. % % Local, integer CONTIG, the number of consecutive values to set. % % Local, integer SKIP, the distance from the current value of START % to the next location of a block of values to set. % % Local, integer REP, the number of blocks of values to set. % persistent contig; persistent rep; persistent skip; if ( factor_index == 1 ) contig = 1; skip = 1; rep = point_num; w(1:point_num) = 1.0; end rep = rep / factor_order; skip = skip * factor_order; for j = 1 : factor_order start = 1 + ( j - 1 ) * contig; for k = 1 : rep w(start:start+contig-1) = w(start:start+contig-1) * factor_value(j); start = start + skip; end end contig = contig * factor_order; return end function r8vec_print_some ( n, a, i_lo, i_hi, title ) %*****************************************************************************80 % %% R8VEC_PRINT_SOME prints "some" of an R8VEC. % % Discussion: % % An R8VEC is a vector of R8 values. % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 16 October 2006 % % Author: % % John Burkardt % % Parameters: % % Input, integer N, the dimension of the vector. % % Input, real A(N), the vector to be printed. % % Input, integer MAX_PRINT, the maximum number of lines to print. % % Input, string TITLE, a title. % fprintf ( 1, '\n' ); fprintf ( 1, '%s\n', title ); fprintf ( 1, '\n' ); for i = max ( 1, i_lo ) : min ( n, i_hi ) fprintf ( 1, ' %8d %12f\n', i, a(i) ); end return end function [ grid_weight, grid_point ] = sparse_grid_herm ( dim_num, ... level_max, point_num ) %*****************************************************************************80 % %% SPARSE_GRID_HERM computes a sparse grid of Gauss-Hermite points. % % Discussion: % % The quadrature rule is associated with a sparse grid derived from % a Smolyak construction using a 1D Gauss-Hermite quadrature rule. % % The user specifies: % * the spatial dimension of the quadrature region, % * the level that defines the Smolyak grid. % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 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, integer DIM_NUM, the spatial dimension. % % Input, integer LEVEL_MAX, controls the size of the % sparse grid. % % Input, integer POINT_NUM, the number of points in the grid, % as determined by SPARSE_GRID_HERM_SIZE. % % Output, real GRID_WEIGHT(1,POINT_NUM), the weights. % % Output, real GRID_POINT(DIM_NUM,POINT_NUM), the points. % grid_point = zeros ( dim_num, point_num ); grid_weight = zeros ( 1, point_num ); % % The outer loop generates LEVELs from LEVEL_MIN to LEVEL_MAX. % point_num2 = 0; level_min = max ( 0, level_max + 1 - dim_num ); for level = level_min : level_max % % The middle loop generates the next partition LEVEL_1D(1:DIM_NUM) % that adds up to LEVEL. % level_1d = []; more = 0; h = 0; t = 0; while ( 1 ) [ level_1d, more, h, t ] = comp_next ( level, dim_num, level_1d, more, h, t ); % % Transform each 1D level to a corresponding 1D order. % The relationship is the same as for other OPEN rules. % The Gauss-Hermite rule differs from the other OPEN rules % only in the nesting behavior. % order_1d = level_to_order_open ( dim_num, level_1d ); grid_base2(1:dim_num) = round ( ( order_1d(1:dim_num) - 1 ) / 2 ); % % The product of the 1D orders gives us the number of points in this subgrid. % order_nd = prod ( order_1d(1:dim_num) ); % % Compute the weights for this product grid. % grid_weight2 = product_weight_herm ( dim_num, order_1d, order_nd ); % % Now determine the coefficient of the weight. % coeff = (-1)^( level_max - level ) ... * i4_choose ( dim_num - 1, level_max - level ); % % The inner (hidden) loop generates all points corresponding to given grid. % The grid indices will be between -M to +M, where 2*M + 1 = ORDER_1D(DIM). % grid_index2 = multigrid_index_z ( dim_num, order_1d, order_nd ); % % Determine the first level of appearance of each of the points. % This allows us to flag certain points as being repeats of points % generated on a grid of lower level. % % This is SLIGHTLY tricky. % grid_level = index_level_herm ( level, level_max, dim_num, order_nd, ... grid_index2, grid_base2 ); % % Only keep those points which first appear on this level. % for point = 1 : order_nd % % Either a "new" point (increase count, create point, create weight) % if ( grid_level(point) == level ) point_num2 = point_num2 + 1; grid_point(1:dim_num,point_num2) = hermite_abscissa ( dim_num, 1, ... grid_index2(1:dim_num,point), grid_base2(1:dim_num) ); grid_weight(point_num2) = coeff * grid_weight2(point); % % or an already existing point (create point temporarily, find match, % add weight to matched point's weight). % else grid_point_temp(1:dim_num) = hermite_abscissa ( dim_num, 1, ... grid_index2(1:dim_num,point), grid_base2(1:dim_num) ); point3 = -1; for point2 = 1 : point_num2 if ( grid_point(1:dim_num,point2) == grid_point_temp(1:dim_num)' ) point3 = point2; break end end if ( point3 == -1 ) fprintf ( 1, '\n' ); fprintf ( 1, 'SPARSE_GRID_HERMITE - Fatal error!\n' ); fprintf ( 1, ' Could not match point.\n' ); error ( 'SPARSE_GRID_HERMITE - Fatal error!' ); end grid_weight(point3) = grid_weight(point3) ... + coeff * grid_weight2(point); end end if ( ~more ) break end end end return end function point_num = sparse_grid_herm_size ( dim_num, level_max ) %*****************************************************************************80 % %% SPARSE_GRID_HERM_SIZE sizes a sparse grid of Gauss-Hermite points. % % Discussion: % % The grid is defined as the sum of the product rules whose LEVEL % satisfies: % % LEVEL_MIN <= LEVEL <= LEVEL_MAX. % % where LEVEL_MAX is user specified, and % % LEVEL_MIN = max ( 0, LEVEL_MAX + 1 - DIM_NUM ). % % The grids are only very weakly nested, since Gauss-Hermite rules % only have the origin in common. % % 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, integer DIM_NUM, the spatial dimension. % % Input, integer LEVEL_MAX, the maximum value of LEVEL. % % Output, integer POINT_NUM, the number of points in the grid. % % % Special case. % if ( level_max == 0 ) point_num = 1; return; end % % The outer loop generates LEVELs from LEVEL_MIN to LEVEL_MAX. % point_num = 0; level_min = max ( 0, level_max + 1 - dim_num ); for level = level_min : level_max % % The middle loop generates the next partition that adds up to LEVEL. % level_1d = []; more = false; h = 0; t = 0; while ( 1 ) [ level_1d, more, h, t ] = comp_next ( level, dim_num, level_1d, more, h, t ); % % Transform each 1D level to a corresponding 1D order. % order_1d = level_to_order_open ( dim_num, level_1d ); for dim = 1 : dim_num % % If we can reduce the level in this dimension by 1 and % still not go below LEVEL_MIN. % if ( level_min < level && 1 < order_1d(dim) ) order_1d(dim) = order_1d(dim) - 1; end end point_num = point_num + prod ( order_1d(1:dim_num) ); if ( ~more ) break end end end return end function timestamp ( ) %*****************************************************************************80 % %% TIMESTAMP prints the current YMDHMS date as a timestamp. % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 14 February 2003 % % Author: % % John Burkardt % t = now; c = datevec ( t ); s = datestr ( c, 0 ); fprintf ( 1, '%s\n', s ); return end function [ a, more ] = vec_colex_next2 ( dim_num, base, a, more ) %*****************************************************************************80 % %% 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,BASE(2)-1,...,BASE(DIM_NUM)-1). % % Example: % % 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: % % 11 November 2007 % % Author: % % John Burkardt % % Reference: % % Dennis Stanton, Dennis White, % Constructive Combinatorics, % Springer, 1986, % ISBN: 0387963472, % LC: QA164.S79. % % Parameters: % % Input, integer DIM_NUM, the spatial dimension. % % Input, integer BASE(DIM_NUM), the base to be used in each dimension. % % Input, integer A(DIM_NUM), except on the first call, this should % be the output value of A on the last call. % % Input, logical MORE, should be FALSE on the first call, and % thereafter should be the output value of MORE from the previous call. % % Output, integer A(DIM_NUM), the next vector. % % Output, logical MORE, is TRUE if another vector was computed. % If MORE is FALSE on return, then ignore the output value A, and % stop calling the routine. % if ( ~more ) a(1:dim_num) = 0; more = 1; else for i = 1 : dim_num a(i) = a(i) + 1; if ( a(i) < base(i) ) return end a(i) = 0; end more = 0; end return end