Newer
Older
Spencer Sherwin
committed
#include <cstdio>
#include <cstdlib>
#include <LibUtilities/Polylib/Polylib.h>
#include <LibUtilities/Foundations/ManagerAccess.h>
using namespace Nektar;
using namespace std;
#define WITHSOLUTION 1
int main(int, char **)
Spencer Sherwin
committed
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
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
{
cout << "==========================================================="<< endl;
cout << "| INTEGRATION ON 2D ELEMENT in Local Region |"<< endl;
cout << "==========================================================="<< endl;
cout << endl;
cout << "Integrate the function f(x1,x2) = x1^12 * x2^14 " << endl;
cout << "on a local quadrilateral element:" << endl;
// Specify the number of quadrature points in both directions
int nQuadPointsDir1 = 6;
int nQuadPointsDir2 = 7;
// Specify the type of quadrature points in both directions
LibUtilities::PointsType quadPointsTypeDir1 =
LibUtilities::eGaussLobattoLegendre;
LibUtilities::PointsType quadPointsTypeDir2 =
LibUtilities::eGaussLobattoLegendre;
// Declare variables (of type Array) to hold the quadrature
// zeros and weights
Array<OneD, NekDouble> quadZerosDir1 (nQuadPointsDir1);
Array<OneD, NekDouble> quadWeightsDir1(nQuadPointsDir1);
Array<OneD, NekDouble> quadZerosDir2 (nQuadPointsDir2);
Array<OneD, NekDouble> quadWeightsDir2(nQuadPointsDir2);
// Calculate the GLL-quadrature zeros and weights in both
// directions. This is done in 2 steps.
//
// Step 1: Declare the PointsKeys which uniquely defines the
// quadrature points
const LibUtilities::PointsKey quadPointsKeyDir1(nQuadPointsDir1,
quadPointsTypeDir1);
const LibUtilities::PointsKey quadPointsKeyDir2(nQuadPointsDir2,
quadPointsTypeDir2);
// Step 2: Using this key, the quadrature zeros and weights
// can now be retrieved through the PointsManager
quadZerosDir1 = LibUtilities::PointsManager()[quadPointsKeyDir1]->GetZ();
quadWeightsDir1 = LibUtilities::PointsManager()[quadPointsKeyDir1]->GetW();
quadZerosDir2 = LibUtilities::PointsManager()[quadPointsKeyDir2]->GetZ();
quadWeightsDir2 = LibUtilities::PointsManager()[quadPointsKeyDir2]->GetW();
// The local (straight-sided) quadrilateral element has the
// following vertices:
//
// - Vertex A: (x1_A,x2_A) = (0,-1)
// - Vertex B: (x1_A,x2_A) = (1,-1)
// - Vertex C: (x1_A,x2_A) = (1,1)
// - Vertex D: (x1_A,x2_A) = (0,0)
//
NekDouble x1_A = 0.0;
NekDouble x2_A = -1.0;
NekDouble x1_B = 1.0;
NekDouble x2_B = -1.0;
NekDouble x1_C = 1.0;
NekDouble x2_C = 1.0;
NekDouble x1_D = 0.0;
NekDouble x2_D = 0.0;
// Integrate the function f(x1,x2) = x1^12 * x2^14 on a local
// quadrilateral element (defined above). Use 7th order
// Gauss-Lobatto-Legendre quadrature in both direction.
//
// Your code can be based on the previous exercise. However,
// as we are calculating the integral of a function defined on
// a local element rather than on a reference element, we have
// to take into account the geometry of the element.
//
// Therefore, the implementation should be altered in two ways:
//
// (1) The quadrature zeros should be transformed to local
// coordinates to evaluate the integrand f(x1,x2)
//
// (2) Take into account the Jacobian of the transformation
// between local and reference coordinates when evaluating
// the integral. (Evaluate the expression for the Jacobian
// analytically rather than using numerical
// differentiation)
//
// Store the solution in the variable 'result'
NekDouble result = 0.0;
// Apply the Gaussian quadrature technique to integrate the
// function f(x1,x2) = x1^12 * x2^14 on the standard
// quadrilateral. To do so, edit the (double) loop which performs
// the summation.
double x1;
double x2;
double jacobian;
for(int i = 0; i < nQuadPointsDir1; i++)
{
for(int j = 0; j < nQuadPointsDir2; j++)
{
// Calculate the local coordinates of the quadrature zeros
// using the mapping from reference to local element
x1 = x1_A * 0.25 * (1-quadZerosDir1[i]) * (1-quadZerosDir2[j]) +
x1_B * 0.25 * (1+quadZerosDir1[i]) * (1-quadZerosDir2[j]) +
x1_D * 0.25 * (1-quadZerosDir1[i]) * (1+quadZerosDir2[j]) +
x1_C * 0.25 * (1+quadZerosDir1[i]) * (1+quadZerosDir2[j]);
x2 = x2_A * 0.25 * (1-quadZerosDir1[i]) * (1-quadZerosDir2[j]) +
x2_B * 0.25 * (1+quadZerosDir1[i]) * (1-quadZerosDir2[j]) +
x2_D * 0.25 * (1-quadZerosDir1[i]) * (1+quadZerosDir2[j]) +
x2_C * 0.25 * (1+quadZerosDir1[i]) * (1+quadZerosDir2[j]);
#if WITHSOLUTION
// Analytically evaluate the Jacobian
jacobian = (0.25 * (1-quadZerosDir2[j]) * (x1_B-x1_A) +
0.25 * (1+quadZerosDir2[j]) * (x1_C-x1_D)) *
(0.25 * (1-quadZerosDir1[i]) * (x2_D-x2_A) +
0.25 * (1+quadZerosDir1[i]) * (x2_C-x2_B)) -
(0.25 * (1-quadZerosDir1[i]) * (x1_D-x1_A) +
0.25 * (1+quadZerosDir1[i]) * (x1_C-x1_B)) *
(0.25 * (1-quadZerosDir2[j]) * (x2_B-x2_A) +
0.25 * (1+quadZerosDir2[j]) * (x2_C-x2_D));
#endif
result += pow(x1,12) * pow(x2,14) *
quadWeightsDir1[i] * quadWeightsDir2[j] * fabs(jacobian);
}
}
// Display the output
NekDouble exactResult = 1.0/195.0+1.0/420.0;
cout << "\t Error = "<< fabs(result-exactResult) << endl;
}