dune-localfunctions  2.3beta2
qklocalcoefficients.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 
4 #ifndef DUNE_LOCALFUNCTIONS_QK_LOCALCOEFFICIENTS_HH
5 #define DUNE_LOCALFUNCTIONS_QK_LOCALCOEFFICIENTS_HH
6 
7 #include <cassert>
8 #include <vector>
9 
10 #include <dune/common/exceptions.hh>
11 #include <dune/common/power.hh>
12 
14 
15 namespace Dune
16 {
22  template<int k, int d>
24 
25  // Return i as a d-digit number in the (k+1)-nary system
26  static array<unsigned int,d> multiindex (unsigned int i)
27  {
28  array<unsigned int,d> alpha;
29  for (int j=0; j<d; j++)
30  {
31  alpha[j] = i % (k+1);
32  i = i/(k+1);
33  }
34  return alpha;
35  }
36 
37 
38  void setup2d(std::vector<unsigned int>& subEntity)
39  {
40  assert(k>0);
41  unsigned lastIndex=0;
42 
43  // LocalKey: entity number , entity codim, dof indices within each entity
44  /* edge and vertex numbering
45  2----3----3
46  | |
47  | |
48  0 1
49  | |
50  | |
51  0----2----1
52  */
53 
54  // lower edge (2)
55  subEntity[lastIndex++] = 0; // corner 0
56  for (unsigned i = 0; i < k - 1; ++i)
57  subEntity[lastIndex++] = 2; // inner dofs of lower edge (2)
58 
59  subEntity[lastIndex++] = 1; // corner 1
60 
61  // iterate from bottom to top over inner edge dofs
62  for (unsigned e = 0; e < k - 1; ++e) {
63  subEntity[lastIndex++] = 0; // left edge (0)
64  for (unsigned i = 0; i < k - 1; ++i)
65  subEntity[lastIndex++] = 0; // face dofs
66  subEntity[lastIndex++] = 1; // right edge (1)
67  }
68 
69  // upper edge (3)
70  subEntity[lastIndex++] = 2; // corner 2
71  for (unsigned i = 0; i < k - 1; ++i)
72  subEntity[lastIndex++] = 3; // inner dofs of upper edge (3)
73 
74  subEntity[lastIndex++] = 3; // corner 3
75 
76  assert((StaticPower<k+1,d>::power==lastIndex));
77  }
78 
79 
80  void setup3d(std::vector<unsigned int>& subEntity)
81  {
82  assert(k>0);
83  unsigned lastIndex=0;
84 #ifndef NDEBUG
85  const unsigned numIndices = StaticPower<k+1,d>::power;
86  const unsigned numFaceIndices = StaticPower<k+1,d-1>::power;
87 #endif
88  const unsigned numInnerEdgeDofs = k-1;
89 
90  // LocalKey: entity number , entity codim, dof indices within each entity
91  /* edge and vertex numbering
92  2----3----3
93  | |
94  | |
95  0 1
96  | |
97  | |
98  0----2----1
99  */
100 
101  // bottom face (4)
102  lastIndex=0;
103  // lower edge (2)
104  subEntity[lastIndex++] = 0; // corner 0
105  for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
106  subEntity[lastIndex++] = 6; // inner dofs of lower edge (2)
107 
108  subEntity[lastIndex++] = 1; // corner 1
109 
110  // iterate from bottom to top over inner edge dofs
111  for (unsigned e = 0; e < numInnerEdgeDofs; ++e) {
112  subEntity[lastIndex++] = 4; // left edge (4)
113  for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
114  subEntity[lastIndex++] = 4; // inner face dofs
115  subEntity[lastIndex++] = 5; // right edge (5)
116  }
117 
118  // upper edge (7)
119  subEntity[lastIndex++] = 2; // corner 2
120  for (unsigned i = 0; i < k - 1; ++i)
121  subEntity[lastIndex++] = 7; // inner dofs of upper edge (7)
122  subEntity[lastIndex++] = 3; // corner 3
123 
124  assert(numFaceIndices==lastIndex); // added 1 face so far
126 
128  for(unsigned f = 0; f < numInnerEdgeDofs; ++f) {
129 
130  // lower edge (connecting edges 0 and 1)
131  subEntity[lastIndex++] = 0; // dof on edge 0
132  for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
133  subEntity[lastIndex++] = 2; // dof in front face
134  subEntity[lastIndex++] = 1; // dof on edge 1
135 
136  // iterate from bottom to top over inner edge dofs
137  for (unsigned e = 0; e < numInnerEdgeDofs; ++e) {
138  subEntity[lastIndex++] = 0; // on left face (0)
139  for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
140  subEntity[lastIndex++] = 0; // volume dofs
141  subEntity[lastIndex++] = 1; // right edge (5)
142  }
143 
144  // upper edge (connecting edges 0 and 1)
145  subEntity[lastIndex++] = 2; // dof on edge 2
146  for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
147  subEntity[lastIndex++] = 3; // dof on rear face (3)
148  subEntity[lastIndex++] = 3; // dof on edge 3
149 
150  assert(lastIndex==(f+1+1)*numFaceIndices);
151  }
152 
154  // lower edge (10)
155  subEntity[lastIndex++] = 4; // corner 4
156  for (unsigned i = 0; i < k - 1; ++i)
157  subEntity[lastIndex++] = 10; // inner dofs on lower edge (10)
158  subEntity[lastIndex++] = 5; // corner 5
159 
160  // iterate from bottom to top over inner edge dofs
161  for (unsigned e = 0; e < k - 1; ++e) {
162  subEntity[lastIndex++] = 8; // left edge (8)
163  for (unsigned i = 0; i < k - 1; ++i)
164  subEntity[lastIndex++] = 5; // face dofs
165  subEntity[lastIndex++] = 9; // right edge (9)
166  }
167 
168  // upper edge (11)
169  subEntity[lastIndex++] = 6; // corner 6
170  for (unsigned i = 0; i < k - 1; ++i)
171  subEntity[lastIndex++] = 11; // inner dofs of upper edge (11)
172  subEntity[lastIndex++] = 7; // corner 7
173 
174  assert(numIndices==lastIndex);
175  }
176 
177  public:
179  QkLocalCoefficients () : li(StaticPower<k+1,d>::power)
180  {
181  // Set up array of codimension-per-dof-number
182  std::vector<unsigned int> codim(li.size());
183 
184  for (std::size_t i=0; i<codim.size(); i++) {
185  codim[i] = 0;
186  // Codimension gets increased by 1 for each coordinate direction
187  // where dof is on boundary
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)
191  codim[i]++;
192  }
193 
194  // Set up index vector (the index of the dof in the set of dofs of a given subentity)
195  // Algorithm: the 'index' has the same ordering as the dof number 'i'.
196  // To make it consecutive we interpret 'i' in the (k+1)-adic system, omit all digits
197  // that correspond to axes where the dof is on the element boundary, and transform the
198  // rest to the (k-1)-adic system.
199  std::vector<unsigned int> index(size());
200 
201  for (std::size_t i=0; i<size(); i++) {
202 
203  index[i] = 0;
204 
205  array<unsigned int,d> mIdx = multiindex(i);
206 
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);
210 
211  }
212 
213  // Set up entity and dof numbers for each (supported) dimension separately
214  std::vector<unsigned int> subEntity(li.size());
215 
216  if (k==1) {
217 
218  for (std::size_t i=0; i<size(); i++)
219  subEntity[i] = i;
220 
221  } else if (d==2) {
222 
223  setup2d(subEntity);
224 
225  } else if (d==3) {
226 
227  setup3d(subEntity);
228 
229  } else
230  DUNE_THROW(Dune::NotImplemented, "QkLocalCoefficients for k==" << k << " and d==" << d);
231 
232  for (size_t i=0; i<li.size(); i++)
233  li[i] = LocalKey(subEntity[i], codim[i], index[i]);
234  }
235 
237  std::size_t size () const
238  {
239  return StaticPower<k+1,d>::power;
240  }
241 
243  const LocalKey& localKey (std::size_t i) const
244  {
245  return li[i];
246  }
247 
248  private:
249  std::vector<LocalKey> li;
250  };
251 
252 }
253 
254 #endif