Loading [MathJax]/extensions/tex2jax.js

dune-composites (unstable)

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