createProductRuleGrid {SparseGrid}R Documentation

Create a multidimensional grid of nodes and weights for integration

Description

Creates nodes and weights according to the product rule, combining 1D nodes and weights. Sparse grids can be created with the function createSparseGrid.

Usage

createProductRuleGrid(type, dimension, k, sym = FALSE)

Arguments

type

String or function for type of 1D integration rule, can take on values

"KPU"

Nested rule for unweighted integral over [0,1]

"KPN"

Nested rule for integral with Gaussian weight

"GQU"

Gaussian quadrature for unweighted integral over [0,1] (Gauss-Legendre)

"GQN"

Gaussian quadrature for integral with Gaussian weight (Gauss-Hermite)

func

any function. Function must accept level k and return a list with two elements nodes and weights for univariate quadrature rule with polynomial exactness 2k-1.

dimension

dimension of the integration problem.

k

Accuracy level. The rule will be exact for polynomial up to total order 2k-1.

sym

(optional) only used for own 1D quadrature rule (type not "KPU",...). If sym is supplied and not FALSE, the code will run faster but will produce incorrect results if 1D quadrature rule is asymmetric.

Value

The return value contains a list with nodes and weights

nodes

matrix with a node in each row

weights

vector with corresponding weights

Author(s)

Jelmer Ypma

References

Florian Heiss, Viktor Winschel, Likelihood approximation by numerical integration on sparse grids, Journal of Econometrics, Volume 144, Issue 1, May 2008, Pages 62-80, http://www.sparse-grids.de

See Also

createSparseGrid createMonteCarloGrid createIntegrationGrid integrate pmvnorm

Examples

# load library
library('SparseGrid')

# define function to be integrated
# g(x) = x[1] * x[2] * ... * x[n]
g <- function( x ) {
	return( prod( x ) )
}

#
# Create sparse integration grid to approximate integral of a function with uniform weights
#
sp.grid <- createSparseGrid( 'KPU', dimension=3, k=5 )

# number of nodes and weights
length( sp.grid$weights )

# evaluate function g in nodes
gx.sp <- apply( sp.grid$nodes, 1, g )

# take weighted sum to get approximation for the integral
val.sp <- gx.sp %*% sp.grid$weights

#
# Create integration grid to approximate integral of a function with uniform weights
#
pr.grid <- createProductRuleGrid( 'KPU', dimension=3, k=5 )

# number of nodes and weights
length( pr.grid$weights )

# evaluate function g in nodes
gx.pr <- apply( pr.grid$nodes, 1, g )

# take weighted sum to get approximation for the integral
val.pr <- gx.pr %*% pr.grid$weights

#
# Create integration grid to approximation integral using Monte Carlo simulation
#
set.seed( 3141 )
mc.grid <- createMonteCarloGrid( runif, dimension=3, num.sim=1000 )

# number of nodes and weights
length( mc.grid$weights )

# evaluate function g in MC nodes
gx.mc   <- apply( mc.grid$nodes, 1, g )

# take weighted sum to get approximation for the integral
# the weights are all equal to 1/1000 in this case
val.mc <- gx.mc %*% mc.grid$weights

val.sp
val.pr
val.mc

[Package SparseGrid version 0.8.2 Index]