DUNE PDELab (unstable)

lagrangesimplex.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 © 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_LAGRANGE_LAGRANGESIMPLEX_HH
6#define DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGESIMPLEX_HH
7
8#include <array>
9#include <cstddef>
10#include <iterator>
11#include <numeric>
12#include <algorithm>
13#include <span>
14#include <type_traits>
15
19#include <dune/common/hybridutilities.hh>
20#include <dune/common/math.hh>
22#include <dune/common/std/mdarray.hh>
23#include <dune/common/std/mdspan.hh>
24
25#include <dune/geometry/referenceelements.hh>
26
27#include <dune/localfunctions/common/localbasis.hh>
28#include <dune/localfunctions/common/localfiniteelementtraits.hh>
29#include <dune/localfunctions/common/localkey.hh>
30
31namespace Dune { namespace Impl
32{
33 // The traits provide static or dynamic order and size information
34 template<unsigned int dim, int compileTimeOrder>
35 struct LagrangeSimplexTraits;
36
37 // The traits provide static order and size information
38 template<unsigned int dim, int compileTimeOrder>
39 requires (compileTimeOrder >= 0)
40 struct LagrangeSimplexTraits<dim,compileTimeOrder>
41 {
42 static constexpr bool is_static_order = true;
43
44 constexpr LagrangeSimplexTraits(int /*runTimeOrder*/ = compileTimeOrder) {}
45
47 static constexpr unsigned int order() { return compileTimeOrder; }
48
53 static constexpr unsigned int size() { return binomial(compileTimeOrder+dim,dim); }
54 };
55
56
57 // this specialization is for the dynamic order case
58 template<unsigned int dim>
59 struct LagrangeSimplexTraits<dim,-1>
60 {
61 unsigned int order_;
62 std::size_t size_;
63
64 static constexpr bool is_static_order = false;
65
66 constexpr explicit LagrangeSimplexTraits(int runTimeOrder)
67 : order_(runTimeOrder >= 0 ? (unsigned int)(runTimeOrder) : 0u)
68 , size_(binomial(order_+dim,dim))
69 {
70 if (runTimeOrder < 0)
71 DUNE_THROW(Dune::InvalidStateException, "LagrangeSimplex: order must be non-negative");
72 }
73
75 constexpr unsigned int order() const
76 {
77 return order_;
78 }
79
81 constexpr std::size_t size() const
82 {
83 return size_;
84 }
85 };
86
87
96 template<class R, unsigned int dim, int k>
97 struct LagrangeSimplexLocalBasisBuffers
98 {
99 LagrangeSimplexLocalBasisBuffers(int /*order*/) {}
100
101 // Buffer for the evaluation of shape function, used in evaluateFunction()
102 static auto makeValueBuffer(int /*order*/)
103 {
104 using E = Std::extents<int,dim+1,k+1>;
105 using C = std::array<R,(dim+1)*(k+1)>;
106 return Std::mdarray<R,E,Std::layout_right,C>{};
107 }
108
109 // Buffer for the evaluation of shape function Jacobians, used in evaluateJacobian()
110 static auto makeJacobianBuffer(int /*order*/)
111 {
112 using E = Std::extents<int,dim+1,2,k+1>;
113 using C = std::array<R,(dim+1)*2*(k+1)>;
114 return Std::mdarray<R,E,Std::layout_right,C>{};
115 }
116
117 // Buffer for the evaluation of shape function partial derivatives, used in partial()
118 template <class T, T partialOrders>
119 static auto makePartialBuffer(std::integral_constant<T,partialOrders>, int /*order*/)
120 {
121 using E = Std::extents<int,dim+1,partialOrders+1,k+1>;
122 using C = std::array<R,(dim+1)*(partialOrders+1)*(k+1)>;
123 return Std::mdarray<R,E,Std::layout_right,C>{};
124 }
125 };
126
127
135 template<class R, unsigned int dim>
136 struct LagrangeSimplexLocalBasisBuffers<R,dim,-1>
137 {
138 LagrangeSimplexLocalBasisBuffers(int order)
139 : buffer_((dim+1)*std::max<int>(2,order+1)*(order+1))
140 {}
141
142 mutable std::vector<R> buffer_{};
143
144 // Buffer for the evaluation of shape function, used in evaluateFunction()
145 auto makeValueBuffer(int order) const
146 {
147 using E = Std::extents<int,dim+1,std::dynamic_extent>;
148 return Std::mdspan<R,E>{buffer_.data(), E(order+1)};
149 }
150
151 // Buffer for the evaluation of shape function Jacobians, used in evaluateJacobian()
152 auto makeJacobianBuffer(int order) const
153 {
154 using E = Std::extents<int,dim+1,2,std::dynamic_extent>;
155 return Std::mdspan<R,E>{buffer_.data(), E(order+1)};
156 }
157
158 // Buffer for the evaluation of shape function partial derivatives, used in partial()
159 auto makePartialBuffer(int partialOrders, int order) const
160 {
161 using E = Std::extents<int,dim+1,std::dynamic_extent,std::dynamic_extent>;
162 return Std::mdspan<R,E>{buffer_.data(), E(partialOrders+1, order+1)};
163 }
164 };
165
166
177 template<class D, class R, unsigned int dim, int compileTimeOrder>
178 class LagrangeSimplexLocalBasis
179 : private LagrangeSimplexTraits<dim,compileTimeOrder>
180 , private LagrangeSimplexLocalBasisBuffers<R,dim,compileTimeOrder>
181 {
182 using OrderTraits = LagrangeSimplexTraits<dim,compileTimeOrder>;
183 using Buffers = LagrangeSimplexLocalBasisBuffers<R,dim,compileTimeOrder>;
184 static constexpr bool is_static_order = OrderTraits::is_static_order;
185
186 public:
188 constexpr LagrangeSimplexLocalBasis () requires (OrderTraits::is_static_order)
189 : LagrangeSimplexLocalBasis(OrderTraits())
190 {}
191
193 explicit constexpr LagrangeSimplexLocalBasis(OrderTraits orderTraits)
194 : OrderTraits(orderTraits)
195 , Buffers(orderTraits.order())
196 {}
197
198 using OrderTraits::order;
199 using OrderTraits::size;
200
201 private:
202
203 // Fix the first index in a multi-dimensional array to `i`. For matrices this
204 // corresponds to the ith row of the matrix.
205 // Equivalent to std::submdspan(L, i, std::full_extent...)
206 template <class I, std::size_t e0, std::size_t... ee, class A>
207 static auto slice(Std::mdspan<R,Std::extents<I,e0,ee...>,Std::layout_right,A> L, int i)
208 {
209 return Dune::unpackIntegerSequence([&](auto... ii) {
210 return Std::mdspan<R,Std::extents<I,ee...>,Std::layout_right>{
211 L.data_handle() + i * L.stride(0), L.extent(ii+1)...};
212 }, std::make_index_sequence<sizeof...(ee)>{});
213 }
214
215 // Overload for mdarray type: forward to the mdspan implementation
216 template <class I, std::size_t e0, std::size_t... ee, class C>
217 static auto slice(Std::mdarray<R,Std::extents<I,e0,ee...>,Std::layout_right,C>& L, int i)
218 {
219 return slice(L.to_mdspan(), i);
220 }
221
222 // Compute the rescaled barycentric coordinates of x.
223 // We rescale the simplex by k and then compute the
224 // barycentric coordinates with respect to the points
225 // p_i = e_i (for i=0,...,dim-1) and p_dim=0.
226 // Notice that then the Lagrange points have the barycentric
227 // coordinates (i_0,...,i_d) where i_j are all non-negative
228 // integers satisfying the constraint sum i_j = k.
229 constexpr auto barycentric(const auto& x) const
230 {
231 auto b = std::array<R,dim+1>{};
232 b[dim] = order();
233 for (auto i : Dune::range(dim))
234 {
235 b[i] = order() * x[i];
236 b[dim] -= b[i];
237 }
238 return b;
239 }
240
241 // Evaluate the univariate Lagrange polynomials L_i(t) for i=0,...,k where
242 //
243 // L_i(t) = (t-0)/(i-0) * ... * (t-(i-1))/(i-(i-1))
244 // = (t-0)*...*(t-(i-1))/(i!);
245 static constexpr void evaluateLagrangePolynomials(const R& t, auto&& L, int k)
246 {
247 L[0] = 1;
248 for (auto i : Dune::range(k))
249 L[i+1] = L[i] * (t - i) / (i+1);
250 }
251
252 // Evaluate the univariate Lagrange polynomial derivatives L_i(t) for i=0,...,k
253 // up to given maxDerivativeOrder.
254 static constexpr void evaluateLagrangePolynomialDerivative(const R& t, auto&& LL, unsigned int maxDerivativeOrder, int k)
255 {
256 auto L = slice(LL,0);
257 L[0] = 1;
258 for (auto i : Dune::range(k))
259 L[i+1] = L[i] * (t - i) / (i+1);
260 for(auto j : Dune::range(maxDerivativeOrder))
261 {
262 auto F = slice(LL,j);
263 auto DF = slice(LL,j+1);
264 DF[0] = 0;
265 for (auto i : Dune::range(k))
266 DF[i+1] = (DF[i] * (t - i) + (j+1)*F[i]) / (i+1);
267 }
268 }
269
270 using BarycentricMultiIndex = std::array<int,dim+1>;
271
272
273 // This computed the required partial derivative given by the multi-index
274 // beta of a product of a function given as a product of dim+1 derivatives
275 // of univariate Lagrange polynomials of the dim+1 barycentric coordinates.
276 // The polynomials in the product are specified as follows:
277 //
278 // The table L contains all required derivatives of all univariate
279 // polynomials evaluated at all barycentric coordinates. The two
280 // multi-indices i and alpha encode that the polynomial for the
281 // j-th barycentric coordinate is the alpha-j-th derivative of
282 // the i_j-the Lagrange polynomial.
283 //
284 // Hence this method computes D_beta f(x) where f(x) is the product
285 // \f$f(x) = \prod_{j=0}^{d} L_{i_j}^{(alpha_j)}(x_j) \f$.
286 static constexpr R barycentricDerivative(
287 BarycentricMultiIndex beta,
288 const auto&L,
289 int k,
290 const BarycentricMultiIndex& i,
291 const BarycentricMultiIndex& alpha = {})
292 {
293 // If there are unprocessed derivatives left we search the first unprocessed
294 // partial derivative direction j and compute it using the product and chain rule.
295 // The remaining partial derivatives are forwarded to the recursion.
296 // Notice that the partial derivative of the last barycentric component
297 // wrt to the j-th component is -1 which is responsible for the second
298 // term in the sum. Furthermore we get the factor k due to the scaling of
299 // the simplex.
300 for(auto j : Dune::range(dim))
301 {
302 if (beta[j] > 0)
303 {
304 auto leftDerivatives = alpha;
305 auto rightDerivatives = alpha;
306 leftDerivatives[j]++;
307 rightDerivatives.back()++;
308 beta[j]--;
309 return (barycentricDerivative(beta, L, k, i, leftDerivatives) -
310 barycentricDerivative(beta, L, k, i, rightDerivatives)) * k;
311 }
312 }
313 // If there are no unprocessed derivatives we can simply evaluate
314 // the product of the derivatives of the Lagrange polynomials with
315 // given indices i_j and derivative orders alpha_j
316 // Evaluate the product of the univariate Lagrange polynomials
317 // with given indices and orders.
318 auto y = R(1);
319 for(auto j : Dune::range(dim+1))
320 y *= L(j,alpha[j],i[j]);
321 return y;
322 }
323
324
325 public:
326 using Traits = LocalBasisTraits<D,dim,FieldVector<D,dim>,R,1,FieldVector<R,1>,FieldMatrix<R,1,dim> >;
327
329 void evaluateFunction(const typename Traits::DomainType& x,
330 std::vector<typename Traits::RangeType>& out) const
331 {
332 const int k = order();
333 out.resize(size());
334
335 // Specialization for zero-order case
336 if (k==0)
337 {
338 out[0] = 1;
339 return;
340 }
341
342 // Specialization for first-order case
343 if (k==1)
344 {
345 out[0] = 1.0;
346 for (size_t i=0; i<dim; i++)
347 {
348 out[0] -= x[i];
349 out[i+1] = x[i];
350 }
351 return;
352 }
353
354 // Compute rescaled barycentric coordinates of x
355 auto z = barycentric(x);
356
357 auto L = Buffers::makeValueBuffer(k);
358 for (auto j : Dune::range(dim+1))
359 evaluateLagrangePolynomials(z[j], slice(L,j), k);
360
361 if (dim==1)
362 {
363 unsigned int n = 0;
364 for (auto i0 : Dune::range(k + 1))
365 for (auto i1 : std::array{k - i0})
366 out[n++] = L(0,i0) * L(1,i1);
367 return;
368 }
369 if (dim==2)
370 {
371 unsigned int n=0;
372 for (auto i1 : Dune::range(k + 1))
373 for (auto i0 : Dune::range(k - i1 + 1))
374 for (auto i2 : std::array{k - i1 - i0})
375 out[n++] = L(0,i0) * L(1,i1) * L(2,i2);
376 return;
377 }
378 if (dim==3)
379 {
380 unsigned int n = 0;
381 for (auto i2 : Dune::range(k + 1))
382 for (auto i1 : Dune::range(k - i2 + 1))
383 for (auto i0 : Dune::range(k - i2 - i1 + 1))
384 for (auto i3 : std::array{k - i2 - i1 - i0})
385 out[n++] = L(0,i0) * L(1,i1) * L(2,i2) * L(3,i3);
386 return;
387 }
388
389 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalBasis for k>=2 only implemented for dim<=3");
390 }
391
397 void evaluateJacobian(const typename Traits::DomainType& x,
398 std::vector<typename Traits::JacobianType>& out) const
399 {
400 const int k = order();
401 out.resize(size());
402
403 // Specialization for k==0
404 if (k==0)
405 {
406 std::fill(out[0][0].begin(), out[0][0].end(), 0);
407 return;
408 }
409
410 // Specialization for k==1
411 if (k==1)
412 {
413 std::fill(out[0][0].begin(), out[0][0].end(), -1);
414
415 for (unsigned int i=0; i<dim; i++)
416 for (unsigned int j=0; j<dim; j++)
417 out[i+1][0][j] = (i==j);
418
419 return;
420 }
421
422 // Compute rescaled barycentric coordinates of x
423 auto z = barycentric(x);
424
425 // L(j,m,i) is the m-th derivative of the i-th Lagrange polynomial at z[j]
426 auto L = Buffers::makeJacobianBuffer(k);
427 for (auto j : Dune::range(dim+1))
428 evaluateLagrangePolynomialDerivative(z[j], slice(L,j), 1, k);
429
430 if (dim==1)
431 {
432 unsigned int n = 0;
433 for (auto i0 : Dune::range(k + 1))
434 {
435 for (auto i1 : std::array{k-i0})
436 {
437 out[n][0][0] = (L(0,1,i0) * L(1,0,i1) - L(0,0,i0) * L(1,1,i1))*k;
438 ++n;
439 }
440 }
441 return;
442 }
443 if (dim==2)
444 {
445 unsigned int n=0;
446 for (auto i1 : Dune::range(k + 1))
447 {
448 for (auto i0 : Dune::range(k - i1 + 1))
449 {
450 for (auto i2 : std::array{k - i1 - i0})
451 {
452 out[n][0][0] = (L(0,1,i0) * L(1,0,i1) * L(2,0,i2) - L(0,0,i0) * L(1,0,i1) * L(2,1,i2))*k;
453 out[n][0][1] = (L(0,0,i0) * L(1,1,i1) * L(2,0,i2) - L(0,0,i0) * L(1,0,i1) * L(2,1,i2))*k;
454 ++n;
455 }
456 }
457 }
458 return;
459 }
460 if (dim==3)
461 {
462 unsigned int n = 0;
463 for (auto i2 : Dune::range(k + 1))
464 {
465 for (auto i1 : Dune::range(k - i2 + 1))
466 {
467 for (auto i0 : Dune::range(k - i2 - i1 + 1))
468 {
469 for (auto i3 : std::array{k - i2 - i1 - i0})
470 {
471 out[n][0][0] = (L(0,1,i0) * L(1,0,i1) * L(2,0,i2) * L(3,0,i3) - L(0,0,i0) * L(1,0,i1) * L(2,0,i2) * L(3,1,i3))*k;
472 out[n][0][1] = (L(0,0,i0) * L(1,1,i1) * L(2,0,i2) * L(3,0,i3) - L(0,0,i0) * L(1,0,i1) * L(2,0,i2) * L(3,1,i3))*k;
473 out[n][0][2] = (L(0,0,i0) * L(1,0,i1) * L(2,1,i2) * L(3,0,i3) - L(0,0,i0) * L(1,0,i1) * L(2,0,i2) * L(3,1,i3))*k;
474 ++n;
475 }
476 }
477 }
478 }
479
480 return;
481 }
482
483 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalBasis for k>=2 only implemented for dim<=3");
484 }
485
492 void partial(const std::array<unsigned int,dim>& partialOrders,
493 const typename Traits::DomainType& in,
494 std::vector<typename Traits::RangeType>& out) const
495 {
496 const int k = order();
497 int totalOrder = std::accumulate(partialOrders.begin(), partialOrders.end(), 0u);
498
499 out.resize(size());
500
501 // Derivative order zero corresponds to the function evaluation.
502 if (totalOrder == 0)
503 {
504 evaluateFunction(in, out);
505 return;
506 }
507
508 // Derivatives of order >k are all zero.
509 if (totalOrder > k)
510 {
511 for(auto& out_i : out)
512 out_i = 0;
513 return;
514 }
515
516 // It remains to cover the cases 0 < totalOrder<= k.
517
518 if (k==1)
519 {
520 if (totalOrder==1)
521 {
522 auto direction = std::find(partialOrders.begin(), partialOrders.end(), 1);
523 out[0] = -1;
524 for (unsigned int i=0; i<dim; i++)
525 out[i+1] = (i==(direction-partialOrders.begin()));
526 }
527 return;
528 }
529
530 // Since the required stack storage depends on the dynamic total order,
531 // we need to do a dynamic to static dispatch by enumerating all supported
532 // static orders.
533 auto supportedOrders = [&]{
534 if constexpr(is_static_order)
536 else
537 return Dune::range(1, k+1);
538 }();
539 return Dune::Hybrid::switchCases(supportedOrders, totalOrder, [&](auto staticTotalOrder) {
540
541 // Compute rescaled barycentric coordinates of x
542 auto z = barycentric(in);
543
544 // L(j,m,i) is the m-th derivative of the i-th Lagrange polynomial at z[j]
545 auto L = Buffers::makePartialBuffer(staticTotalOrder,k);
546 for (auto j : Dune::range(dim))
547 evaluateLagrangePolynomialDerivative(z[j], slice(L,j), partialOrders[j], k);
548 evaluateLagrangePolynomialDerivative(z[dim], slice(L,dim), totalOrder, k);
549
550 auto barycentricOrder = BarycentricMultiIndex{};
551 for (auto j : Dune::range(dim))
552 barycentricOrder[j] = partialOrders[j];
553 barycentricOrder[dim] = 0;
554
555 if constexpr (dim==1)
556 {
557 unsigned int n = 0;
558 for (auto i0 : Dune::range(k + 1))
559 for (auto i1 : std::array{k - i0})
560 out[n++] = barycentricDerivative(barycentricOrder, L, k, BarycentricMultiIndex{i0, i1});
561 }
562 if constexpr (dim==2)
563 {
564 unsigned int n=0;
565 for (auto i1 : Dune::range(k + 1))
566 for (auto i0 : Dune::range(k - i1 + 1))
567 for (auto i2 : std::array{k - i1 - i0})
568 out[n++] = barycentricDerivative(barycentricOrder, L, k, BarycentricMultiIndex{i0, i1, i2});
569 }
570 if constexpr (dim==3)
571 {
572 unsigned int n = 0;
573 for (auto i2 : Dune::range(k + 1))
574 for (auto i1 : Dune::range(k - i2 + 1))
575 for (auto i0 : Dune::range(k - i2 - i1 + 1))
576 for (auto i3 : std::array{k - i2 - i1 - i0})
577 out[n++] = barycentricDerivative(barycentricOrder, L, k, BarycentricMultiIndex{i0, i1, i2, i3});
578 }
579 });
580 }
581 };
582
588 template<unsigned int dim, int compileTimeOrder>
589 class LagrangeSimplexLocalCoefficients
590 : private LagrangeSimplexTraits<dim,compileTimeOrder>
591 {
592 using OrderTraits = LagrangeSimplexTraits<dim,compileTimeOrder>;
593
594 public:
596 LagrangeSimplexLocalCoefficients () requires (OrderTraits::is_static_order)
597 : LagrangeSimplexLocalCoefficients(OrderTraits())
598 {}
599
601 LagrangeSimplexLocalCoefficients (OrderTraits orderTraits)
602 : OrderTraits(orderTraits)
603 , localKeys_(OrderTraits::size())
604 {
605 const int k = OrderTraits::order();
606 if (k==0)
607 {
608 localKeys_[0] = LocalKey(0,0,0);
609 return;
610 }
611
612 if (k==1)
613 {
614 for (std::size_t i=0; i<OrderTraits::size(); i++)
615 localKeys_[i] = LocalKey(i,dim,0);
616 return;
617 }
618
619 if (dim==1)
620 {
621 // Order is at least 2 here
622 localKeys_[0] = LocalKey(0,1,0); // vertex dof
623 for (int i=1; i<k; i++)
624 localKeys_[i] = LocalKey(0,0,i-1); // element dofs
625 localKeys_.back() = LocalKey(1,1,0); // vertex dof
626 return;
627 }
628
629 if (dim==2)
630 {
631 int n=0;
632 int c=0;
633 for (int j=0; j<=k; j++)
634 for (int i=0; i<=k-j; i++)
635 {
636 if (i==0 && j==0)
637 {
638 localKeys_[n++] = LocalKey(0,2,0);
639 continue;
640 }
641 if (i==k && j==0)
642 {
643 localKeys_[n++] = LocalKey(1,2,0);
644 continue;
645 }
646 if (i==0 && j==k)
647 {
648 localKeys_[n++] = LocalKey(2,2,0);
649 continue;
650 }
651 if (j==0)
652 {
653 localKeys_[n++] = LocalKey(0,1,i-1);
654 continue;
655 }
656 if (i==0)
657 {
658 localKeys_[n++] = LocalKey(1,1,j-1);
659 continue;
660 }
661 if (i+j==k)
662 {
663 localKeys_[n++] = LocalKey(2,1,j-1);
664 continue;
665 }
666 localKeys_[n++] = LocalKey(0,0,c++);
667 }
668 return;
669 }
670
671 if (dim==3)
672 {
673 std::array<unsigned int, dim+1> vertexMap;
674 for (unsigned int i=0; i<=dim; i++)
675 vertexMap[i] = i;
676 generateLocalKeys(vertexMap);
677 return;
678 }
679 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalCoefficients only implemented for k<=1 or dim<=3!");
680 }
681
688 constexpr LagrangeSimplexLocalCoefficients (const std::array<unsigned int, dim+1> vertexMap) requires (OrderTraits::is_static_order)
689 : localKeys_(OrderTraits::size())
690 {
691 if (dim!=2 && dim!=3)
692 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalCoefficients only implemented for dim==2 and dim==3!");
693
694 generateLocalKeys(vertexMap);
695 }
696
697
698 template<std::input_iterator VertexMap>
699 constexpr LagrangeSimplexLocalCoefficients(const VertexMap &vertexmap) requires (OrderTraits::is_static_order)
700 : localKeys_(OrderTraits::size())
701 {
702 if (dim!=2 && dim!=3)
703 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalCoefficients only implemented for dim==2 and dim==3!");
704
705 std::array<unsigned int, dim+1> vertexmap_array;
706 std::copy(vertexmap, vertexmap + dim + 1, vertexmap_array.begin());
707 generateLocalKeys(vertexmap_array);
708 }
709
711 constexpr const LocalKey& localKey (std::size_t i) const
712 {
713 return localKeys_[i];
714 }
715
716 using OrderTraits::order;
717 using OrderTraits::size;
718
719 private:
720 std::vector<LocalKey> localKeys_;
721
722 constexpr void generateLocalKeys(const std::array<unsigned int, dim+1> vertexMap)
723 {
724 const int k = OrderTraits::order();
725 if (k==0)
726 {
727 localKeys_[0] = LocalKey(0,0,0);
728 return;
729 }
730
731 if (dim==2)
732 {
733 // Create default assignment
734 int n=0;
735 int c=0;
736 for (int j=0; j<=k; j++)
737 for (int i=0; i<=k-j; i++)
738 {
739 if (i==0 && j==0)
740 {
741 localKeys_[n++] = LocalKey(0,2,0);
742 continue;
743 }
744 if (i==k && j==0)
745 {
746 localKeys_[n++] = LocalKey(1,2,0);
747 continue;
748 }
749 if (i==0 && j==k)
750 {
751 localKeys_[n++] = LocalKey(2,2,0);
752 continue;
753 }
754 if (j==0)
755 {
756 localKeys_[n++] = LocalKey(0,1,i-1);
757 continue;
758 }
759 if (i==0)
760 {
761 localKeys_[n++] = LocalKey(1,1,j-1);
762 continue;
763 }
764 if (i+j==k)
765 {
766 localKeys_[n++] = LocalKey(2,1,j-1);
767 continue;
768 }
769 localKeys_[n++] = LocalKey(0,0,c++);
770 }
771
772 // Flip edge orientations, if requested
773 bool flip[3];
774 flip[0] = vertexMap[0] > vertexMap[1];
775 flip[1] = vertexMap[0] > vertexMap[2];
776 flip[2] = vertexMap[1] > vertexMap[2];
777 for (std::size_t i=0; i<OrderTraits::size(); i++)
778 if (localKeys_[i].codim()==1 && flip[localKeys_[i].subEntity()])
779 localKeys_[i].index(k-2-localKeys_[i].index());
780
781 return;
782 }
783
784 if (dim!=3)
785 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalCoefficients only implemented for dim==3!");
786
787 unsigned int subindex[16];
788 unsigned int codim_count[4] = {0};
789 for (unsigned int m = 1; m < 16; ++m)
790 {
791 unsigned int codim = !(m&1) + !(m&2) + !(m&4) + !(m&8);
792 subindex[m] = codim_count[codim]++;
793 }
794
795 int a1 = (3*k + 12)*k + 11;
796 int a2 = -3*k - 6;
797 unsigned int dof_count[16] = {0};
798 int i[4];
799 for (i[3] = 0; i[3] <= k; ++i[3])
800 for (i[2] = 0; i[2] <= k - i[3]; ++i[2])
801 for (i[1] = 0; i[1] <= k - i[2] - i[3]; ++i[1])
802 {
803 i[0] = k - i[1] - i[2] - i[3];
804 int j[4];
805 unsigned int entity = 0;
806 unsigned int codim = 0;
807 for (unsigned int m = 0; m < 4; ++m)
808 {
809 j[m] = i[vertexMap[m]];
810 entity += !!j[m] << m;
811 codim += !j[m];
812 }
813 int local_index = j[3]*(a1 + (a2 + j[3])*j[3])/6
814 + j[2]*(2*(k - j[3]) + 3 - j[2])/2 + j[1];
815 localKeys_[local_index] = LocalKey(subindex[entity], codim, dof_count[entity]++);
816 }
817 }
818 };
819
827 template<class D, class R, int dim, int compileTimeOrder>
828 class LagrangeSimplexLocalInterpolation
829 : private LagrangeSimplexTraits<dim, compileTimeOrder>
830 {
831 using OrderTraits = LagrangeSimplexTraits<dim, compileTimeOrder>;
832 using Traits = typename LagrangeSimplexLocalBasis<D, R, dim, compileTimeOrder>::Traits;
833
834 using OrderTraits::size;
835
836 public:
837 constexpr LagrangeSimplexLocalInterpolation () requires (OrderTraits::is_static_order)
838 : OrderTraits()
839 {}
840
842 explicit constexpr LagrangeSimplexLocalInterpolation (OrderTraits orderTraits)
843 : OrderTraits(orderTraits)
844 {}
845
853 template<typename F, typename C>
854 constexpr void interpolate (const F& f, std::vector<C>& out) const
855 {
856 const int k = OrderTraits::order();
857 out.resize(OrderTraits::size());
858
859 // Specialization for zero-order case
860 if (k==0)
861 {
862 auto center = ReferenceElements<D,dim>::simplex().position(0,0);
863 out[0] = f(center);
864 return;
865 }
866
867 typename Traits::DomainType x;
868
869 // Specialization for first-order case
870 if (k==1)
871 {
872 // vertex 0
873 std::fill(x.begin(), x.end(), 0);
874 out[0] = f(x);
875
876 // remaining vertices
877 for (unsigned int i=0; i<dim; i++)
878 {
879 for (unsigned int j=0; j<dim; j++)
880 x[j] = (i==j);
881
882 out[i+1] = f(x);
883 }
884 return;
885 }
886
887 if (dim==1)
888 {
889 for (int i=0; i<k+1; i++)
890 {
891 x[0] = ((D)i)/k;
892 out[i] = f(x);
893 }
894 return;
895 }
896
897 if (dim==2)
898 {
899 int n=0;
900 for (int j=0; j<=k; j++)
901 for (int i=0; i<=k-j; i++)
902 {
903 x = { ((D)i)/k, ((D)j)/k };
904 out[n] = f(x);
905 n++;
906 }
907 return;
908 }
909
910 if (dim!=3)
911 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalInterpolation only implemented for dim<=3!");
912
913 const int kdiv = (k==0 ? 1 : k);
914
915 int n=0;
916 for (int i2 = 0; i2 <= k; i2++)
917 for (int i1 = 0; i1 <= k-i2; i1++)
918 for (int i0 = 0; i0 <= k-i1-i2; i0++)
919 {
920 x[0] = ((D)i0)/((D)kdiv);
921 x[1] = ((D)i1)/((D)kdiv);
922 x[2] = ((D)i2)/((D)kdiv);
923 out[n] = f(x);
924 n++;
925 }
926 }
927
928 };
929
930} } // namespace Dune::Impl
931
932namespace Dune
933{
987 template<class D, class R, int dim, int compileTimeOrder = -1>
989 : private Impl::LagrangeSimplexTraits<dim, compileTimeOrder>
990 {
991 using OrderTraits = Impl::LagrangeSimplexTraits<dim, compileTimeOrder>;
992
993 public:
997 Impl::LagrangeSimplexLocalBasis<D,R,dim,compileTimeOrder>,
998 Impl::LagrangeSimplexLocalCoefficients<dim,compileTimeOrder>,
999 Impl::LagrangeSimplexLocalInterpolation<D,R,dim,compileTimeOrder> >;
1000
1003 : OrderTraits()
1004 , basis_(*this)
1005 , coefficients_(*this)
1006 , interpolation_(*this)
1007 {
1008 static_assert(OrderTraits::is_static_order, "Default constructor only allowed for compile-time order >= 0");
1009 static_assert(compileTimeOrder >= 0, "Default constructor only allowed for compile-time order >= 0");
1010 }
1011
1013 explicit constexpr LagrangeSimplexLocalFiniteElement(int runTimeOrder)
1014 : OrderTraits(runTimeOrder)
1015 , basis_(*this)
1016 , coefficients_(*this)
1017 , interpolation_(*this)
1018 {
1019 if (runTimeOrder < 0)
1020 DUNE_THROW(Dune::InvalidStateException, "LagrangeSimplexLocalFiniteElement: run-time order must be non-negative!");
1021 if (compileTimeOrder >= 0 && compileTimeOrder != runTimeOrder)
1022 DUNE_THROW(Dune::InvalidStateException, "LagrangeSimplexLocalFiniteElement: Compile-time order must be identical to run-time order!");
1023 }
1024
1027 template<typename VertexMap>
1028 explicit constexpr LagrangeSimplexLocalFiniteElement(const VertexMap& vertexmap) requires(compileTimeOrder >= 0)
1029 : OrderTraits()
1030 , basis_(*this)
1031 , coefficients_(vertexmap)
1032 , interpolation_(*this)
1033 {}
1034
1037 constexpr const typename Traits::LocalBasisType& localBasis () const
1038 {
1039 return basis_;
1040 }
1041
1044 constexpr const typename Traits::LocalCoefficientsType& localCoefficients () const
1045 {
1046 return coefficients_;
1047 }
1048
1051 constexpr const typename Traits::LocalInterpolationType& localInterpolation () const
1052 {
1053 return interpolation_;
1054 }
1055
1057 using OrderTraits::size;
1058
1061 static constexpr GeometryType type ()
1062 {
1063 return GeometryTypes::simplex(dim);
1064 }
1065
1066 private:
1067 typename Traits::LocalBasisType basis_;
1068 typename Traits::LocalCoefficientsType coefficients_;
1069 typename Traits::LocalInterpolationType interpolation_;
1070 };
1071
1072} // namespace Dune
1073
1074#endif // DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGESIMPLEX_HH
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:375
Lagrange finite element for simplices with arbitrary compile-time dimension and static or dynamic pol...
Definition: lagrangesimplex.hh:990
constexpr const Traits::LocalCoefficientsType & localCoefficients() const
Returns the assignment of the degrees of freedom to the element subentities.
Definition: lagrangesimplex.hh:1044
constexpr LagrangeSimplexLocalFiniteElement(const VertexMap &vertexmap)
Constructs a finite element given a vertex reordering.
Definition: lagrangesimplex.hh:1028
static constexpr GeometryType type()
The reference element that the local finite element is defined on.
Definition: lagrangesimplex.hh:1061
constexpr const Traits::LocalBasisType & localBasis() const
Returns the local basis, i.e., the set of shape functions.
Definition: lagrangesimplex.hh:1037
constexpr const Traits::LocalInterpolationType & localInterpolation() const
Returns object that evaluates degrees of freedom.
Definition: lagrangesimplex.hh:1051
constexpr LagrangeSimplexLocalFiniteElement(int runTimeOrder)
Constructor for run-time order.
Definition: lagrangesimplex.hh:1013
constexpr LagrangeSimplexLocalFiniteElement()
Constructor for compile-time order.
Definition: lagrangesimplex.hh:1002
A few common exception classes.
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.
decltype(auto) constexpr unpackIntegerSequence(F &&f, std::integer_sequence< I, i... > sequence)
Unpack an std::integer_sequence<I,i...> to std::integral_constant<I,i>...
Definition: indices.hh:124
std::integral_constant< std::size_t, i > index_constant
An index constant with value i.
Definition: indices.hh:29
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:453
void interpolate(const F &f, const GFS &gfs, XG &xg)
interpolation from a given grid function
Definition: interpolate.hh:177
constexpr decltype(auto) switchCases(const Cases &cases, const Value &value, Branches &&branches, ElseBranch &&elseBranch)
Switch statement.
Definition: hybridutilities.hh:661
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:489
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:284
static constexpr IntegralRange< std::decay_t< T > > range(T &&from, U &&to) noexcept
free standing function for setting up a range based for loop over an integer range for (auto i: range...
Definition: rangeutilities.hh:288
Some useful basic math stuff.
Dune namespace
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
static constexpr T binomial(const T &n, const T &k) noexcept
calculate the binomial coefficient n over k as a constexpr
Definition: math.hh:113
STL namespace.
Utilities for reduction like operations on ranges.
static const ReferenceElement & simplex()
get simplex reference elements
Definition: referenceelements.hh:162
D DomainType
domain type
Definition: localbasis.hh:43
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (May 19, 22:31, 2026)