4 #ifndef DUNE_LOCALFUNCTIONS_QK_LOCALCOEFFICIENTS_HH
5 #define DUNE_LOCALFUNCTIONS_QK_LOCALCOEFFICIENTS_HH
10 #include <dune/common/exceptions.hh>
11 #include <dune/common/power.hh>
22 template<
int k,
int d>
26 static array<unsigned int,d> multiindex (
unsigned int i)
28 array<unsigned int,d> alpha;
29 for (
int j=0; j<d; j++)
38 void setup2d(std::vector<unsigned int>& subEntity)
55 subEntity[lastIndex++] = 0;
56 for (
unsigned i = 0; i < k - 1; ++i)
57 subEntity[lastIndex++] = 2;
59 subEntity[lastIndex++] = 1;
62 for (
unsigned e = 0; e < k - 1; ++e) {
63 subEntity[lastIndex++] = 0;
64 for (
unsigned i = 0; i < k - 1; ++i)
65 subEntity[lastIndex++] = 0;
66 subEntity[lastIndex++] = 1;
70 subEntity[lastIndex++] = 2;
71 for (
unsigned i = 0; i < k - 1; ++i)
72 subEntity[lastIndex++] = 3;
74 subEntity[lastIndex++] = 3;
76 assert((StaticPower<k+1,d>::power==lastIndex));
80 void setup3d(std::vector<unsigned int>& subEntity)
85 const unsigned numIndices = StaticPower<k+1,d>::power;
86 const unsigned numFaceIndices = StaticPower<k+1,d-1>::power;
88 const unsigned numInnerEdgeDofs = k-1;
104 subEntity[lastIndex++] = 0;
105 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
106 subEntity[lastIndex++] = 6;
108 subEntity[lastIndex++] = 1;
111 for (
unsigned e = 0; e < numInnerEdgeDofs; ++e) {
112 subEntity[lastIndex++] = 4;
113 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
114 subEntity[lastIndex++] = 4;
115 subEntity[lastIndex++] = 5;
119 subEntity[lastIndex++] = 2;
120 for (
unsigned i = 0; i < k - 1; ++i)
121 subEntity[lastIndex++] = 7;
122 subEntity[lastIndex++] = 3;
124 assert(numFaceIndices==lastIndex);
128 for(
unsigned f = 0; f < numInnerEdgeDofs; ++f) {
131 subEntity[lastIndex++] = 0;
132 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
133 subEntity[lastIndex++] = 2;
134 subEntity[lastIndex++] = 1;
137 for (
unsigned e = 0; e < numInnerEdgeDofs; ++e) {
138 subEntity[lastIndex++] = 0;
139 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
140 subEntity[lastIndex++] = 0;
141 subEntity[lastIndex++] = 1;
145 subEntity[lastIndex++] = 2;
146 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
147 subEntity[lastIndex++] = 3;
148 subEntity[lastIndex++] = 3;
150 assert(lastIndex==(f+1+1)*numFaceIndices);
155 subEntity[lastIndex++] = 4;
156 for (
unsigned i = 0; i < k - 1; ++i)
157 subEntity[lastIndex++] = 10;
158 subEntity[lastIndex++] = 5;
161 for (
unsigned e = 0; e < k - 1; ++e) {
162 subEntity[lastIndex++] = 8;
163 for (
unsigned i = 0; i < k - 1; ++i)
164 subEntity[lastIndex++] = 5;
165 subEntity[lastIndex++] = 9;
169 subEntity[lastIndex++] = 6;
170 for (
unsigned i = 0; i < k - 1; ++i)
171 subEntity[lastIndex++] = 11;
172 subEntity[lastIndex++] = 7;
174 assert(numIndices==lastIndex);
182 std::vector<unsigned int> codim(li.size());
184 for (std::size_t i=0; i<codim.size(); i++) {
188 array<unsigned int,d> mIdx = multiindex(i);
189 for (
int j=0; j<d; j++)
190 if (mIdx[j]==0 or mIdx[j]==k)
199 std::vector<unsigned int> index(
size());
201 for (std::size_t i=0; i<
size(); i++) {
205 array<unsigned int,d> mIdx = multiindex(i);
207 for (
int j=d-1; j>=0; j--)
208 if (mIdx[j]>0 and mIdx[j]<k)
209 index[i] = (k-1)*index[i] + (mIdx[j]-1);
214 std::vector<unsigned int> subEntity(li.size());
218 for (std::size_t i=0; i<
size(); i++)
230 DUNE_THROW(Dune::NotImplemented,
"QkLocalCoefficients for k==" << k <<
" and d==" << d);
232 for (
size_t i=0; i<li.size(); i++)
233 li[i] =
LocalKey(subEntity[i], codim[i], index[i]);
239 return StaticPower<k+1,d>::power;
249 std::vector<LocalKey> li;