4#ifndef DUNE_LOCALFUNCTIONS_SERENDIPITY_LOCALCOEFFICIENTS_HH
5#define DUNE_LOCALFUNCTIONS_SERENDIPITY_LOCALCOEFFICIENTS_HH
11#include <dune/common/exceptions.hh>
12#include <dune/common/power.hh>
14#include <dune/localfunctions/common/localkey.hh>
23 template<
int k,
int d>
26 typedef typename Dune::FieldVector<int,d> DuneVector;
28 static void multiindex (DuneVector& alpha)
31 for (
int j=0; j<d; j++){
32 i+=alpha[j]*pow(k+1,j);
35 for (
int j=0; j<d; j++)
41 for (
unsigned int j=0; j<d; j++){
42 if(alpha[j] > 0 && alpha[j] < k){ sum++; }
50 void setup1d(std::vector<unsigned int>& subEntity)
66 subEntity[lastIndex++] = 0;
67 for (
unsigned i = 0; i < k - 1; ++i)
68 subEntity[lastIndex++] = 0;
70 subEntity[lastIndex++] = 1;
72 assert((
size()==lastIndex));
75 void setup2d(std::vector<unsigned int>& subEntity)
98 subEntity[lastIndex++] = 0;
99 for (
unsigned i = 0; i < k - 1; ++i)
100 subEntity[lastIndex++] = 2;
101 subEntity[lastIndex++] = 1;
104 for (
unsigned e = 0; e < k - 1; ++e) {
105 subEntity[lastIndex++] = 0;
106 subEntity[lastIndex++] = 1;
110 subEntity[lastIndex++] = 2;
111 for (
unsigned i = 0; i < k - 1; ++i)
112 subEntity[lastIndex++] = 3;
113 subEntity[lastIndex++] = 3;
115 assert((
size()==lastIndex));
120 void setup3d(std::vector<unsigned int>& subEntity)
130 const unsigned numIndices =
size();
132 const unsigned numInnerEdgeDofs = k-1;
151 unsigned lastIndex=0;
153 subEntity[lastIndex++] = 0;
154 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
155 subEntity[lastIndex++] = 6;
156 subEntity[lastIndex++] = 1;
159 for (
unsigned e = 0; e < numInnerEdgeDofs; ++e) {
160 subEntity[lastIndex++] = 4;
161 subEntity[lastIndex++] = 5;
165 subEntity[lastIndex++] = 2;
166 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
167 subEntity[lastIndex++] = 7;
168 subEntity[lastIndex++] = 3;
173 for(
unsigned f = 0; f < numInnerEdgeDofs; ++f) {
175 subEntity[lastIndex++] = 0;
176 subEntity[lastIndex++] = 1;
179 subEntity[lastIndex++] = 2;
180 subEntity[lastIndex++] = 3;
186 subEntity[lastIndex++] = 4;
187 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
188 subEntity[lastIndex++] = 10;
189 subEntity[lastIndex++] = 5;
192 for (
unsigned e = 0; e < numInnerEdgeDofs; ++e) {
193 subEntity[lastIndex++] = 8;
194 subEntity[lastIndex++] = 9;
198 subEntity[lastIndex++] = 6;
199 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
200 subEntity[lastIndex++] = 11;
201 subEntity[lastIndex++] = 7;
204 assert(numIndices==lastIndex);
212 std::vector<unsigned int> codim(li.size());
214 for (std::size_t i=0; i <
size(); i++) {
220 if(i>0){ multiindex(mIdx); }
221 for (
int j=0; j<d; j++)
222 if (mIdx[j]==0 or mIdx[j]==k)
231 std::vector<unsigned int> index(
size());
233 for (std::size_t i=0; i<
size(); i++) {
235 if(i>0){ multiindex(mIdx); }
237 for (
int j=d-1; j>=0; j--)
238 if (mIdx[j]>0 and mIdx[j]<k)
239 index[i] = (k-1)*index[i] + (mIdx[j]-1);
243 std::vector<unsigned int> subEntity(li.size());
246 for (std::size_t i=0; i<
size(); i++)
259 DUNE_THROW(Dune::NotImplemented,
"SerendipityLocalCoefficients for k==" << k <<
" and d==" << d);
261 for (
size_t i=0; i<li.size(); i++){
262 li[i] = LocalKey(subEntity[i], codim[i], index[i]);
273 return 2*(k+1)+2*(k-1);
276 return 4*(k+1)+8*(k-1);
287 std::vector<LocalKey> li;
Attaches a shape function to an entity.
Definition: serendipitylocalcoefficients.hh:24
SerendipityLocalCoefficients()
Default constructor.
Definition: serendipitylocalcoefficients.hh:209
const LocalKey & localKey(std::size_t i) const
get i'th index
Definition: serendipitylocalcoefficients.hh:281
std::size_t size() const
number of coefficients
Definition: serendipitylocalcoefficients.hh:267