Doxygen 1.9.1
Toolkit for Adaptive Stochastic Modeling and Non-Intrusive ApproximatioN: Tasmanian v8.1
Tasmanian Sparse Grids module, example 2
Collaboration diagram for Tasmanian Sparse Grids module, example 2:

Functions

void sparse_grids_example_02 ()
 Sparse Grids Example 2: integrate a simple function over a canonical domain. More...
 

Detailed Description

Example 2
Compute $ \int_{-5}^{+5} \int_{-2}^{+3} \exp(-x^2) \cos(y) dy dx $ using a sparse grid with Gauss-Patterson nodes and weights. The grid is constructed over a non-canononical domain and the points are selected based on a total degree polynomial space.

Function Documentation

◆ sparse_grids_example_02()

void sparse_grids_example_02 ( )

Sparse Grids Example 2: integrate a simple function over a canonical domain.

#endif
cout << "\n---------------------------------------------------------------------------------------------------\n";
cout << std::scientific; cout.precision(17);
cout << "Example 2: integrate f(x,y) = exp(-x^2) * cos(y) over [-5,5] x [-2,3]\n"
<< " using Gauss-Patterson nodes and total degree polynomial space\n";
int dimension = 2;
int exactness = 20;
// the type_qptotal will guarantee exact integral for all polynomials with degree 20 or less
auto grid = TasGrid::makeGlobalGrid(dimension, 0, exactness,
grid.setDomainTransform({-5.0, -2.0}, {5.0, 3.0}); // set the non-canonical domain
auto points = grid.getPoints();
auto weights = grid.getQuadratureWeights();
double I = 0.0;
for(int i=0; i<grid.getNumPoints(); i++){
double x = points[i*dimension];
double y = points[i*dimension+1];
I += weights[i] * std::exp(-x*x) * std::cos(y);
}
double exact = 1.861816427518323e+00;
double E = std::abs(exact - I);
cout << " at polynomial exactness: " << exactness
<< "\n the grid has: " << weights.size()
<< "\n integral: " << I
<< "\n error: " << E << "\n\n";
exactness = 40;
grid = TasGrid::makeGlobalGrid(dimension, 0, exactness,
grid.setDomainTransform({-5.0, -2.0}, {5.0, 3.0}); // set the non-canonical domain
points = grid.getPoints();
weights = grid.getQuadratureWeights();
I = 0.0;
for(int i=0; i<grid.getNumPoints(); i++){
double x = points[i*dimension];
double y = points[i*dimension+1];
I += weights[i] * std::exp(-x*x) * std::cos(y);
}
E = std::abs(exact - I);
cout << " at polynomial exactness: " << exactness
<< "\n the grid has: " << weights.size()
<< "\n integral: " << I
<< "\n error: " << E << endl;
#ifndef __TASMANIAN_DOXYGEN_SKIP
@ rule_gausspatterson
Nested rule that is optimized for integration, probably the best integration rule in more than 2 dime...
Definition: tsgEnumerates.hpp:333
@ type_qptotal
Total degree polynomial space for quadrature/integration.
Definition: tsgEnumerates.hpp:224
TasmanianSparseGrid makeGlobalGrid(int dimensions, int outputs, int depth, TypeDepth type, TypeOneDRule rule, std::vector< int > const &anisotropic_weights=std::vector< int >(), double alpha=0.0, double beta=0.0, const char *custom_filename=nullptr, std::vector< int > const &level_limits=std::vector< int >())
Factory method, creates a new grid and calls TasmanianSparseGrid::makeGlobalGrid().
Definition: TasmanianSparseGrid.hpp:2272