DUNE-FEM (unstable)

pyramidpoints.hh
1#ifndef DUNE_FEM_PYRAMIDPOINTS_HH
2#define DUNE_FEM_PYRAMIDPOINTS_HH
3
6
7namespace Dune
8{
9
10 namespace Fem
11 {
12
17 {
18 public:
19 enum { numQuads = 2 };
20 enum { MAXP = 8 };
21 enum { highest_order = 2 };
22
24 const FieldVector<double, 3>& point(int m, int i) const
25 {
26 return G[m][i];
27 }
28
30 double weight (int m, int i) const
31 {
32 return W[m][i];
33 }
34
36 int order (int m) const
37 {
38 return O[m];
39 }
40
42 int numPoints(int m) const
43 {
44 return N[m];
45 }
46
49 {
50 int m = 0;
51 O[m] = 0;
52 N[m] = 0;
53
54 // polynom degree 2 ???
55 m = 1;
56 G[m][0][0] =0.58541020;
57 G[m][0][1] =0.72819660;
58 G[m][0][2] =0.13819660;
59
60 G[m][1][0] =0.13819660;
61 G[m][1][1] =0.72819660;
62 G[m][1][2] =0.13819660;
63
64 G[m][2][0] =0.13819660;
65 G[m][2][1] =0.27630920;
66 G[m][2][2] =0.58541020;
67
68 G[m][3][0] =0.13819660;
69 G[m][3][1] =0.27630920;
70 G[m][3][2] =0.13819660;
71
72 G[m][4][0] =0.72819660;
73 G[m][4][1] =0.13819660;
74 G[m][4][2] =0.13819660;
75
76 G[m][5][0] =0.72819660;
77 G[m][5][1] =0.58541020;
78 G[m][5][2] =0.13819660;
79
80 G[m][6][0] =0.27630920;
81 G[m][6][1] =0.13819660;
82 G[m][6][2] =0.58541020;
83
84 G[m][7][0] =0.27630920;
85 G[m][7][1] =0.13819660;
86 G[m][7][2] =0.13819660;
87
88 W[m][0] = 0.041666667;
89 W[m][1] = 0.041666667;
90 W[m][2] = 0.041666667;
91 W[m][3] = 0.041666667;
92 W[m][4] = 0.041666667;
93 W[m][5] = 0.041666667;
94 W[m][6] = 0.041666667;
95 W[m][7] = 0.041666667;
96
97 O[m] = 2;// verify ????????
98 N[m] = 8;
99 }
100
101 private:
102 FieldVector<double, 3> G[numQuads][MAXP];
103 double W[numQuads][MAXP]; // weights associated with points
104 int O[numQuads]; // order of the rule
105 int N[numQuads]; // number of points in quadrature rule
106 };
107
108 } // namespace Fem
109
110} // namespace Dune
111
112#endif // #ifndef DUNE_FEM_PYRAMIDPOINTS_HH
Definition: pyramidpoints.hh:17
const FieldVector< double, 3 > & point(int m, int i) const
Access to the ith point of quadrature rule m.
Definition: pyramidpoints.hh:24
PyramidPoints()
initialize quadrature points on the interval for all orders
Definition: pyramidpoints.hh:48
int numPoints(int m) const
Number of points in the quadrature rule m.
Definition: pyramidpoints.hh:42
int order(int m) const
Actual order of quadrature rule m.
Definition: pyramidpoints.hh:36
double weight(int m, int i) const
Access to the ith weight of quadrature rule m.
Definition: pyramidpoints.hh:30
Implements a vector constructed from a given type representing a field and a compile-time given size.
Dune namespace.
Definition: alignedallocator.hh:13
Definition of macros controlling symbol visibility at the ABI level.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 13, 23:29, 2024)