Dune Core Modules (2.6.0)

qklocalcoefficients.hh
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 <array>
8#include <cassert>
9#include <vector>
10
12#include <dune/common/power.hh>
13
14#include <dune/localfunctions/common/localkey.hh>
15
16namespace Dune
17{
23 template<int k, int d>
25
26 static const unsigned unsignedK = k;
27
28 // Return i as a d-digit number in the (k+1)-nary system
29 static std::array<unsigned int,d> multiindex (unsigned int i)
30 {
31 std::array<unsigned int,d> alpha;
32 for (int j=0; j<d; j++)
33 {
34 alpha[j] = i % (k+1);
35 i = i/(k+1);
36 }
37 return alpha;
38 }
39
41 void setup1d(std::vector<unsigned int>& subEntity)
42 {
43 // Special-handling for piecewise constant elements
44 if (k==0)
45 {
46 subEntity[0] = 0;
47 return;
48 }
49
50 unsigned lastIndex=0;
51
52 /* edge and vertex numbering
53 0----0----1
54 */
55
56 // edge (0)
57 subEntity[lastIndex++] = 0; // corner 0
58 for (unsigned i = 0; i < unsignedK - 1; ++i)
59 subEntity[lastIndex++] = 0; // inner dofs of element (0)
60
61 subEntity[lastIndex++] = 1; // corner 1
62
63 assert((StaticPower<k+1,d>::power==lastIndex));
64 }
65
66 void setup2d(std::vector<unsigned int>& subEntity)
67 {
68 // Special-handling for piecewise constant elements
69 if (k==0)
70 {
71 subEntity[0] = 0;
72 return;
73 }
74
75 unsigned lastIndex=0;
76
77 // LocalKey: entity number, entity codim, dof indices within each entity
78 /* edge and vertex numbering
79 2----3----3
80 | |
81 | |
82 0 1
83 | |
84 | |
85 0----2----1
86 */
87
88 // lower edge (2)
89 subEntity[lastIndex++] = 0; // corner 0
90 for (unsigned i = 0; i < unsignedK - 1; ++i)
91 subEntity[lastIndex++] = 2; // inner dofs of lower edge (2)
92
93 subEntity[lastIndex++] = 1; // corner 1
94
95 // iterate from bottom to top over inner edge dofs
96 for (unsigned e = 0; e < unsignedK - 1; ++e) {
97 subEntity[lastIndex++] = 0; // left edge (0)
98 for (unsigned i = 0; i < unsignedK - 1; ++i)
99 subEntity[lastIndex++] = 0; // face dofs
100 subEntity[lastIndex++] = 1; // right edge (1)
101 }
102
103 // upper edge (3)
104 subEntity[lastIndex++] = 2; // corner 2
105 for (unsigned i = 0; i < unsignedK - 1; ++i)
106 subEntity[lastIndex++] = 3; // inner dofs of upper edge (3)
107
108 subEntity[lastIndex++] = 3; // corner 3
109
110 assert((StaticPower<k+1,d>::power==lastIndex));
111 }
112
113
114
115 void setup3d(std::vector<unsigned int>& subEntity)
116 {
117 // Special-handling for piecewise constant elements
118 if (k==0)
119 {
120 subEntity[0] = 0;
121 return;
122 }
123
124 unsigned lastIndex=0;
125#ifndef NDEBUG
126 const unsigned numIndices = StaticPower<k+1,d>::power;
127 const unsigned numFaceIndices = StaticPower<k+1,d-1>::power;
128#endif
129 const unsigned numInnerEdgeDofs = k-1;
130
131 // LocalKey: entity number, entity codim, dof indices within each entity
132 /* edge and vertex numbering
133
134 6---(11)--7 6---------7
135 /| /| /| (5) /|
136 (8)| (9)| / | top / |
137 / (2) / (3) / |(3)bac/k |
138 4---(10)--5 | 4---------5 |
139 | | | | left|(0)| |(1)|right
140 | 2--(7)|---3 | 2-----|---3
141 (0) / (1) / |(2)front | /
142 |(4) |(5) | / (4) | /
143 |/ |/ |/ bottom |/
144 0---(6)---1 0---------1
145 */
146
147 // bottom face (4)
148 lastIndex=0;
149 // lower edge (6)
150 subEntity[lastIndex++] = 0; // corner 0
151 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
152 subEntity[lastIndex++] = 6; // inner dofs of lower edge (6)
153
154 subEntity[lastIndex++] = 1; // corner 1
155
156 // iterate from bottom to top over inner edge dofs
157 for (unsigned e = 0; e < numInnerEdgeDofs; ++e) {
158 subEntity[lastIndex++] = 4; // left edge (4)
159 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
160 subEntity[lastIndex++] = 4; // inner face dofs
161 subEntity[lastIndex++] = 5; // right edge (5)
162 }
163
164 // upper edge (7)
165 subEntity[lastIndex++] = 2; // corner 2
166 for (unsigned i = 0; i < unsignedK - 1; ++i)
167 subEntity[lastIndex++] = 7; // inner dofs of upper edge (7)
168 subEntity[lastIndex++] = 3; // corner 3
169
170 assert(numFaceIndices==lastIndex); // added 1 face so far
172
174 for(unsigned f = 0; f < numInnerEdgeDofs; ++f) {
175
176 // lower edge (connecting edges 0 and 1)
177 subEntity[lastIndex++] = 0; // dof on edge 0
178 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
179 subEntity[lastIndex++] = 2; // dof in front face
180 subEntity[lastIndex++] = 1; // dof on edge 1
181
182 // iterate from bottom to top over inner edge dofs
183 for (unsigned e = 0; e < numInnerEdgeDofs; ++e) {
184 subEntity[lastIndex++] = 0; // on left face (0)
185 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
186 subEntity[lastIndex++] = 0; // volume dofs
187 subEntity[lastIndex++] = 1; // right face (1)
188 }
189
190 // upper edge (connecting edges 0 and 1)
191 subEntity[lastIndex++] = 2; // dof on edge 2
192 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
193 subEntity[lastIndex++] = 3; // dof on rear face (3)
194 subEntity[lastIndex++] = 3; // dof on edge 3
195
196 assert(lastIndex==(f+1+1)*numFaceIndices);
197 }
198
200 // lower edge (10)
201 subEntity[lastIndex++] = 4; // corner 4
202 for (unsigned i = 0; i < unsignedK - 1; ++i)
203 subEntity[lastIndex++] = 10; // inner dofs on lower edge (10)
204 subEntity[lastIndex++] = 5; // corner 5
205
206 // iterate from bottom to top over inner edge dofs
207 for (unsigned e = 0; e < unsignedK - 1; ++e) {
208 subEntity[lastIndex++] = 8; // left edge (8)
209 for (unsigned i = 0; i < unsignedK - 1; ++i)
210 subEntity[lastIndex++] = 5; // face dofs
211 subEntity[lastIndex++] = 9; // right edge (9)
212 }
213
214 // upper edge (11)
215 subEntity[lastIndex++] = 6; // corner 6
216 for (unsigned i = 0; i < unsignedK - 1; ++i)
217 subEntity[lastIndex++] = 11; // inner dofs of upper edge (11)
218 subEntity[lastIndex++] = 7; // corner 7
219
220 assert(numIndices==lastIndex);
221 }
222
223 public:
225 QkLocalCoefficients () : li(StaticPower<k+1,d>::power)
226 {
227 // Set up array of codimension-per-dof-number
228 std::vector<unsigned int> codim(li.size());
229
230 for (std::size_t i=0; i<codim.size(); i++) {
231 codim[i] = 0;
232 if (k==0)
233 continue;
234 // Codimension gets increased by 1 for each coordinate direction
235 // where dof is on boundary
236 std::array<unsigned int,d> mIdx = multiindex(i);
237 for (int j=0; j<d; j++)
238 if (mIdx[j]==0 or mIdx[j]==k)
239 codim[i]++;
240 }
241
242 // Set up index vector (the index of the dof in the set of dofs of a given subentity)
243 // Algorithm: the 'index' has the same ordering as the dof number 'i'.
244 // To make it consecutive we interpret 'i' in the (k+1)-adic system, omit all digits
245 // that correspond to axes where the dof is on the element boundary, and transform the
246 // rest to the (k-1)-adic system.
247 std::vector<unsigned int> index(size());
248
249 for (std::size_t i=0; i<size(); i++) {
250
251 index[i] = 0;
252
253 std::array<unsigned int,d> mIdx = multiindex(i);
254
255 for (int j=d-1; j>=0; j--)
256 if (mIdx[j]>0 and mIdx[j]<unsignedK)
257 index[i] = (k-1)*index[i] + (mIdx[j]-1);
258
259 }
260
261 // Set up entity and dof numbers for each (supported) dimension separately
262 std::vector<unsigned int> subEntity(li.size());
263
264 if (k==1) { // We can handle the first-order case in any dimension
265
266 for (std::size_t i=0; i<size(); i++)
267 subEntity[i] = i;
268
269 } else if (d==1) {
270
271 setup1d(subEntity);
272
273 } else if (d==2) {
274
275 setup2d(subEntity);
276
277 } else if (d==3) {
278
279 setup3d(subEntity);
280
281 } else
282 DUNE_THROW(Dune::NotImplemented, "QkLocalCoefficients for k==" << k << " and d==" << d);
283
284 for (size_t i=0; i<li.size(); i++)
285 li[i] = LocalKey(subEntity[i], codim[i], index[i]);
286 }
287
289 std::size_t size () const
290 {
292 }
293
295 const LocalKey& localKey (std::size_t i) const
296 {
297 return li[i];
298 }
299
300 private:
301 std::vector<LocalKey> li;
302 };
303
304}
305
306#endif
Describe position of one degree of freedom.
Definition: localkey.hh:21
Default exception for dummy implementations.
Definition: exceptions.hh:261
Attaches a shape function to an entity.
Definition: qklocalcoefficients.hh:24
const LocalKey & localKey(std::size_t i) const
get i'th index
Definition: qklocalcoefficients.hh:295
QkLocalCoefficients()
Default constructor.
Definition: qklocalcoefficients.hh:225
std::size_t size() const
number of coefficients
Definition: qklocalcoefficients.hh:289
A few common exception classes.
Various implementations of the power function for run-time and static arguments.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
Dune namespace.
Definition: alignedallocator.hh:10
Calculates m^p at compile time.
Definition: power.hh:20
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 24, 23:30, 2024)