Dune Core Modules (2.9.0)

nedelec1stkindcube.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3// SPDX-FileCopyrightInfo: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_LOCALFUNCTIONS_NEDELEC_NEDELEC1STKINDCUBE_HH
6#define DUNE_LOCALFUNCTIONS_NEDELEC_NEDELEC1STKINDCUBE_HH
7
8#include <numeric>
9
12#include <dune/common/math.hh>
13
14#include <dune/geometry/referenceelements.hh>
15#include <dune/geometry/type.hh>
16
17#include <dune/localfunctions/common/localbasis.hh>
18#include <dune/localfunctions/common/localfiniteelementtraits.hh>
19#include <dune/localfunctions/common/localinterpolation.hh> // For deprecated makeFunctionWithCallOperator
20#include <dune/localfunctions/common/localkey.hh>
21
22namespace Dune
23{
24namespace Impl
25{
36 template<class D, class R, int dim, int k>
37 class Nedelec1stKindCubeLocalBasis
38 {
39 // Number of edges of the reference cube
40 constexpr static std::size_t numberOfEdges = power(2,dim-1)*dim;
41
42 public:
43 using Traits = LocalBasisTraits<D,dim,FieldVector<D,dim>,
44 R,dim,FieldVector<R,dim>,
45 FieldMatrix<R,dim,dim> >;
46
53 Nedelec1stKindCubeLocalBasis()
54 {
55 std::fill(edgeOrientation_.begin(), edgeOrientation_.end(), 1.0);
56 }
57
60 Nedelec1stKindCubeLocalBasis(std::bitset<numberOfEdges> edgeOrientation)
61 : Nedelec1stKindCubeLocalBasis()
62 {
63 for (std::size_t i=0; i<edgeOrientation_.size(); i++)
64 edgeOrientation_[i] *= edgeOrientation[i] ? -1.0 : 1.0;
65 }
66
68 static constexpr unsigned int size()
69 {
70 static_assert(dim==2 || dim==3, "Nedelec shape functions are implemented only for 2d and 3d cubes.");
71 if (dim==2)
72 return 2*k * (k+1);
73 if (dim==3)
74 return 3*k * (k+1) * (k+1);
75 }
76
82 void evaluateFunction(const typename Traits::DomainType& in,
83 std::vector<typename Traits::RangeType>& out) const
84 {
85 static_assert(k==1, "Evaluating Nédélec shape functions is implemented only for first order.");
86 out.resize(size());
87
88 if (dim==2)
89 {
90 // First-order Nédélec shape functions on a square are of the form
91 //
92 // (a, b)^T + (c y, d x)^T, a, b, c, d \in R
93 //
94 // The following coefficients create the four basis vectors
95 // that are dual to the edge degrees of freedom:
96 //
97 // a[0] = 0 b[0] = 1 c[0] = 0 d[0] = -1
98 // a[1] = 0 b[1] = 0 c[1] = 0 d[1] = 1
99 // a[2] = 1 b[2] = 0 c[2] = 0 d[2] = -1
100 // a[3] = 0 b[3] = 0 c[3] = 0 d[3] = 1
101
102 out[0] = { 0, D(1) - in[0]};
103 out[1] = { 0, in[0]};
104 out[2] = { D(1) - in[1], 0};
105 out[3] = { in[1], 0};
106 }
107
108 if constexpr (dim==3)
109 {
110 // First-order Nédélec shape functions on a cube are of the form
111 //
112 // (e1 yz)
113 // a + b x + c y + d z + (e2 xz) , a, b, c, d \in R^3 and b[0]=c[1]=d[2]=0
114 // (e3 xy)
115 //
116 // The following coefficients create the twelve basis vectors
117 // that are dual to the edge degrees of freedom:
118 //
119 // a[0] = { 0, 0, 1} b[0] = { 0, 0, -1} c[0] = { 0, 0, -1} d[0] = { 0, 0, 0} e[0] = { 0, 0, 1}
120 // a[1] = { 0, 0, 0} b[1] = { 0, 0, 1} c[1] = { 0, 0, 0} d[1] = { 0, 0, 0} e[1] = { 0, 0, -1}
121 // a[2] = { 0, 0, 0} b[2] = { 0, 0, 0} c[2] = { 0, 0, 1} d[2] = { 0, 0, 0} e[2] = { 0, 0, -1}
122 // a[3] = { 0, 0, 0} b[3] = { 0, 0, 0} c[3] = { 0, 0, 0} d[3] = { 0, 0, 0} e[3] = { 0, 0, 1}
123 //
124 // The following implementation uses these values, and simply
125 // skips all the zeros.
126
127 for (std::size_t i=0; i<out.size(); i++)
128 out[i] = {0,0,0};
129
130 out[0][2] = { 1 - in[0] - in[1] + in[0]*in[1]};
131 out[1][2] = { in[0] - in[0]*in[1]};
132 out[2][2] = { in[1] - in[0]*in[1]};
133 out[3][2] = { in[0]*in[1]};
134
135 out[4][1] = { 1 - in[0] - in[2] + in[0]*in[2]};
136 out[5][1] = { in[0] - in[0]*in[2]};
137 out[8][1] = { in[2] - in[0]*in[2]};
138 out[9][1] = { in[0]*in[2]};
139
140 out[6][0] = { 1 - in[1] - in[2] + in[1]*in[2]};
141 out[7][0] = { in[1] - in[1]*in[2]};
142 out[10][0] = { in[2] - in[1]*in[2]};
143 out[11][0] = { in[1]*in[2]};
144 }
145
146 for (std::size_t i=0; i<out.size(); i++)
147 out[i] *= edgeOrientation_[i];
148 }
149
155 void evaluateJacobian(const typename Traits::DomainType& in,
156 std::vector<typename Traits::JacobianType>& out) const
157 {
158 out.resize(size());
159 if (dim==2)
160 {
161 for (std::size_t i=0; i<out.size(); i++)
162 for (std::size_t j=0; j<dim; j++)
163 out[i][j] = { 0, 0};
164
165 out[0][1] = { -1, 0};
166 out[1][1] = { 1, 0};
167
168 out[2][0] = { 0, -1};
169 out[3][0] = { 0, 1};
170 }
171 if (dim==3)
172 {
173 for (std::size_t i=0; i<out.size(); i++)
174 for(std::size_t j=0;j<dim; j++)
175 out[i][j] = {0,0,0};
176
177
178 out[0][2] = {-1 +in[1], -1 + in[0], 0};
179 out[1][2] = { 1 -in[1], - in[0], 0};
180 out[2][2] = { -in[1], 1 - in[0], 0};
181 out[3][2] = { in[1], in[0], 0};
182
183 out[4][1] = {-1 +in[2], 0, -1 + in[0]};
184 out[5][1] = { 1 -in[2], 0, - in[0]};
185 out[8][1] = { -in[2], 0, 1 - in[0]};
186 out[9][1] = { in[2], 0, in[0]};
187
188 out[6][0] = { 0, -1 + in[2], -1 + in[1]};
189 out[7][0] = { 0, 1 - in[2], - in[1]};
190 out[10][0] = { 0, - in[2], 1 - in[1]};
191 out[11][0] = { 0, in[2], in[1]};
192
193 }
194
195 for (std::size_t i=0; i<out.size(); i++)
196 out[i] *= edgeOrientation_[i];
197
198 }
199
206 void partial(const std::array<unsigned int, dim>& order,
207 const typename Traits::DomainType& in,
208 std::vector<typename Traits::RangeType>& out) const
209 {
210 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
211 if (totalOrder == 0) {
212 evaluateFunction(in, out);
213 } else if (totalOrder == 1) {
214 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
215 out.resize(size());
216
217 if (dim==2)
218 {
219 if (direction==0)
220 {
221 out[0] = { 0, -1};
222 out[1] = { 0, 1};
223 out[2] = { 0, 0};
224 out[3] = { 0, 0};
225 }
226 else
227 {
228 out[0] = { 0, 0};
229 out[1] = { 0, 0};
230 out[2] = { -1, 0};
231 out[3] = { 1, 0};
232 }
233 }
234
235 if (dim==3)
236 {
237 switch (direction)
238 {
239 case 0:
240 out[0] = { 0, 0, -1 +in[1]};
241 out[1] = { 0, 0, 1 -in[1]};
242 out[2] = { 0, 0, -in[1]};
243 out[3] = { 0, 0, in[1]};
244
245 out[4] = { 0, -1 +in[2], 0};
246 out[5] = { 0, 1 -in[2], 0};
247 out[8] = { 0, -in[2], 0};
248 out[9] = { 0, in[2], 0};
249
250 out[6] = {0,0,0};
251 out[7] = {0,0,0};
252 out[10] = {0,0,0};
253 out[11] = {0,0,0};
254 break;
255
256 case 1:
257 out[0] = { 0, 0, -1 + in[0]};
258 out[1] = { 0, 0, - in[0]};
259 out[2] = { 0, 0, 1 - in[0]};
260 out[3] = { 0, 0, in[0]};
261
262 out[4] = {0,0,0};
263 out[5] = {0,0,0};
264 out[8] = {0,0,0};
265 out[9] = {0,0,0};
266
267 out[6] = { -1 + in[2], 0, 0};
268 out[7] = { 1 - in[2], 0, 0};
269 out[10] = { - in[2], 0, 0};
270 out[11] = { in[2], 0, 0};
271 break;
272
273 case 2:
274 out[0] = {0,0,0};
275 out[1] = {0,0,0};
276 out[2] = {0,0,0};
277 out[3] = {0,0,0};
278
279 out[4] = { 0, -1 + in[0], 0};
280 out[5] = { 0, - in[0], 0};
281 out[8] = { 0, 1 - in[0], 0};
282 out[9] = { 0, in[0], 0};
283
284 out[6] = { -1 + in[1], 0, 0};
285 out[7] = { - in[1], 0, 0};
286 out[10] = { 1 - in[1], 0, 0};
287 out[11] = { in[1], 0, 0};
288 break;
289 }
290 }
291
292 for (std::size_t i=0; i<out.size(); i++)
293 out[i] *= edgeOrientation_[i];
294
295 } else if (totalOrder == 2) {
296 out.resize(size());
297
298 if (dim==2)
299 for (std::size_t i = 0; i < size(); ++i)
300 for (std::size_t j = 0; j < dim; ++j)
301 out[i][j] = 0;
302
303 if (dim==3)
304 {
305 for(size_t i=0; i<out.size(); i++)
306 out[i] = { 0, 0, 0};
307
308 //case (1,1,0):
309 if( order[0] == 1 and order[1]==1)
310 {
311 out[0] = { 0, 0, 1};
312 out[1] = { 0, 0, -1};
313 out[2] = { 0, 0, -1};
314 out[3] = { 0, 0, 1};
315 }
316
317 //case (1,0,1):
318 if( order[0] == 1 and order[2]==1)
319 {
320 out[4] = { 0, 1, 0};
321 out[5] = { 0, -1, 0};
322 out[8] = { 0, -1, 0};
323 out[9] = { 0, 1, 0};
324 }
325
326 //case (0,1,1):
327 if( order[1] == 1 and order[2]==1)
328 {
329 out[6] = { 1, 0, 0};
330 out[7] = { -1, 0, 0};
331 out[10] = { -1, 0, 0};
332 out[11] = { 1, 0, 0};
333 }
334
335 for (std::size_t i=0; i<out.size(); i++)
336 out[i] *= edgeOrientation_[i];
337 }
338
339
340 }else {
341 out.resize(size());
342 for (std::size_t i = 0; i < size(); ++i)
343 for (std::size_t j = 0; j < dim; ++j)
344 out[i][j] = 0;
345 }
346
347 }
348
350 unsigned int order() const
351 {
352 if (dim==2)
353 return 2*k-1;
354 if (dim==3)
355 return 3*k-1;
356 }
357
358 private:
359
360 // Orientations of the cube edges
361 std::array<R,numberOfEdges> edgeOrientation_;
362 };
363
364
369 template <int dim, int k>
370 class Nedelec1stKindCubeLocalCoefficients
371 {
372 public:
374 Nedelec1stKindCubeLocalCoefficients ()
375 : localKey_(size())
376 {
377 static_assert(k==1, "Only first-order Nédélec local coefficients are implemented.");
378 // Assign all degrees of freedom to edges
379 // TODO: This is correct only for first-order Nédélec elements
380 for (std::size_t i=0; i<size(); i++)
381 localKey_[i] = LocalKey(i,dim-1,0);
382 }
383
385 std::size_t size() const
386 {
387 static_assert(dim==2 || dim==3, "Nédélec shape functions are implemented only for 2d and 3d cubes.");
388 return (dim==2) ? 2*k * (k+1)
389 : 3*k * (k+1) * (k+1);
390 }
391
394 const LocalKey& localKey (std::size_t i) const
395 {
396 return localKey_[i];
397 }
398
399 private:
400 std::vector<LocalKey> localKey_;
401 };
402
407 template<class LB>
408 class Nedelec1stKindCubeLocalInterpolation
409 {
410 static constexpr auto dim = LB::Traits::dimDomain;
411 static constexpr auto size = LB::size();
412
413 // Number of edges of the reference cube
414 constexpr static std::size_t numberOfEdges = power(2,dim-1)*dim;
415
416 public:
417
419 Nedelec1stKindCubeLocalInterpolation (std::bitset<numberOfEdges> s = 0)
420 {
421 auto refElement = Dune::referenceElement<double,dim>(GeometryTypes::cube(dim));
422
423 for (std::size_t i=0; i<numberOfEdges; i++)
424 m_[i] = refElement.position(i,dim-1);
425
426 for (std::size_t i=0; i<numberOfEdges; i++)
427 {
428 auto vertexIterator = refElement.subEntities(i,dim-1,dim).begin();
429 auto v0 = *vertexIterator;
430 auto v1 = *(++vertexIterator);
431
432 // By default, edges point from the vertex with the smaller index
433 // to the vertex with the larger index.
434 if (v0>v1)
435 std::swap(v0,v1);
436 edge_[i] = refElement.position(v1,dim) - refElement.position(v0,dim);
437 edge_[i] *= (s[i]) ? -1.0 : 1.0;
438 }
439 }
440
446 template<typename F, typename C>
447 void interpolate (const F& ff, std::vector<C>& out) const
448 {
449 out.resize(size);
450 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
451
452 for (std::size_t i=0; i<size; i++)
453 {
454 auto y = f(m_[i]);
455 out[i] = 0.0;
456 for (int j=0; j<dim; j++)
457 out[i] += y[j]*edge_[i][j];
458 }
459 }
460
461 private:
462 // Edge midpoints of the reference cube
463 std::array<typename LB::Traits::DomainType, numberOfEdges> m_;
464 // Edges of the reference cube
465 std::array<typename LB::Traits::DomainType, numberOfEdges> edge_;
466 };
467
468}
469
470
494 template<class D, class R, int dim, int k>
496 {
497 public:
499 Impl::Nedelec1stKindCubeLocalCoefficients<dim,k>,
500 Impl::Nedelec1stKindCubeLocalInterpolation<Impl::Nedelec1stKindCubeLocalBasis<D,R,dim,k> > >;
501
502 static_assert(dim==2 || dim==3, "Nedelec elements are only implemented for 2d and 3d elements.");
503 static_assert(k==1, "Nedelec elements of the first kind are currently only implemented for order k==1.");
504
508
514 Nedelec1stKindCubeLocalFiniteElement (std::bitset<power(2,dim-1)*dim> s) :
515 basis_(s),
516 interpolation_(s)
517 {}
518
519 const typename Traits::LocalBasisType& localBasis () const
520 {
521 return basis_;
522 }
523
524 const typename Traits::LocalCoefficientsType& localCoefficients () const
525 {
526 return coefficients_;
527 }
528
529 const typename Traits::LocalInterpolationType& localInterpolation () const
530 {
531 return interpolation_;
532 }
533
534 static constexpr unsigned int size ()
535 {
536 return Traits::LocalBasisType::size();
537 }
538
539 static constexpr GeometryType type ()
540 {
541 return GeometryTypes::cube(dim);
542 }
543
544 private:
545 typename Traits::LocalBasisType basis_;
546 typename Traits::LocalCoefficientsType coefficients_;
547 typename Traits::LocalInterpolationType interpolation_;
548 };
549
550}
551
552#endif
Nédélec elements of the first kind for cube elements.
Definition: nedelec1stkindcube.hh:496
Nedelec1stKindCubeLocalFiniteElement()=default
Default constructor.
Nedelec1stKindCubeLocalFiniteElement(std::bitset< power(2, dim-1) *dim > s)
Constructor with explicitly given edge orientations.
Definition: nedelec1stkindcube.hh:514
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:472
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:291
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:13
constexpr Base power(Base m, Exponent p)
Power method for integer exponents.
Definition: math.hh:75
D DomainType
domain type
Definition: localbasis.hh:42
traits helper struct
Definition: localfiniteelementtraits.hh:13
LB LocalBasisType
Definition: localfiniteelementtraits.hh:16
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:20
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:24
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)