dune-composites (2.5.1)

serendipitylocalcoefficients.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_SERENDIPITY_LOCALCOEFFICIENTS_HH
5#define DUNE_LOCALFUNCTIONS_SERENDIPITY_LOCALCOEFFICIENTS_HH
6
7#include <array>
8#include <cassert>
9#include <vector>
10
11#include <dune/common/exceptions.hh>
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 typedef typename Dune::FieldVector<int,d> DuneVector;
27
28 static void multiindex (DuneVector& alpha)
29 {
30 int i = 0;
31 for (int j=0; j<d; j++){
32 i+=alpha[j]*pow(k+1,j);
33 }
34 i++;
35 for (int j=0; j<d; j++)
36 {
37 alpha[j] = i % (k+1);
38 i = i/(k+1);
39 }
40 unsigned int sum = 0;
41 for (unsigned int j=0; j<d; j++){
42 if(alpha[j] > 0 && alpha[j] < k){ sum++; }
43 }
44 if(sum > 1){
45 multiindex(alpha);
46 }
47 }
48
50 void setup1d(std::vector<unsigned int>& subEntity)
51 {
52 // Special-handling for piecewise constant elements
53 if (k==0)
54 {
55 subEntity[0] = 0;
56 return;
57 }
58
59 unsigned lastIndex=0;
60
61 /* edge and vertex numbering
62 0----0----1
63 */
64
65 // lower edge (2)
66 subEntity[lastIndex++] = 0; // corner 0
67 for (unsigned i = 0; i < k - 1; ++i)
68 subEntity[lastIndex++] = 0; // inner dofs of element (0)
69
70 subEntity[lastIndex++] = 1; // corner 1
71
72 assert((size()==lastIndex));
73 }
74
75 void setup2d(std::vector<unsigned int>& subEntity)
76 {
77 // Special-handling for piecewise constant elements
78 if (k==0)
79 {
80 subEntity[0] = 0;
81 return;
82 }
83
84 unsigned lastIndex=0;
85
86 // LocalKey: entity number , entity codim, dof indices within each entity
87 /* edge and vertex numbering
88 2----3----3
89 | |
90 | |
91 0 1
92 | |
93 | |
94 0----2----1
95 */
96
97 // lower edge (2)
98 subEntity[lastIndex++] = 0; // corner 0
99 for (unsigned i = 0; i < k - 1; ++i)
100 subEntity[lastIndex++] = 2; // inner dofs of lower edge (2)
101 subEntity[lastIndex++] = 1; // corner 1
102
103 // iterate from bottom to top over inner edge dofs
104 for (unsigned e = 0; e < k - 1; ++e) {
105 subEntity[lastIndex++] = 0; // left edge (0)
106 subEntity[lastIndex++] = 1; // right edge (1)
107 }
108
109 // upper edge (3)
110 subEntity[lastIndex++] = 2; // corner 2
111 for (unsigned i = 0; i < k - 1; ++i)
112 subEntity[lastIndex++] = 3; // inner dofs of upper edge (3)
113 subEntity[lastIndex++] = 3; // corner 3
114
115 assert((size()==lastIndex));
116 }
117
118
119
120 void setup3d(std::vector<unsigned int>& subEntity)
121 {
122 // Special-handling for piecewise constant elements
123 if (k==0)
124 {
125 subEntity[0] = 0;
126 return;
127 }
128
129#ifndef NDEBUG
130 const unsigned numIndices = size();
131#endif
132 const unsigned numInnerEdgeDofs = k-1;
133
134 // LocalKey: entity number , entity codim, dof indices within each entity
135 /* edge and vertex numbering
136
137 6---(11)--7 6---------7
138 /| /| /| (5) /|
139 (8)| (9)| / | top / |
140 / (2) / (3) / |(3)bac/k |
141 4---(10)--5 | 4---------5 |
142 | | | | left|(0)| |(1)|right
143 | 2--(7)|---3 | 2-----|---3
144 (0) / (1) / |(2)front | /
145 |(4) |(5) | / (4) | /
146 |/ |/ |/ bottom |/
147 0---(6)---1 0---------1
148 */
149
150 //============================== bottom face (4) =================================
151 unsigned lastIndex=0;
152 // lower edge (6)
153 subEntity[lastIndex++] = 0; // corner 0
154 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
155 subEntity[lastIndex++] = 6; // inner dofs of lower edge (6)
156 subEntity[lastIndex++] = 1; // corner 1
157
158 // iterate from bottom to top over inner edge dofs
159 for (unsigned e = 0; e < numInnerEdgeDofs; ++e) {
160 subEntity[lastIndex++] = 4; // left edge (4)
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 < numInnerEdgeDofs; ++i)
167 subEntity[lastIndex++] = 7; // inner dofs of upper edge (7)
168 subEntity[lastIndex++] = 3; // corner 3
169
170 //=============================== end bottom face (4) ============================
171
172 //=============================== inner faces ====================================
173 for(unsigned f = 0; f < numInnerEdgeDofs; ++f) {
174 // lower edge (connecting edges 0 and 1)
175 subEntity[lastIndex++] = 0; // dof on edge 0
176 subEntity[lastIndex++] = 1; // dof on edge 1
177
178 // upper edge (connecting edges 0 and 1)
179 subEntity[lastIndex++] = 2; // dof on edge 2
180 subEntity[lastIndex++] = 3; // dof on edge 3
181 }
182 //============================= end inner faces =================================
183
184 //============================= top face (5) ====================================
185 // lower edge (10)
186 subEntity[lastIndex++] = 4; // corner 4
187 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
188 subEntity[lastIndex++] = 10; // inner dofs on lower edge (10)
189 subEntity[lastIndex++] = 5; // corner 5
190
191 // iterate from bottom to top over inner edge dofs
192 for (unsigned e = 0; e < numInnerEdgeDofs; ++e) {
193 subEntity[lastIndex++] = 8; // left edge (8)
194 subEntity[lastIndex++] = 9; // right edge (9)
195 }
196
197 // upper edge (11)
198 subEntity[lastIndex++] = 6; // corner 6
199 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
200 subEntity[lastIndex++] = 11; // inner dofs of upper edge (11)
201 subEntity[lastIndex++] = 7; // corner 7
202 //============================ end top face ====================================
203
204 assert(numIndices==lastIndex);
205 }
206
207 public:
210 {
211 // Set up array of codimension-per-dof-number
212 std::vector<unsigned int> codim(li.size());
213 DuneVector mIdx(0);
214 for (std::size_t i=0; i < size(); i++) {
215 codim[i] = 0;
216 if (k==0)
217 continue;
218 // Codimension gets increased by 1 for each coordinate direction
219 // where dof is on boundary
220 if(i>0){ multiindex(mIdx); }
221 for (int j=0; j<d; j++)
222 if (mIdx[j]==0 or mIdx[j]==k)
223 codim[i]++;
224 }
225
226 // Set up index vector (the index of the dof in the set of dofs of a given subentity)
227 // Algorithm: the 'index' has the same ordering as the dof number 'i'.
228 // To make it consecutive we interpret 'i' in the (k+1)-adic system, omit all digits
229 // that correspond to axes where the dof is on the element boundary, and transform the
230 // rest to the (k-1)-adic system.
231 std::vector<unsigned int> index(size());
232 mIdx = 0;
233 for (std::size_t i=0; i<size(); i++) {
234 index[i] = 0;
235 if(i>0){ multiindex(mIdx); }
236
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);
240 }
241
242 // Set up entity and dof numbers for each (supported) dimension separately
243 std::vector<unsigned int> subEntity(li.size());
244
245 if (k==1) { // We can handle the first-order case in any dimension
246 for (std::size_t i=0; i<size(); i++)
247 subEntity[i] = i;
248 }
249 else if (d==1) {
250 setup1d(subEntity);
251 }
252 else if (d==2) {
253 setup2d(subEntity);
254 }
255 else if (d==3) {
256 setup3d(subEntity);
257 }
258 else
259 DUNE_THROW(Dune::NotImplemented, "SerendipityLocalCoefficients for k==" << k << " and d==" << d);
260
261 for (size_t i=0; i<li.size(); i++){
262 li[i] = LocalKey(subEntity[i], codim[i], index[i]);
263 }
264 }
265
267 std::size_t size () const
268 {
269 if(d == 1){
270 return k+1;
271 }
272 if(d == 2){
273 return 2*(k+1)+2*(k-1);
274 }
275 if(d == 3){
276 return 4*(k+1)+8*(k-1);
277 }
278 }
279
281 const LocalKey& localKey (std::size_t i) const
282 {
283 return li[i];
284 }
285
286 private:
287 std::vector<LocalKey> li;
288 };
289
290}
291
292#endif
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Sep 4, 22:38, 2025)