Loading [MathJax]/extensions/tex2jax.js

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 <numeric>
10#include <algorithm>
11
15#include <dune/common/hybridutilities.hh>
16#include <dune/common/math.hh>
18
19#include <dune/geometry/referenceelements.hh>
20
21#include <dune/localfunctions/common/localbasis.hh>
22#include <dune/localfunctions/common/localfiniteelementtraits.hh>
23#include <dune/localfunctions/common/localkey.hh>
24
25namespace Dune { namespace Impl
26{
37 template<class D, class R, unsigned int dim, unsigned int k>
38 class LagrangeSimplexLocalBasis
39 {
40
41 // Compute the rescaled barycentric coordinates of x.
42 // We rescale the simplex by k and then compute the
43 // barycentric coordinates with respect to the points
44 // p_i = e_i (for i=0,...,dim-1) and p_dim=0.
45 // Notice that then the Lagrange points have the barycentric
46 // coordinates (i_0,...,i_d) where i_j are all non-negative
47 // integers satisfying the constraint sum i_j = k.
48 constexpr auto barycentric(const auto& x) const
49 {
50 auto b = std::array<R,dim+1>{};
51 b[dim] = k;
52 for(auto i : Dune::range(dim))
53 {
54 b[i] = k*x[i];
55 b[dim] -= b[i];
56 }
57 return b;
58 }
59
60 // Evaluate the univariate Lagrange polynomials L_i(t) for i=0,...,k where
61 //
62 // L_i(t) = (t-0)/(i-0) * ... * (t-(i-1))/(i-(i-1))
63 // = (t-0)*...*(t-(i-1))/(i!);
64 static constexpr void evaluateLagrangePolynomials(const R& t, auto& L)
65 {
66 L[0] = 1;
67 for (auto i : Dune::range(k))
68 L[i+1] = L[i] * (t - i) / (i+1);
69 }
70
71 // Evaluate the univariate Lagrange polynomial derivatives L_i(t) for i=0,...,k
72 // up to given maxDerivativeOrder.
73 static constexpr void evaluateLagrangePolynomialDerivative(const R& t, auto& LL, unsigned int maxDerivativeOrder)
74 {
75 auto& L = LL[0];
76 L[0] = 1;
77 for (auto i : Dune::range(k))
78 L[i+1] = L[i] * (t - i) / (i+1);
79 for(auto j : Dune::range(maxDerivativeOrder))
80 {
81 auto& F = LL[j];
82 auto& DF = LL[j+1];
83 DF[0] = 0;
84 for (auto i : Dune::range(k))
85 DF[i+1] = (DF[i] * (t - i) + (j+1)*F[i]) / (i+1);
86 }
87 }
88
89 using BarycentricMultiIndex = std::array<unsigned int,dim+1>;
90
91
92 // This computed the required partial derivative given by the multi-index
93 // beta of a product of a function given as a product of dim+1 derivatives
94 // of univariate Lagrange polynomials of the dim+1 barycentric coordinates.
95 // The polynomials in the product are specified as follows:
96 //
97 // The table L contains all required derivatives of all univariate
98 // polynomials evaluated at all barycentric coordinates. The two
99 // multi-indices i and alpha encode that the polynomial for the
100 // j-th barycentric coordinate is the alpha-j-th derivative of
101 // the i_j-the Lagrange polynomial.
102 //
103 // Hence this method computes D_beta f(x) where f(x) is the product
104 // \f$f(x) = \prod_{j=0}^{d} L_{i_j}^{(alpha_j)}(x_j) \f$.
105 static constexpr R barycentricDerivative(
106 BarycentricMultiIndex beta,
107 const auto&L,
108 const BarycentricMultiIndex& i,
109 const BarycentricMultiIndex& alpha = {})
110 {
111 // If there are unprocessed derivatives left we search the first unprocessed
112 // partial derivative direction j and compute it using the product and chain rule.
113 // The remaining partial derivatives are forwarded to the recursion.
114 // Notice that the partial derivative of the last barycentric component
115 // wrt to the j-th component is -1 which is responsible for the second
116 // term in the sum. Furthermore we get the factor k due to the scaling of
117 // the simplex.
118 for(auto j : Dune::range(dim))
119 {
120 if (beta[j] > 0)
121 {
122 auto leftDerivatives = alpha;
123 auto rightDerivatives = alpha;
124 leftDerivatives[j]++;
125 rightDerivatives.back()++;
126 beta[j]--;
127 return (barycentricDerivative(beta, L, i, leftDerivatives) - barycentricDerivative(beta, L, i, rightDerivatives))*k;
128 }
129 }
130 // If there are no unprocessed derivatives we can simply evaluate
131 // the product of the derivatives of the Lagrange polynomials with
132 // given indices i_j and derivative orders alpha_j
133 // Evaluate the product of the univariate Lagrange polynomials
134 // with given indices and orders.
135 auto y = R(1);
136 for(auto j : Dune::range(dim+1))
137 y *= L[j][alpha[j]][i[j]];
138 return y;
139 }
140
141
142 public:
143 using Traits = LocalBasisTraits<D,dim,FieldVector<D,dim>,R,1,FieldVector<R,1>,FieldMatrix<R,1,dim> >;
144
149 static constexpr unsigned int size ()
150 {
151 return binomial(k+dim,dim);
152 }
153
155 void evaluateFunction(const typename Traits::DomainType& x,
156 std::vector<typename Traits::RangeType>& out) const
157 {
158 out.resize(size());
159
160 // Specialization for zero-order case
161 if (k==0)
162 {
163 out[0] = 1;
164 return;
165 }
166
167 // Specialization for first-order case
168 if (k==1)
169 {
170 out[0] = 1.0;
171 for (size_t i=0; i<dim; i++)
172 {
173 out[0] -= x[i];
174 out[i+1] = x[i];
175 }
176 return;
177 }
178
179 // Compute rescaled barycentric coordinates of x
180 auto z = barycentric(x);
181
182 auto L = std::array<std::array<R,k+1>, dim+1>();
183 for (auto j : Dune::range(dim+1))
184 evaluateLagrangePolynomials(z[j], L[j]);
185
186 if (dim==1)
187 {
188 unsigned int n = 0;
189 for (auto i0 : Dune::range(k + 1))
190 for (auto i1 : std::array{k - i0})
191 out[n++] = L[0][i0] * L[1][i1];
192 return;
193 }
194 if (dim==2)
195 {
196 unsigned int n=0;
197 for (auto i1 : Dune::range(k + 1))
198 for (auto i0 : Dune::range(k - i1 + 1))
199 for (auto i2 : std::array{k - i1 - i0})
200 out[n++] = L[0][i0] * L[1][i1] * L[2][i2];
201 return;
202 }
203 if (dim==3)
204 {
205 unsigned int n = 0;
206 for (auto i2 : Dune::range(k + 1))
207 for (auto i1 : Dune::range(k - i2 + 1))
208 for (auto i0 : Dune::range(k - i2 - i1 + 1))
209 for (auto i3 : std::array{k - i2 - i1 - i0})
210 out[n++] = L[0][i0] * L[1][i1] * L[2][i2] * L[3][i3];
211 return;
212 }
213
214 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalBasis for k>=2 only implemented for dim<=3");
215 }
216
222 void evaluateJacobian(const typename Traits::DomainType& x,
223 std::vector<typename Traits::JacobianType>& out) const
224 {
225 out.resize(size());
226
227 // Specialization for k==0
228 if (k==0)
229 {
230 std::fill(out[0][0].begin(), out[0][0].end(), 0);
231 return;
232 }
233
234 // Specialization for k==1
235 if (k==1)
236 {
237 std::fill(out[0][0].begin(), out[0][0].end(), -1);
238
239 for (unsigned int i=0; i<dim; i++)
240 for (unsigned int j=0; j<dim; j++)
241 out[i+1][0][j] = (i==j);
242
243 return;
244 }
245
246 // Compute rescaled barycentric coordinates of x
247 auto z = barycentric(x);
248
249 // L[j][m][i] is the m-th derivative of the i-th Lagrange polynomial at z[j]
250 auto L = std::array<std::array<std::array<R,k+1>, 2>, dim+1>();
251 for (auto j : Dune::range(dim+1))
252 evaluateLagrangePolynomialDerivative(z[j], L[j], 1);
253
254 if (dim==1)
255 {
256 unsigned int n = 0;
257 for (auto i0 : Dune::range(k + 1))
258 {
259 for (auto i1 : std::array{k-i0})
260 {
261 out[n][0][0] = (L[0][1][i0] * L[1][0][i1] - L[0][0][i0] * L[1][1][i1])*k;
262 ++n;
263 }
264 }
265 return;
266 }
267 if (dim==2)
268 {
269 unsigned int n=0;
270 for (auto i1 : Dune::range(k + 1))
271 {
272 for (auto i0 : Dune::range(k - i1 + 1))
273 {
274 for (auto i2 : std::array{k - i1 - i0})
275 {
276 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;
277 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;
278 ++n;
279 }
280 }
281 }
282 return;
283 }
284 if (dim==3)
285 {
286 unsigned int n = 0;
287 for (auto i2 : Dune::range(k + 1))
288 {
289 for (auto i1 : Dune::range(k - i2 + 1))
290 {
291 for (auto i0 : Dune::range(k - i2 - i1 + 1))
292 {
293 for (auto i3 : std::array{k - i2 - i1 - i0})
294 {
295 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;
296 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;
297 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;
298 ++n;
299 }
300 }
301 }
302 }
303
304 return;
305 }
306
307 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalBasis for k>=2 only implemented for dim<=3");
308 }
309
316 void partial(const std::array<unsigned int,dim>& order,
317 const typename Traits::DomainType& in,
318 std::vector<typename Traits::RangeType>& out) const
319 {
320 auto totalOrder = std::accumulate(order.begin(), order.end(), 0u);
321
322 out.resize(size());
323
324 // Derivative order zero corresponds to the function evaluation.
325 if (totalOrder == 0)
326 {
327 evaluateFunction(in, out);
328 return;
329 }
330
331 // Derivatives of order >k are all zero.
332 if (totalOrder > k)
333 {
334 for(auto& out_i : out)
335 out_i = 0;
336 return;
337 }
338
339 // It remains to cover the cases 0 < totalOrder<= k.
340
341 if (k==1)
342 {
343 if (totalOrder==1)
344 {
345 auto direction = std::find(order.begin(), order.end(), 1);
346 out[0] = -1;
347 for (unsigned int i=0; i<dim; i++)
348 out[i+1] = (i==(direction-order.begin()));
349 }
350 return;
351 }
352
353 // Since the required stack storage depends on the dynamic total order,
354 // we need to do a dynamic to static dispatch by enumerating all supported
355 // static orders.
356 auto supportedStaticOrders = Dune::range(Dune::index_constant<1>{}, Dune::index_constant<k+1>{});
357 return Dune::Hybrid::switchCases(supportedStaticOrders, totalOrder, [&](auto staticTotalOrder) {
358
359 // Compute rescaled barycentric coordinates of x
360 auto z = barycentric(in);
361
362 // L[j][m][i] is the m-th derivative of the i-th Lagrange polynomial at z[j]
363 auto L = std::array<std::array<std::array<R, k+1>, staticTotalOrder+1>, dim+1>();
364 for (auto j : Dune::range(dim))
365 evaluateLagrangePolynomialDerivative(z[j], L[j], order[j]);
366 evaluateLagrangePolynomialDerivative(z[dim], L[dim], totalOrder);
367
368 auto barycentricOrder = BarycentricMultiIndex{};
369 for (auto j : Dune::range(dim))
370 barycentricOrder[j] = order[j];
371 barycentricOrder[dim] = 0;
372
373 if constexpr (dim==1)
374 {
375 unsigned int n = 0;
376 for (auto i0 : Dune::range(k + 1))
377 for (auto i1 : std::array{k - i0})
378 out[n++] = barycentricDerivative(barycentricOrder, L, BarycentricMultiIndex{i0, i1});
379 }
380 if constexpr (dim==2)
381 {
382 unsigned int n=0;
383 for (auto i1 : Dune::range(k + 1))
384 for (auto i0 : Dune::range(k - i1 + 1))
385 for (auto i2 : std::array{k - i1 - i0})
386 out[n++] = barycentricDerivative(barycentricOrder, L, BarycentricMultiIndex{i0, i1, i2});
387 }
388 if constexpr (dim==3)
389 {
390 unsigned int n = 0;
391 for (auto i2 : Dune::range(k + 1))
392 for (auto i1 : Dune::range(k - i2 + 1))
393 for (auto i0 : Dune::range(k - i2 - i1 + 1))
394 for (auto i3 : std::array{k - i2 - i1 - i0})
395 out[n++] = barycentricDerivative(barycentricOrder, L, BarycentricMultiIndex{i0, i1, i2, i3});
396 }
397 });
398 }
399
401 static constexpr unsigned int order ()
402 {
403 return k;
404 }
405 };
406
412 template<unsigned int dim, unsigned int k>
413 class LagrangeSimplexLocalCoefficients
414 {
415 public:
417 LagrangeSimplexLocalCoefficients ()
418 : localKeys_(size())
419 {
420 if (k==0)
421 {
422 localKeys_[0] = LocalKey(0,0,0);
423 return;
424 }
425
426 if (k==1)
427 {
428 for (std::size_t i=0; i<size(); i++)
429 localKeys_[i] = LocalKey(i,dim,0);
430 return;
431 }
432
433 if (dim==1)
434 {
435 // Order is at least 2 here
436 localKeys_[0] = LocalKey(0,1,0); // vertex dof
437 for (unsigned int i=1; i<k; i++)
438 localKeys_[i] = LocalKey(0,0,i-1); // element dofs
439 localKeys_.back() = LocalKey(1,1,0); // vertex dof
440 return;
441 }
442
443 if (dim==2)
444 {
445 int n=0;
446 int c=0;
447 for (unsigned int j=0; j<=k; j++)
448 for (unsigned int i=0; i<=k-j; i++)
449 {
450 if (i==0 && j==0)
451 {
452 localKeys_[n++] = LocalKey(0,2,0);
453 continue;
454 }
455 if (i==k && j==0)
456 {
457 localKeys_[n++] = LocalKey(1,2,0);
458 continue;
459 }
460 if (i==0 && j==k)
461 {
462 localKeys_[n++] = LocalKey(2,2,0);
463 continue;
464 }
465 if (j==0)
466 {
467 localKeys_[n++] = LocalKey(0,1,i-1);
468 continue;
469 }
470 if (i==0)
471 {
472 localKeys_[n++] = LocalKey(1,1,j-1);
473 continue;
474 }
475 if (i+j==k)
476 {
477 localKeys_[n++] = LocalKey(2,1,j-1);
478 continue;
479 }
480 localKeys_[n++] = LocalKey(0,0,c++);
481 }
482 return;
483 }
484
485 if (dim==3)
486 {
487 std::array<unsigned int, dim+1> vertexMap;
488 for (unsigned int i=0; i<=dim; i++)
489 vertexMap[i] = i;
490 generateLocalKeys(vertexMap);
491 return;
492 }
493 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalCoefficients only implemented for k<=1 or dim<=3!");
494 }
495
502 LagrangeSimplexLocalCoefficients (const std::array<unsigned int, dim+1> vertexMap)
503 : localKeys_(size())
504 {
505 if (dim!=2 && dim!=3)
506 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalCoefficients only implemented for dim==2 and dim==3!");
507
508 generateLocalKeys(vertexMap);
509 }
510
511
512 template<class VertexMap>
513 LagrangeSimplexLocalCoefficients(const VertexMap &vertexmap)
514 : localKeys_(size())
515 {
516 if (dim!=2 && dim!=3)
517 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalCoefficients only implemented for dim==2 and dim==3!");
518
519 std::array<unsigned int, dim+1> vertexmap_array;
520 std::copy(vertexmap, vertexmap + dim + 1, vertexmap_array.begin());
521 generateLocalKeys(vertexmap_array);
522 }
523
525 static constexpr std::size_t size ()
526 {
527 return binomial(k+dim,dim);
528 }
529
531 const LocalKey& localKey (std::size_t i) const
532 {
533 return localKeys_[i];
534 }
535
536 private:
537 std::vector<LocalKey> localKeys_;
538
539 void generateLocalKeys(const std::array<unsigned int, dim+1> vertexMap)
540 {
541 if (k==0)
542 {
543 localKeys_[0] = LocalKey(0,0,0);
544 return;
545 }
546
547 if (dim==2)
548 {
549 // Create default assignment
550 int n=0;
551 int c=0;
552 for (unsigned int j=0; j<=k; j++)
553 for (unsigned int i=0; i<=k-j; i++)
554 {
555 if (i==0 && j==0)
556 {
557 localKeys_[n++] = LocalKey(0,2,0);
558 continue;
559 }
560 if (i==k && j==0)
561 {
562 localKeys_[n++] = LocalKey(1,2,0);
563 continue;
564 }
565 if (i==0 && j==k)
566 {
567 localKeys_[n++] = LocalKey(2,2,0);
568 continue;
569 }
570 if (j==0)
571 {
572 localKeys_[n++] = LocalKey(0,1,i-1);
573 continue;
574 }
575 if (i==0)
576 {
577 localKeys_[n++] = LocalKey(1,1,j-1);
578 continue;
579 }
580 if (i+j==k)
581 {
582 localKeys_[n++] = LocalKey(2,1,j-1);
583 continue;
584 }
585 localKeys_[n++] = LocalKey(0,0,c++);
586 }
587
588 // Flip edge orientations, if requested
589 bool flip[3];
590 flip[0] = vertexMap[0] > vertexMap[1];
591 flip[1] = vertexMap[0] > vertexMap[2];
592 flip[2] = vertexMap[1] > vertexMap[2];
593 for (std::size_t i=0; i<size(); i++)
594 if (localKeys_[i].codim()==1 && flip[localKeys_[i].subEntity()])
595 localKeys_[i].index(k-2-localKeys_[i].index());
596
597 return;
598 }
599
600 if (dim!=3)
601 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalCoefficients only implemented for dim==3!");
602
603 unsigned int subindex[16];
604 unsigned int codim_count[4] = {0};
605 for (unsigned int m = 1; m < 16; ++m)
606 {
607 unsigned int codim = !(m&1) + !(m&2) + !(m&4) + !(m&8);
608 subindex[m] = codim_count[codim]++;
609 }
610
611 int a1 = (3*k + 12)*k + 11;
612 int a2 = -3*k - 6;
613 unsigned int dof_count[16] = {0};
614 unsigned int i[4];
615 for (i[3] = 0; i[3] <= k; ++i[3])
616 for (i[2] = 0; i[2] <= k - i[3]; ++i[2])
617 for (i[1] = 0; i[1] <= k - i[2] - i[3]; ++i[1])
618 {
619 i[0] = k - i[1] - i[2] - i[3];
620 unsigned int j[4];
621 unsigned int entity = 0;
622 unsigned int codim = 0;
623 for (unsigned int m = 0; m < 4; ++m)
624 {
625 j[m] = i[vertexMap[m]];
626 entity += !!j[m] << m;
627 codim += !j[m];
628 }
629 int local_index = j[3]*(a1 + (a2 + j[3])*j[3])/6
630 + j[2]*(2*(k - j[3]) + 3 - j[2])/2 + j[1];
631 localKeys_[local_index] = LocalKey(subindex[entity], codim, dof_count[entity]++);
632 }
633 }
634 };
635
640 template<class LocalBasis>
641 class LagrangeSimplexLocalInterpolation
642 {
643 static const int kdiv = (LocalBasis::order() == 0 ? 1 : LocalBasis::order());
644 public:
645
653 template<typename F, typename C>
654 void interpolate (const F& f, std::vector<C>& out) const
655 {
656 constexpr auto dim = LocalBasis::Traits::dimDomain;
657 constexpr auto k = LocalBasis::order();
658 using D = typename LocalBasis::Traits::DomainFieldType;
659
660 typename LocalBasis::Traits::DomainType x;
661
662 out.resize(LocalBasis::size());
663
664 // Specialization for zero-order case
665 if (k==0)
666 {
667 auto center = ReferenceElements<D,dim>::simplex().position(0,0);
668 out[0] = f(center);
669 return;
670 }
671
672 // Specialization for first-order case
673 if (k==1)
674 {
675 // vertex 0
676 std::fill(x.begin(), x.end(), 0);
677 out[0] = f(x);
678
679 // remaining vertices
680 for (int i=0; i<dim; i++)
681 {
682 for (int j=0; j<dim; j++)
683 x[j] = (i==j);
684
685 out[i+1] = f(x);
686 }
687 return;
688 }
689
690 if (dim==1)
691 {
692 for (unsigned int i=0; i<k+1; i++)
693 {
694 x[0] = ((D)i)/k;
695 out[i] = f(x);
696 }
697 return;
698 }
699
700 if (dim==2)
701 {
702 int n=0;
703 for (unsigned int j=0; j<=k; j++)
704 for (unsigned int i=0; i<=k-j; i++)
705 {
706 x = { ((D)i)/k, ((D)j)/k };
707 out[n] = f(x);
708 n++;
709 }
710 return;
711 }
712
713 if (dim!=3)
714 DUNE_THROW(NotImplemented, "LagrangeSimplexLocalInterpolation only implemented for dim<=3!");
715
716 int n=0;
717 for (int i2 = 0; i2 <= (int)k; i2++)
718 for (int i1 = 0; i1 <= (int)k-i2; i1++)
719 for (int i0 = 0; i0 <= (int)k-i1-i2; i0++)
720 {
721 x[0] = ((D)i0)/((D)kdiv);
722 x[1] = ((D)i1)/((D)kdiv);
723 x[2] = ((D)i2)/((D)kdiv);
724 out[n] = f(x);
725 n++;
726 }
727 }
728
729 };
730
731} } // namespace Dune::Impl
732
733namespace Dune
734{
788 template<class D, class R, int d, int k>
790 {
791 public:
795 Impl::LagrangeSimplexLocalCoefficients<d,k>,
796 Impl::LagrangeSimplexLocalInterpolation<Impl::LagrangeSimplexLocalBasis<D,R,d,k> > >;
797
800
803 template<typename VertexMap>
804 LagrangeSimplexLocalFiniteElement(const VertexMap& vertexmap)
805 : coefficients_(vertexmap)
806 {}
807
810 const typename Traits::LocalBasisType& localBasis () const
811 {
812 return basis_;
813 }
814
818 {
819 return coefficients_;
820 }
821
825 {
826 return interpolation_;
827 }
828
830 static constexpr std::size_t size ()
831 {
832 return Impl::LagrangeSimplexLocalBasis<D,R,d,k>::size();
833 }
834
837 static constexpr GeometryType type ()
838 {
839 return GeometryTypes::simplex(d);
840 }
841
842 private:
843 Impl::LagrangeSimplexLocalBasis<D,R,d,k> basis_;
844 Impl::LagrangeSimplexLocalCoefficients<d,k> coefficients_;
845 Impl::LagrangeSimplexLocalInterpolation<Impl::LagrangeSimplexLocalBasis<D,R,d,k> > interpolation_;
846 };
847
848} // namespace Dune
849
850#endif // DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGESIMPLEX_HH
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
Lagrange finite element for simplices with arbitrary compile-time dimension and polynomial order.
Definition: lagrangesimplex.hh:790
const Traits::LocalInterpolationType & localInterpolation() const
Returns object that evaluates degrees of freedom.
Definition: lagrangesimplex.hh:824
const Traits::LocalBasisType & localBasis() const
Returns the local basis, i.e., the set of shape functions.
Definition: lagrangesimplex.hh:810
LagrangeSimplexLocalFiniteElement()
Definition: lagrangesimplex.hh:799
static constexpr std::size_t size()
The number of shape functions.
Definition: lagrangesimplex.hh:830
LagrangeSimplexLocalFiniteElement(const VertexMap &vertexmap)
Constructs a finite element given a vertex reordering.
Definition: lagrangesimplex.hh:804
const Traits::LocalCoefficientsType & localCoefficients() const
Returns the assignment of the degrees of freedom to the element subentities.
Definition: lagrangesimplex.hh:817
static constexpr GeometryType type()
The reference element that the local finite element is defined on.
Definition: lagrangesimplex.hh:837
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.
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:674
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:280
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:294
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:131
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 (Apr 5, 23:02, 2025)