DUNE PDELab (unstable)

lagrangecube.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_LAGRANGECUBE_HH
6#define DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGECUBE_HH
7
8#include <array>
9#include <numeric>
10
13#include <dune/common/math.hh>
15
16#include <dune/geometry/referenceelements.hh>
17
18#include <dune/localfunctions/common/localbasis.hh>
19#include <dune/localfunctions/common/localfiniteelementtraits.hh>
20#include <dune/localfunctions/common/localkey.hh>
21
22namespace Dune { namespace Impl
23{
24 // The traits provide static or dynamic order and size information
25 template<unsigned int dim, int compileTimeOrder>
26 struct LagrangeCubeOrderTraits;
27
28 // The traits provide static order and size information
29 template<unsigned int dim, int compileTimeOrder>
30 requires (compileTimeOrder >= 0)
31 struct LagrangeCubeOrderTraits<dim,compileTimeOrder>
32 {
33 static constexpr bool is_static_order = true;
34
35 constexpr LagrangeCubeOrderTraits(int /*runTimeOrder*/ = compileTimeOrder) {}
36
40 static constexpr std::size_t size ()
41 {
42 return power(compileTimeOrder+1, dim);
43 }
44
48 static constexpr unsigned int order ()
49 {
50 return compileTimeOrder;
51 }
52
53 // Return i as a d-digit number in the (k+1)-nary system
54 static constexpr std::array<unsigned int,dim> multiindex (unsigned int i)
55 {
56 std::array<unsigned int,dim> alpha;
57 for (unsigned int j=0; j<dim; j++)
58 {
59 alpha[j] = i % (order()+1);
60 i = i/(order()+1);
61 }
62 return alpha;
63 }
64 };
65
66
67 template<unsigned int dim>
68 struct LagrangeCubeOrderTraits<dim,-1>
69 {
70 unsigned int order_;
71 std::size_t size_;
72
73 static constexpr bool is_static_order = false;
74
78 constexpr LagrangeCubeOrderTraits (int runTimeOrder)
79 : order_(runTimeOrder >= 0 ? (unsigned int)(runTimeOrder) : 0u)
80 , size_(power(order_+1, dim))
81 {
82 if (runTimeOrder < 0)
83 DUNE_THROW(Dune::InvalidStateException, "LagrangeCube: run-time order must be >= 0");
84 }
85
89 constexpr std::size_t size () const
90 {
91 return size_;
92 }
93
97 constexpr unsigned int order () const
98 {
99 return order_;
100 }
101
102 // Return i as a d-digit number in the (k+1)-nary system
103 constexpr std::array<unsigned int,dim> multiindex (unsigned int i) const
104 {
105 std::array<unsigned int,dim> alpha;
106 for (unsigned int j=0; j<dim; j++)
107 {
108 alpha[j] = i % (order()+1);
109 i = i/(order()+1);
110 }
111 return alpha;
112 }
113 };
114
125 template<class D, class R, unsigned int dim, int compileTimeOrder>
126 class LagrangeCubeLocalBasis
127 : private LagrangeCubeOrderTraits<dim,compileTimeOrder>
128 {
129 using OrderTraits = LagrangeCubeOrderTraits<dim,compileTimeOrder>;
130
131 // i-th Lagrange polynomial of degree k in one dimension
132 constexpr R p (unsigned int i, D x) const
133 {
134 R result(1.0);
135 for (unsigned int j=0; j<=order(); j++)
136 if (j!=i) result *= (order()*x-j)/((int)i-(int)j);
137 return result;
138 }
139
140 // derivative of ith Lagrange polynomial of degree k in one dimension
141 constexpr R dp (unsigned int i, D x) const
142 {
143 R result(0.0);
144
145 for (unsigned int j=0; j<=order(); j++)
146 {
147 if (j!=i)
148 {
149 R prod( (order()*1.0)/((int)i-(int)j) );
150 for (unsigned int l=0; l<=order(); l++)
151 if (l!=i && l!=j)
152 prod *= (order()*x-l)/((int)i-(int)l);
153 result += prod;
154 }
155 }
156 return result;
157 }
158
159 // Second derivative of j-th Lagrange polynomial of degree k in one dimension
160 // Formula and notation taken from https://en.wikipedia.org/wiki/Lagrange_polynomial#Derivatives
161 constexpr R ddp(unsigned int j, D x) const
162 {
163 R result(0.0);
164
165 for (unsigned int i=0; i<=order(); i++)
166 {
167 if (i==j)
168 continue;
169
170 R sum(0);
171
172 for (unsigned int m=0; m<=order(); m++)
173 {
174 if (m==i || m==j)
175 continue;
176
177 R prod( (order()*1.0)/((int)j-(int)m) );
178 for (unsigned int l=0; l<=order(); l++)
179 if (l!=i && l!=j && l!=m)
180 prod *= (order()*x-l)/((int)j-(int)l);
181 sum += prod;
182 }
183
184 result += sum * ( (order()*1.0)/((int)j-(int)i) );
185 }
186
187 return result;
188 }
189
190 using OrderTraits::multiindex;
191
192 public:
193 using Traits = LocalBasisTraits<D,dim,FieldVector<D,dim>,R,1,FieldVector<R,1>,FieldMatrix<R,1,dim> >;
194
195 constexpr LagrangeCubeLocalBasis () requires (OrderTraits::is_static_order)
196 : OrderTraits()
197 {}
198
201 explicit constexpr LagrangeCubeLocalBasis (OrderTraits orderTraits)
202 : OrderTraits(orderTraits)
203 {}
204
205 using OrderTraits::size;
206 using OrderTraits::order;
207
209 constexpr void evaluateFunction (const typename Traits::DomainType& x,
210 std::vector<typename Traits::RangeType>& out) const
211 {
212 out.resize(size());
213
214 // Specialization for zero-order case
215 if (order()==0)
216 {
217 out[0] = 1;
218 return;
219 }
220
221 if (order()==1)
222 {
223 for (size_t i=0; i<size(); i++)
224 {
225 out[i] = 1;
226
227 for (unsigned int j=0; j<dim; j++)
228 // if j-th bit of i is set multiply with x[j], else with 1-x[j]
229 out[i] *= (i & (1<<j)) ? x[j] : 1-x[j];
230 }
231 return;
232 }
233
234 // General case
235 for (size_t i=0; i<size(); i++)
236 {
237 // convert index i to multiindex
238 std::array<unsigned int,dim> alpha(multiindex(i));
239
240 // initialize product
241 out[i] = 1.0;
242
243 // dimension by dimension
244 for (unsigned int j=0; j<dim; j++)
245 out[i] *= p(alpha[j],x[j]);
246 }
247 }
248
254 constexpr void evaluateJacobian (const typename Traits::DomainType& x,
255 std::vector<typename Traits::JacobianType>& out) const
256 {
257 out.resize(size());
258
259 // Specialization for common case
260 if (order()==0)
261 {
262 std::fill(out[0][0].begin(), out[0][0].end(), 0);
263 return;
264 }
265
266 // Specialization for common case
267 if (order()==1)
268 {
269 // Loop over all shape functions
270 for (size_t i=0; i<size(); i++)
271 {
272 // Loop over all coordinate directions
273 for (unsigned int j=0; j<dim; j++)
274 {
275 // Initialize: the overall expression is a product
276 // if j-th bit of i is set to 1, else -11
277 out[i][0][j] = (i & (1<<j)) ? 1 : -1;
278
279 for (unsigned int l=0; l<dim; l++)
280 {
281 if (j!=l)
282 // if l-th bit of i is set multiply with x[l], else with 1-x[l]
283 out[i][0][j] *= (i & (1<<l)) ? x[l] : 1-x[l];
284 }
285 }
286 }
287 return;
288 }
289
290 // The general case
291
292 // Loop over all shape functions
293 for (size_t i=0; i<size(); i++)
294 {
295 // convert index i to multiindex
296 std::array<unsigned int,dim> alpha(multiindex(i));
297
298 // Loop over all coordinate directions
299 for (unsigned int j=0; j<dim; j++)
300 {
301 // Initialize: the overall expression is a product
302 // if j-th bit of i is set to -1, else 1
303 out[i][0][j] = dp(alpha[j],x[j]);
304
305 // rest of the product
306 for (unsigned int l=0; l<dim; l++)
307 if (l!=j)
308 out[i][0][j] *= p(alpha[l],x[l]);
309 }
310 }
311 }
312
319 constexpr void partial (const std::array<unsigned int,dim>& partialOrders,
320 const typename Traits::DomainType& in,
321 std::vector<typename Traits::RangeType>& out) const
322 {
323 auto totalOrder = std::accumulate(partialOrders.begin(), partialOrders.end(), 0);
324
325 out.resize(size());
326
327 if (order()==0)
328 {
329 out[0] = (totalOrder==0);
330 return;
331 }
332
333 if (order()==1)
334 {
335 if (totalOrder == 0)
336 {
337 evaluateFunction(in, out);
338 }
339 else if (totalOrder == 1)
340 {
341 out.resize(size());
342
343 auto direction = std::distance(partialOrders.begin(),
344 std::find(partialOrders.begin(), partialOrders.end(), 1));
345 if (direction >= dim)
346 DUNE_THROW(RangeError, "Direction of partial derivative not found!");
347
348 // Loop over all shape functions
349 for (std::size_t i = 0; i < size(); ++i)
350 {
351 // Initialize: the overall expression is a product
352 // if j-th bit of i is set to 1, otherwise to -1
353 out[i] = (i & (1<<direction)) ? 1 : -1;
354
355 for (unsigned int j = 0; j < dim; ++j)
356 {
357 if (direction != j)
358 // if j-th bit of i is set multiply with in[j], else with 1-in[j]
359 out[i] *= (i & (1<<j)) ? in[j] : 1-in[j];
360 }
361 }
362 }
363 else if (totalOrder == 2)
364 {
365
366 for (size_t i=0; i<size(); i++)
367 {
368 // convert index i to multiindex
369 std::array<unsigned int,dim> alpha(multiindex(i));
370
371 // Initialize: the overall expression is a product
372 out[i][0] = 1.0;
373
374 // rest of the product
375 for (std::size_t l=0; l<dim; l++)
376 {
377 switch (partialOrders[l])
378 {
379 case 0:
380 out[i][0] *= p(alpha[l],in[l]);
381 break;
382 case 1:
383 //std::cout << "dp: " << dp(alpha[l],in[l]) << std::endl;
384 out[i][0] *= dp(alpha[l],in[l]);
385 break;
386 case 2:
387 out[i][0] *= 0;
388 break;
389 default:
390 DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
391 }
392 }
393 }
394 }
395 else
396 DUNE_THROW(NotImplemented, "Partial derivative of order " << totalOrder << " is not implemented!");
397
398 return;
399 }
400
401 // The case k>1
402
403 // Loop over all shape functions
404 for (size_t i=0; i<size(); i++)
405 {
406 // convert index i to multiindex
407 std::array<unsigned int,dim> alpha(multiindex(i));
408
409 // Initialize: the overall expression is a product
410 out[i][0] = 1.0;
411
412 // rest of the product
413 for (std::size_t l=0; l<dim; l++)
414 {
415 switch (partialOrders[l])
416 {
417 case 0:
418 out[i][0] *= p(alpha[l],in[l]);
419 break;
420 case 1:
421 out[i][0] *= dp(alpha[l],in[l]);
422 break;
423 case 2:
424 out[i][0] *= ddp(alpha[l],in[l]);
425 break;
426 default:
427 DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
428 }
429 }
430 }
431 }
432 };
433
439 template<unsigned int dim, int compileTimeOrder>
440 class LagrangeCubeLocalCoefficients
441 : private LagrangeCubeOrderTraits<dim,compileTimeOrder>
442 {
443 using OrderTraits = LagrangeCubeOrderTraits<dim,compileTimeOrder>;
444
445 using OrderTraits::multiindex;
446
448 void setup1d (std::vector<unsigned int>& subEntity)
449 {
450 const unsigned int k = order();
451 assert(k>0);
452
453 unsigned lastIndex=0;
454
455 /* edge and vertex numbering
456 0----0----1
457 */
458
459 // edge (0)
460 subEntity[lastIndex++] = 0; // corner 0
461 for (unsigned i = 0; i < k - 1; ++i)
462 subEntity[lastIndex++] = 0; // inner dofs of element (0)
463
464 subEntity[lastIndex++] = 1; // corner 1
465
466 assert(power(k+1,dim)==lastIndex);
467 }
468
469 void setup2d (std::vector<unsigned int>& subEntity)
470 {
471 const unsigned int k = order();
472 assert(k>0);
473
474 unsigned lastIndex=0;
475
476 // LocalKey: entity number, entity codim, dof indices within each entity
477 /* edge and vertex numbering
478 2----3----3
479 | |
480 | |
481 0 1
482 | |
483 | |
484 0----2----1
485 */
486
487 // lower edge (2)
488 subEntity[lastIndex++] = 0; // corner 0
489 for (unsigned i = 0; i < k - 1; ++i)
490 subEntity[lastIndex++] = 2; // inner dofs of lower edge (2)
491
492 subEntity[lastIndex++] = 1; // corner 1
493
494 // iterate from bottom to top over inner edge dofs
495 for (unsigned e = 0; e < k - 1; ++e) {
496 subEntity[lastIndex++] = 0; // left edge (0)
497 for (unsigned i = 0; i < k - 1; ++i)
498 subEntity[lastIndex++] = 0; // face dofs
499 subEntity[lastIndex++] = 1; // right edge (1)
500 }
501
502 // upper edge (3)
503 subEntity[lastIndex++] = 2; // corner 2
504 for (unsigned i = 0; i < k - 1; ++i)
505 subEntity[lastIndex++] = 3; // inner dofs of upper edge (3)
506
507 subEntity[lastIndex++] = 3; // corner 3
508
509 assert(power(k+1,dim)==lastIndex);
510 }
511
512 void setup3d (std::vector<unsigned int>& subEntity)
513 {
514 const unsigned int k = order();
515 assert(k>0);
516
517 unsigned lastIndex=0;
518#ifndef NDEBUG
519 const unsigned numIndices = power(k+1,dim);
520 const unsigned numFaceIndices = power(k+1,dim-1);
521#endif
522 const unsigned numInnerEdgeDofs = k-1;
523
524 // LocalKey: entity number, entity codim, dof indices within each entity
525 /* edge and vertex numbering
526
527 6---(11)--7 6---------7
528 /| /| /| (5) /|
529 (8)| (9)| / | top / |
530 / (2) / (3) / |(3)bac/k |
531 4---(10)--5 | 4---------5 |
532 | | | | left|(0)| |(1)|right
533 | 2--(7)|---3 | 2-----|---3
534 (0) / (1) / |(2)front | /
535 |(4) |(5) | / (4) | /
536 |/ |/ |/ bottom |/
537 0---(6)---1 0---------1
538 */
539
540 // bottom face (4)
541 lastIndex=0;
542 // lower edge (6)
543 subEntity[lastIndex++] = 0; // corner 0
544 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
545 subEntity[lastIndex++] = 6; // inner dofs of lower edge (6)
546
547 subEntity[lastIndex++] = 1; // corner 1
548
549 // iterate from bottom to top over inner edge dofs
550 for (unsigned e = 0; e < numInnerEdgeDofs; ++e) {
551 subEntity[lastIndex++] = 4; // left edge (4)
552 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
553 subEntity[lastIndex++] = 4; // inner face dofs
554 subEntity[lastIndex++] = 5; // right edge (5)
555 }
556
557 // upper edge (7)
558 subEntity[lastIndex++] = 2; // corner 2
559 for (unsigned i = 0; i < k - 1; ++i)
560 subEntity[lastIndex++] = 7; // inner dofs of upper edge (7)
561 subEntity[lastIndex++] = 3; // corner 3
562
563 assert(numFaceIndices==lastIndex); // added 1 face so far
565
567 for(unsigned f = 0; f < numInnerEdgeDofs; ++f) {
568
569 // lower edge (connecting edges 0 and 1)
570 subEntity[lastIndex++] = 0; // dof on edge 0
571 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
572 subEntity[lastIndex++] = 2; // dof in front face
573 subEntity[lastIndex++] = 1; // dof on edge 1
574
575 // iterate from bottom to top over inner edge dofs
576 for (unsigned e = 0; e < numInnerEdgeDofs; ++e) {
577 subEntity[lastIndex++] = 0; // on left face (0)
578 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
579 subEntity[lastIndex++] = 0; // volume dofs
580 subEntity[lastIndex++] = 1; // right face (1)
581 }
582
583 // upper edge (connecting edges 0 and 1)
584 subEntity[lastIndex++] = 2; // dof on edge 2
585 for (unsigned i = 0; i < numInnerEdgeDofs; ++i)
586 subEntity[lastIndex++] = 3; // dof on rear face (3)
587 subEntity[lastIndex++] = 3; // dof on edge 3
588
589 assert(lastIndex==(f+1+1)*numFaceIndices);
590 }
591
593 // lower edge (10)
594 subEntity[lastIndex++] = 4; // corner 4
595 for (unsigned i = 0; i < k - 1; ++i)
596 subEntity[lastIndex++] = 10; // inner dofs on lower edge (10)
597 subEntity[lastIndex++] = 5; // corner 5
598
599 // iterate from bottom to top over inner edge dofs
600 for (unsigned e = 0; e < k - 1; ++e) {
601 subEntity[lastIndex++] = 8; // left edge (8)
602 for (unsigned i = 0; i < k - 1; ++i)
603 subEntity[lastIndex++] = 5; // face dofs
604 subEntity[lastIndex++] = 9; // right edge (9)
605 }
606
607 // upper edge (11)
608 subEntity[lastIndex++] = 6; // corner 6
609 for (unsigned i = 0; i < k - 1; ++i)
610 subEntity[lastIndex++] = 11; // inner dofs of upper edge (11)
611 subEntity[lastIndex++] = 7; // corner 7
612
613 assert(numIndices==lastIndex);
614 }
615
618 void setup ()
619 {
620 const unsigned int k = order();
621 localKeys_.resize(size());
622
623 if (k==0)
624 {
625 localKeys_[0] = LocalKey(0,0,0);
626 return;
627 }
628
629 if (k==1)
630 {
631 for (std::size_t i=0; i<size(); i++)
632 localKeys_[i] = LocalKey(i,dim,0);
633 return;
634 }
635
636 // Now: the general case
637
638 // Set up array of codimension-per-dof-number
639 std::vector<unsigned int> codim(size());
640
641 for (std::size_t i=0; i<codim.size(); i++)
642 {
643 codim[i] = 0;
644
645 // Codimension gets increased by 1 for each coordinate direction
646 // where dof is on boundary
647 std::array<unsigned int,dim> mIdx = multiindex(i);
648 for (unsigned int j=0; j<dim; j++)
649 if (mIdx[j]==0 or mIdx[j]==k)
650 codim[i]++;
651 }
652
653 // Set up index vector (the index of the dof in the set of dofs of a given subentity)
654 // Algorithm: the 'index' has the same ordering as the dof number 'i'.
655 // To make it consecutive we interpret 'i' in the (k+1)-adic system, omit all digits
656 // that correspond to axes where the dof is on the element boundary, and transform the
657 // rest to the (k-1)-adic system.
658 std::vector<unsigned int> index(size());
659
660 for (std::size_t i=0; i<size(); i++)
661 {
662 index[i] = 0;
663
664 std::array<unsigned int,dim> mIdx = multiindex(i);
665
666 for (int j=dim-1; j>=0; j--)
667 if (mIdx[j]>0 && mIdx[j]<k)
668 index[i] = (k-1)*index[i] + (mIdx[j]-1);
669 }
670
671 // Set up entity and dof numbers for each (supported) dimension separately
672 std::vector<unsigned int> subEntity(size());
673
674 if (dim==1) {
675
676 setup1d(subEntity);
677
678 } else if (dim==2) {
679
680 setup2d(subEntity);
681
682 } else if (dim==3) {
683
684 setup3d(subEntity);
685
686 } else
687 DUNE_THROW(Dune::NotImplemented, "LagrangeCubeLocalCoefficients for order " << k << " and dim == " << dim);
688
689 for (size_t i=0; i<size(); i++)
690 localKeys_[i] = LocalKey(subEntity[i], codim[i], index[i]);
691 }
692
693 public:
695 explicit LagrangeCubeLocalCoefficients (OrderTraits orderTraits)
696 : OrderTraits(orderTraits)
697 {
698 setup();
699 }
700
701 using OrderTraits::size;
702 using OrderTraits::order;
703
705 constexpr const LocalKey& localKey (std::size_t i) const
706 {
707 return localKeys_[i];
708 }
709
710 private:
711 std::vector<LocalKey> localKeys_;
712 };
713
721 template<class D, class R, unsigned int dim, int compileTimeOrder>
722 class LagrangeCubeLocalInterpolation
723 : private LagrangeCubeOrderTraits<dim, compileTimeOrder>
724 {
725 using OrderTraits = LagrangeCubeOrderTraits<dim, compileTimeOrder>;
726 using Traits = typename LagrangeCubeLocalBasis<D,R,dim,compileTimeOrder>::Traits;
727
728 using OrderTraits::multiindex;
729 using OrderTraits::size;
730 using OrderTraits::order;
731
732 public:
734 explicit constexpr LagrangeCubeLocalInterpolation (OrderTraits orderTraits)
735 : OrderTraits(orderTraits)
736 {}
737
745 template<typename F, typename C>
746 constexpr void interpolate (const F& f, std::vector<C>& out) const
747 {
748 const unsigned int k = order();
749 out.resize(size());
750
751 // Specialization for zero-order case
752 if (k==0)
753 {
754 auto center = ReferenceElements<D,dim>::cube().position(0,0);
755 out[0] = f(center);
756 return;
757 }
758
759 typename Traits::DomainType x;
760
761 // Specialization for first-order case
762 if (k==1)
763 {
764 for (unsigned int i=0; i<size(); i++)
765 {
766 // Generate coordinate of the i-th corner of the reference cube
767 for (unsigned int j=0; j<dim; j++)
768 x[j] = (i & (1<<j)) ? 1.0 : 0.0;
769
770 out[i] = f(x);
771 }
772 return;
773 }
774
775 // The general case
776 for (unsigned int i=0; i<size(); i++)
777 {
778 // convert index i to multiindex
779 std::array<unsigned int,dim> alpha(multiindex(i));
780
781 // Generate coordinate of the i-th Lagrange point
782 for (unsigned int j=0; j<dim; j++)
783 x[j] = (1.0*alpha[j])/k;
784
785 out[i] = f(x);
786 }
787 }
788 };
789
790} } // namespace Dune::Impl
791
792namespace Dune
793{
802 template<class D, class R, int dim, int compileTimeOrder = -1>
804 : private Impl::LagrangeCubeOrderTraits<dim, compileTimeOrder>
805 {
806 using OrderTraits = Impl::LagrangeCubeOrderTraits<dim, compileTimeOrder>;
807
808 public:
812 Impl::LagrangeCubeLocalBasis<D,R,dim,compileTimeOrder>,
813 Impl::LagrangeCubeLocalCoefficients<dim,compileTimeOrder>,
814 Impl::LagrangeCubeLocalInterpolation<D,R,dim,compileTimeOrder> >;
815
818 : OrderTraits()
819 , basis_(*this)
820 , coefficients_(*this)
821 , interpolation_(*this)
822 {
823 static_assert(OrderTraits::is_static_order, "Default constructor only allowed for compile-time order >= 0");
824 static_assert(compileTimeOrder >= 0, "Default constructor only allowed for compile-time order >= 0");
825 }
826
828 explicit constexpr LagrangeCubeLocalFiniteElement (int runTimeOrder)
829 : OrderTraits(runTimeOrder)
830 , basis_(*this)
831 , coefficients_(*this)
832 , interpolation_(*this)
833 {
834 if (runTimeOrder < 0)
835 DUNE_THROW(Dune::InvalidStateException, "LagrangeCubeLocalFiniteElement: run-time order must be non-negative!");
836 if (compileTimeOrder >= 0 && compileTimeOrder != runTimeOrder)
837 DUNE_THROW(Dune::InvalidStateException, "LagrangeCubeLocalFiniteElement: Compile-time order must be identical to run-time order!");
838 }
839
842 constexpr const typename Traits::LocalBasisType& localBasis () const
843 {
844 return basis_;
845 }
846
849 constexpr const typename Traits::LocalCoefficientsType& localCoefficients () const
850 {
851 return coefficients_;
852 }
853
856 constexpr const typename Traits::LocalInterpolationType& localInterpolation () const
857 {
858 return interpolation_;
859 }
860
862 using OrderTraits::size;
863
866 static constexpr GeometryType type ()
867 {
868 return GeometryTypes::cube(dim);
869 }
870
871 private:
872 typename Traits::LocalBasisType basis_;
873 typename Traits::LocalCoefficientsType coefficients_;
874 typename Traits::LocalInterpolationType interpolation_;
875 };
876
877} // namespace Dune
878
879#endif // DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGECUBE_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 cubes with arbitrary compile-time dimension and polynomial order.
Definition: lagrangecube.hh:805
constexpr const Traits::LocalInterpolationType & localInterpolation() const
Returns object that evaluates degrees of freedom.
Definition: lagrangecube.hh:856
constexpr LagrangeCubeLocalFiniteElement(int runTimeOrder)
Constructor for dynamic order.
Definition: lagrangecube.hh:828
static constexpr GeometryType type()
The reference element that the local finite element is defined on.
Definition: lagrangecube.hh:866
constexpr const Traits::LocalBasisType & localBasis() const
Returns the local basis, i.e., the set of shape functions.
Definition: lagrangecube.hh:842
constexpr const Traits::LocalCoefficientsType & localCoefficients() const
Returns the assignment of the degrees of freedom to the element subentities.
Definition: lagrangecube.hh:849
constexpr LagrangeCubeLocalFiniteElement()
Constructor for static order.
Definition: lagrangecube.hh:817
Default exception for dummy implementations.
Definition: exceptions.hh:357
Traits for type conversions and type information.
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.
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:462
void interpolate(const F &f, const GFS &gfs, XG &xg)
interpolation from a given grid function
Definition: interpolate.hh:177
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:284
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
constexpr Base power(Base m, Exponent p)
Power method for integer exponents.
Definition: math.hh:75
static const ReferenceElement & cube()
get hypercube reference elements
Definition: referenceelements.hh:168
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)