File: example_sparse_grids_02.cpp

package info (click to toggle)
tasmanian 8.2-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 4,852 kB
  • sloc: cpp: 34,523; python: 7,039; f90: 5,080; makefile: 224; sh: 64; ansic: 8
file content (93 lines) | stat: -rw-r--r-- 3,158 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93

#include "Tasmanian.hpp"

using namespace std;

/*!
 * \internal
 * \file example_sparse_grids_02.cpp
 * \brief Examples for the Tasmanian Sparse Grid module.
 * \author Miroslav Stoyanov
 * \ingroup TasmanianSGExamples
 *
 * Tasmanian Sparse Grids Example 2.
 * \endinternal
 */

/*!
 * \ingroup TasmanianSGExamples
 * \addtogroup TasmanianSGExamples2 Tasmanian Sparse Grids module, example 2
 *
 * \par Example 2
 * Compute \f$ \int_{-5}^{+5} \int_{-2}^{+3} \exp(-x^2) \cos(y) dy dx \f$
 * 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.
 */

/*!
 * \ingroup TasmanianSGExamples2
 * \brief Sparse Grids Example 2: integrate a simple function over a canonical domain.
 *
 * \snippet SparseGrids/Examples/example_sparse_grids_02.cpp SG_Example_02 example
 */
void sparse_grids_example_02(){
#ifndef __TASMANIAN_DOXYGEN_SKIP
//! [SG_Example_02 example]
#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,
                                        TasGrid::type_qptotal, TasGrid::rule_gausspatterson);
    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,
                                   TasGrid::type_qptotal, TasGrid::rule_gausspatterson);
    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
//! [SG_Example_02 example]
#endif
}