4#ifndef DUNE_LOCALFUNCTIONS_SERENDIPITY_LOCALCOEFFICIENTS_HH
5#define DUNE_LOCALFUNCTIONS_SERENDIPITY_LOCALCOEFFICIENTS_HH
11#include <dune/common/exceptions.hh>
13#include <dune/localfunctions/common/localkey.hh>
22 template<
int k,
int d>
25 typedef typename Dune::FieldVector<int,d> DuneVector;
27 static void multiindex (DuneVector& alpha)
30 for (
int j=0; j<d; j++){
31 i+=alpha[j]*pow(k+1,j);
34 for (
int j=0; j<d; j++)
40 for (
unsigned int j=0; j<d; j++){
41 if(alpha[j] > 0 && alpha[j] < k){ sum++; }
49 void setup1d(std::vector<unsigned int>& subEntity)
65 subEntity[lastIndex++] = 0;
66 for (
unsigned i = 0; i < k - 1; ++i)
67 subEntity[lastIndex++] = 0;
69 subEntity[lastIndex++] = 1;
71 assert((
size()==lastIndex));
74 void setup2d(std::vector<unsigned int>& subEntity)
97 subEntity[lastIndex++] = 0;
98 for (
unsigned i = 0; i < k - 1; ++i)
99 subEntity[lastIndex++] = 2;
100 subEntity[lastIndex++] = 1;
103 for (
unsigned e = 0; e < k - 1; ++e) {
104 subEntity[lastIndex++] = 0;
105 subEntity[lastIndex++] = 1;
109 subEntity[lastIndex++] = 2;
110 for (
unsigned i = 0; i < k - 1; ++i)
111 subEntity[lastIndex++] = 3;
112 subEntity[lastIndex++] = 3;
114 assert((
size()==lastIndex));
119 void setup3d(std::vector<unsigned int>& subEntity)
129 const unsigned numIndices =
size();
131 const unsigned numInnerEdgeDofs = k-1;
150 unsigned lastIndex=0;
152 subEntity[lastIndex++] = 0;
153 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
154 subEntity[lastIndex++] = 6;
155 subEntity[lastIndex++] = 1;
158 for (
unsigned e = 0; e < numInnerEdgeDofs; ++e) {
159 subEntity[lastIndex++] = 4;
160 subEntity[lastIndex++] = 5;
164 subEntity[lastIndex++] = 2;
165 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
166 subEntity[lastIndex++] = 7;
167 subEntity[lastIndex++] = 3;
172 for(
unsigned f = 0; f < numInnerEdgeDofs; ++f) {
174 subEntity[lastIndex++] = 0;
175 subEntity[lastIndex++] = 1;
178 subEntity[lastIndex++] = 2;
179 subEntity[lastIndex++] = 3;
185 subEntity[lastIndex++] = 4;
186 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
187 subEntity[lastIndex++] = 10;
188 subEntity[lastIndex++] = 5;
191 for (
unsigned e = 0; e < numInnerEdgeDofs; ++e) {
192 subEntity[lastIndex++] = 8;
193 subEntity[lastIndex++] = 9;
197 subEntity[lastIndex++] = 6;
198 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
199 subEntity[lastIndex++] = 11;
200 subEntity[lastIndex++] = 7;
203 assert(numIndices==lastIndex);
211 std::vector<unsigned int> codim(li.size());
213 for (std::size_t i=0; i <
size(); i++) {
219 if(i>0){ multiindex(mIdx); }
220 for (
int j=0; j<d; j++)
221 if (mIdx[j]==0 or mIdx[j]==k)
230 std::vector<unsigned int> index(
size());
232 for (std::size_t i=0; i<
size(); i++) {
234 if(i>0){ multiindex(mIdx); }
236 for (
int j=d-1; j>=0; j--)
237 if (mIdx[j]>0 and mIdx[j]<k)
238 index[i] = (k-1)*index[i] + (mIdx[j]-1);
242 std::vector<unsigned int> subEntity(li.size());
245 for (std::size_t i=0; i<
size(); i++)
258 DUNE_THROW(Dune::NotImplemented,
"SerendipityLocalCoefficients for k==" << k <<
" and d==" << d);
260 for (
size_t i=0; i<li.size(); i++){
261 li[i] = LocalKey(subEntity[i], codim[i], index[i]);
272 return 2*(k+1)+2*(k-1);
275 return 4*(k+1)+8*(k-1);
286 std::vector<LocalKey> li;
Attaches a shape function to an entity.
Definition: serendipitylocalcoefficients.hh:23
SerendipityLocalCoefficients()
Default constructor.
Definition: serendipitylocalcoefficients.hh:208
const LocalKey & localKey(std::size_t i) const
get i'th index
Definition: serendipitylocalcoefficients.hh:280
std::size_t size() const
number of coefficients
Definition: serendipitylocalcoefficients.hh:266