LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
lagr_fe.h
1#ifndef LF_FE
2#define LF_FE
3/***************************************************************************
4 * LehrFEM++ - A simple C++ finite element libray for teaching
5 * Developed from 2018 at the Seminar of Applied Mathematics of ETH Zurich,
6 * lead developers Dr. R. Casagrande and Prof. R. Hiptmair
7 ***************************************************************************/
8
17#include <lf/assemble/dofhandler.h>
18#include <lf/fe/fe_point.h>
19#include <lf/fe/scalar_reference_finite_element.h>
20
21#include <iostream>
22#include <typeinfo>
23
24namespace lf::uscalfe {
37
56template <typename SCALAR>
59 public:
61 FeLagrangeO1Tria(FeLagrangeO1Tria&&) noexcept = default;
62 FeLagrangeO1Tria& operator=(const FeLagrangeO1Tria&) = default;
63 FeLagrangeO1Tria& operator=(FeLagrangeO1Tria&&) noexcept = default;
64 FeLagrangeO1Tria() = default;
65 ~FeLagrangeO1Tria() override = default;
66
67 [[nodiscard]] base::RefEl RefEl() const override {
68 return base::RefEl::kTria();
69 }
70
71 [[nodiscard]] unsigned Degree() const override { return 1; }
72
76 [[nodiscard]] size_type NumRefShapeFunctions() const override { return 3; }
77
78 // clang-format off
83 // clang-format on
84 [[nodiscard]] size_type NumRefShapeFunctions(dim_t codim) const override {
85 LF_ASSERT_MSG(codim <= 2, "Illegal codim " << codim);
86 return (codim == 2) ? 1 : 0;
87 }
88 // clang-format off
93 // clang-format on
95 dim_t codim, sub_idx_t subidx) const override {
96 LF_ASSERT_MSG(codim <= 2, "Illegal codim " << codim);
97 LF_ASSERT_MSG((codim == 2 && subidx < 3) || (codim == 1 && subidx < 3) ||
98 (codim == 0 && subidx == 0),
99 "codim/subidx out of range.");
100 return (codim == 2) ? 1 : 0;
101 }
102
103 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
104 EvalReferenceShapeFunctions(const Eigen::MatrixXd& refcoords) const override {
105 LF_ASSERT_MSG(refcoords.rows() == 2,
106 "Reference coordinates must be 2-vectors");
107 const size_type n_pts(refcoords.cols());
108 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(
109 3, refcoords.cols());
110 result.row(0) = Eigen::RowVectorXd::Ones(refcoords.cols()) -
111 refcoords.row(0) - refcoords.row(1);
112 result.block(1, 0, 2, refcoords.cols()) = refcoords;
113 return result;
114 }
115
116 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
118 const Eigen::MatrixXd& refcoords) const override {
119 LF_ASSERT_MSG(refcoords.rows() == 2,
120 "Reference coordinates must be 2-vectors");
121 const size_type n_pts(refcoords.cols());
122
123 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(
124 3, 2 * refcoords.cols());
125 result.row(0) =
126 Eigen::RowVectorXd::Constant(2 * static_cast<Eigen::Index>(n_pts), -1);
127 result.row(1) = Eigen::RowVector2d(1., 0.).replicate(1, n_pts);
128 result.row(2) = Eigen::RowVector2d(0., 1.).replicate(1, n_pts);
129 return result;
130 }
131
135 [[nodiscard]] Eigen::MatrixXd EvaluationNodes() const override {
136 return RefEl().NodeCoords();
137 }
138
142 [[nodiscard]] size_type NumEvaluationNodes() const override {
143 return RefEl().NumNodes();
144 }
145};
146
164template <typename SCALAR>
166 : public lf::fe::ScalarReferenceFiniteElement<SCALAR> {
167 public:
169 FeLagrangeO1Quad(FeLagrangeO1Quad&&) noexcept = default;
170 FeLagrangeO1Quad& operator=(const FeLagrangeO1Quad&) = default;
171 FeLagrangeO1Quad& operator=(FeLagrangeO1Quad&&) noexcept = default;
172 FeLagrangeO1Quad() = default;
173 ~FeLagrangeO1Quad() override = default;
174
175 [[nodiscard]] base::RefEl RefEl() const override {
176 return base::RefEl::kQuad();
177 }
178
179 [[nodiscard]] unsigned Degree() const override { return 1; }
180
184 [[nodiscard]] size_type NumRefShapeFunctions() const override { return 4; }
185
186 // clang-format off
191 // clang-format on
192 [[nodiscard]] size_type NumRefShapeFunctions(dim_t codim) const override {
193 return (codim == 2) ? 1 : 0;
194 }
195
196 // clang-format off
201 // clang-format on
203 dim_t codim, sub_idx_t subidx) const override {
204 LF_ASSERT_MSG((codim == 2 && subidx < 4) || (codim == 1 && subidx < 4) ||
205 (codim == 0 && subidx == 0),
206 "codim/subidx out of range.");
207 return (codim == 2) ? 1 : 0;
208 }
209
210 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
211 EvalReferenceShapeFunctions(const Eigen::MatrixXd& refcoords) const override {
212 LF_ASSERT_MSG(refcoords.rows() == 2,
213 "Reference coordinates must be 2-vectors");
214
215 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(
216 4, refcoords.cols());
217
218 result.row(0) =
219 ((1 - refcoords.row(0).array()) * (1 - refcoords.row(1).array()))
220 .matrix();
221 result.row(1) =
222 ((refcoords.row(0).array()) * (1 - refcoords.row(1).array())).matrix();
223 result.row(2) =
224 ((refcoords.row(0).array()) * (refcoords.row(1).array())).matrix();
225 result.row(3) =
226 ((1 - refcoords.row(0).array()) * (refcoords.row(1).array())).matrix();
227
228 return result;
229 }
230
231 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
233 const Eigen::MatrixXd& refcoords) const override {
234 LF_ASSERT_MSG(refcoords.rows() == 2,
235 "Reference coordinates must be 2-vectors");
236 const size_type n_pts(refcoords.cols());
237
238 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(4, 2 * n_pts);
239
240 // reshape the result matrix into a 8xn_pts matrix
241 Eigen::Map<Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic>,
242 Eigen::AutoAlign>
243 temp(&result(0, 0), 8, n_pts);
244 temp.row(0) = refcoords.row(1).array() - 1.0;
245 temp.row(1) = 1.0 - refcoords.row(1).array();
246 temp.row(2) = refcoords.row(1).array();
247 temp.row(3) = -refcoords.row(1).array();
248 temp.row(4) = refcoords.row(0).array() - 1.0;
249 temp.row(5) = -refcoords.row(0).array();
250 temp.row(6) = refcoords.row(0).array();
251 temp.row(7) = 1.0 - refcoords.row(0).array();
252
253 return result;
254 }
255
256 [[nodiscard]] Eigen::MatrixXd EvaluationNodes() const override {
257 return RefEl().NodeCoords();
258 }
259
260 [[nodiscard]] size_type NumEvaluationNodes() const override {
261 return RefEl().NumNodes();
262 }
263};
264
283template <typename SCALAR>
285 : public lf::fe::ScalarReferenceFiniteElement<SCALAR> {
286 public:
289 FeLagrangeO1Segment& operator=(const FeLagrangeO1Segment&) = default;
290 FeLagrangeO1Segment& operator=(FeLagrangeO1Segment&&) noexcept = default;
292 ~FeLagrangeO1Segment() override = default;
293
294 [[nodiscard]] base::RefEl RefEl() const override {
295 return base::RefEl::kSegment();
296 }
297
298 [[nodiscard]] unsigned Degree() const override { return 1; }
299
303 [[nodiscard]] size_type NumRefShapeFunctions() const override { return 2; }
304
305 // clang-format off
309 // clang-format on
310 [[nodiscard]] size_type NumRefShapeFunctions(dim_t codim) const override {
311 LF_ASSERT_MSG(codim <= 1, "Illegal codim " << codim);
312 return (codim == 1) ? 1 : 0;
313 }
314
315 // clang-format off
319 // clang-format on
321 dim_t codim, sub_idx_t subidx) const override {
322 LF_ASSERT_MSG(codim <= 1, "Illegal codim " << codim);
323 LF_ASSERT_MSG((codim == 1 && subidx < 2) || (codim == 0 && subidx == 0),
324 "codim/subidx out of range.");
325 return (codim == 1) ? 1 : 0;
326 }
327
328 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
329 EvalReferenceShapeFunctions(const Eigen::MatrixXd& refcoords) const override {
330 LF_ASSERT_MSG(refcoords.rows() == 1,
331 "Reference coordinates must be 1-vectors");
332 const size_type n_pts(refcoords.cols());
333
334 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(
335 2, refcoords.cols());
336 result.row(0) = Eigen::RowVectorXd::Constant(1, n_pts, 1.0) - refcoords;
337 result.row(1) = refcoords;
338
339 return result;
340 }
341
342 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
344 const Eigen::MatrixXd& refcoords) const override {
345 LF_ASSERT_MSG(refcoords.rows() == 1,
346 "Reference coordinates must be 1-vectors");
347 const size_type n_pts(refcoords.cols());
348
349 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(
350 2, refcoords.cols());
351 result.row(0) = Eigen::MatrixXd::Constant(1, n_pts, -1.0);
352 result.row(1) = Eigen::MatrixXd::Constant(1, n_pts, 1.0);
353
354 return result;
355 }
356
360 [[nodiscard]] Eigen::MatrixXd EvaluationNodes() const override {
361 return RefEl().NodeCoords();
362 }
363
367 [[nodiscard]] size_type NumEvaluationNodes() const override {
368 return RefEl().NumNodes();
369 }
370};
371
399template <typename SCALAR>
401 : public lf::fe::ScalarReferenceFiniteElement<SCALAR> {
402 public:
404 FeLagrangeO2Tria(FeLagrangeO2Tria&&) noexcept = default;
405 FeLagrangeO2Tria& operator=(const FeLagrangeO2Tria&) = default;
406 FeLagrangeO2Tria& operator=(FeLagrangeO2Tria&&) noexcept = default;
407 FeLagrangeO2Tria() = default;
408 ~FeLagrangeO2Tria() override = default;
409
412 [[nodiscard]] base::RefEl RefEl() const override {
413 return base::RefEl::kTria();
414 }
415
419 [[nodiscard]] unsigned Degree() const override { return 2; }
420
426 [[nodiscard]] size_type NumRefShapeFunctions() const override { return 6; }
427
434 [[nodiscard]] size_type NumRefShapeFunctions(dim_t codim) const override {
435 LF_ASSERT_MSG(codim <= 2, "Illegal codim" << codim);
436 return (codim == 0) ? 0 : 1;
437 }
438
446 dim_t codim, sub_idx_t /*subidx*/) const override {
447 LF_ASSERT_MSG(codim <= 2, "Illegal codim" << codim);
448 return (codim == 0) ? 0 : 1;
449 }
450
458 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
459 EvalReferenceShapeFunctions(const Eigen::MatrixXd& refcoords) const override {
460 LF_ASSERT_MSG(refcoords.rows() == 2,
461 "Reference coordinates must be 2-vectors");
462 // Number of evaluation points
463 const size_type n_pts(refcoords.cols());
464 // Result returned in 6 x P matrix
465 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(6, n_pts);
466 // Convert to array type for componentwise operations
467 auto x0 = refcoords.row(0).array();
468 auto x1 = refcoords.row(1).array();
469 // Compound evaluation of formulas for reference shape functions
470 result.row(0) = (2.0 * (1 - x0 - x1) * (0.5 - x0 - x1)).matrix();
471 result.row(1) = (2.0 * x0 * (x0 - 0.5)).matrix();
472 result.row(2) = (2.0 * x1 * (x1 - 0.5)).matrix();
473 result.row(3) = (4.0 * (1 - x0 - x1) * x0).matrix();
474 result.row(4) = (4.0 * x0 * x1).matrix();
475 result.row(5) = (4.0 * (1 - x0 - x1) * x1).matrix();
476 return result;
477 }
478
485 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
487 const Eigen::MatrixXd& refcoords) const override {
488 LF_ASSERT_MSG(refcoords.rows() == 2,
489 "Reference coordinates must be 2-vectors");
490 // Number of evaluation points
491 const size_type n_pts(refcoords.cols());
492 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(6, 2 * n_pts);
493 // reshape into a 12xn_pts matrix.
494 Eigen::Map<Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic>,
495 Eigen::AutoAlign>
496 temp(result.data(), 12, n_pts);
497 // Convert to array type for componentwise operations
498 auto x0 = refcoords.row(0).array();
499 auto x1 = refcoords.row(1).array();
500 // d/dx_0 for all reference shape functions
501 temp.row(0) = -3.0 + 4.0 * x0 + 4.0 * x1;
502 temp.row(1) = 4.0 * x0 - 1.0;
503 temp.row(2) = 0.0;
504 temp.row(3) = 4 - 0 - 8.0 * x0 - 4.0 * x1;
505 temp.row(4) = 4.0 * x1;
506 temp.row(5) = -4.0 * x1;
507 // d/dx_1 for the reference shape functions
508 temp.row(6) = -3.0 + 4.0 * x0 + 4.0 * x1;
509 temp.row(7) = 0.0;
510 temp.row(8) = 4.0 * x1 - 1.0;
511 temp.row(9) = -4 * x0;
512 temp.row(10) = 4.0 * x0;
513 temp.row(11) = 4.0 - 8.0 * x1 - 4.0 * x0;
514 return result;
515 }
516
529 [[nodiscard]] Eigen::MatrixXd EvaluationNodes() const override {
530 // clang-format off
531 Eigen::Matrix<double, 2,6> nodes;
532 // Reference coordinates of evaluation nodes in the columns of a 2 x 6 - matrix.
533 nodes << 0.0, 1.0, 0.0, 0.5, 0.5, 0.0,
534 0.0, 0.0, 1.0, 0.0, 0.5, 0.5;
535 // clang-format on
536 return nodes;
537 }
538
543 [[nodiscard]] size_type NumEvaluationNodes() const override {
544 return NumRefShapeFunctions();
545 }
546};
547
562template <typename SCALAR>
564 : public lf::fe::ScalarReferenceFiniteElement<SCALAR> {
565 public:
568 FeLagrangeO2Segment& operator=(const FeLagrangeO2Segment&) = default;
569 FeLagrangeO2Segment& operator=(FeLagrangeO2Segment&&) noexcept = default;
571 ~FeLagrangeO2Segment() override = default;
572
573 [[nodiscard]] base::RefEl RefEl() const override {
574 return base::RefEl::kSegment();
575 }
576
579 [[nodiscard]] unsigned Degree() const override { return 2; }
580
586 [[nodiscard]] size_type NumRefShapeFunctions() const override { return 3; }
587
596 [[nodiscard]] size_type NumRefShapeFunctions(dim_t codim) const override {
597 LF_ASSERT_MSG(codim <= 1, "Illegal codim " << codim);
598 return 1;
599 }
600
607 dim_t codim, sub_idx_t /*subidx*/) const override {
608 LF_ASSERT_MSG(codim <= 1, "Illegal codim " << codim);
609 return 1;
610 }
611
626 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
627 EvalReferenceShapeFunctions(const Eigen::MatrixXd& refcoords) const override {
628 LF_ASSERT_MSG(refcoords.rows() == 1,
629 "Reference coordinates must be 1-vectors");
630 const size_type n_pts(refcoords.cols());
631 // Matrix for returning values of local shape functions
632 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(3, n_pts);
633 // Convert into an array type for componentwise operations
634 auto x = refcoords.row(0).array();
635 // local shape functions belonging to the endpoints
636 result.row(0) = 2.0 * (1.0 - x) * (0.5 - x);
637 result.row(1) = 2.0 * x * (x - 0.5);
638 // local shape function sitting at midpoint (belonging to segment)
639 result.row(2) = 4.0 * (1.0 - x) * x;
640 return result;
641 }
642
652 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
654 const Eigen::MatrixXd& refcoords) const override {
655 LF_ASSERT_MSG(refcoords.rows() == 1,
656 "Reference coordinates must be 1-vectors");
657 const size_type n_pts(refcoords.cols());
658 // Matrix for returning the derivatives of the local shape functions
659 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(3, n_pts);
660 // Convert to array for componentwise operations
661 auto x = refcoords.row(0).array();
662 // LSF at endpoints
663 result.row(0) = (4.0 * x - 3.0).matrix();
664 result.row(1) = (4.0 * x - 1.0).matrix();
665 // LSF at midpoint:
666 result.row(2) = (4.0 - 8.0 * x).matrix();
667 return result;
668 }
669
673 [[nodiscard]] Eigen::MatrixXd EvaluationNodes() const override {
674 Eigen::Matrix<double, 1, 3> nodes;
675 // Reference coordinates of the three evaluation nodes
676 nodes << 0.0, 1.0, 0.5;
677 return nodes;
678 }
679
683 [[nodiscard]] size_type NumEvaluationNodes() const override {
684 return NumRefShapeFunctions();
685 }
686};
687
712template <typename SCALAR>
714 : public lf::fe::ScalarReferenceFiniteElement<SCALAR> {
715 public:
717 FeLagrangeO2Quad(FeLagrangeO2Quad&&) noexcept = default;
718 FeLagrangeO2Quad& operator=(const FeLagrangeO2Quad&) = default;
719 FeLagrangeO2Quad& operator=(FeLagrangeO2Quad&&) noexcept = default;
720 FeLagrangeO2Quad() = default;
721 ~FeLagrangeO2Quad() override = default;
722
724 [[nodiscard]] base::RefEl RefEl() const override {
725 return base::RefEl::kQuad();
726 }
727
730 [[nodiscard]] unsigned Degree() const override { return 2; }
731
738 [[nodiscard]] size_type NumRefShapeFunctions() const override { return 9; }
739
745 [[nodiscard]] size_type NumRefShapeFunctions(dim_t codim) const override {
746 LF_ASSERT_MSG(codim <= 2, "Illegal codim" << codim);
747 return 1;
748 }
749
757 dim_t codim, sub_idx_t /*subidx*/) const override {
758 LF_ASSERT_MSG(codim <= 2, "Illegal codim" << codim);
759 return 1;
760 }
761
779 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
780 EvalReferenceShapeFunctions(const Eigen::MatrixXd& refcoords) const override {
781 LF_ASSERT_MSG(refcoords.rows() == 2,
782 "Reference coordinates must be 2-vectors");
783 // Number of evaluation points
784 const size_type n_pts(refcoords.cols());
785 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(9, n_pts);
786
787 // evaluate "segment" shape functions (b^j and b^l)
788 Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic> segment_x0_eval =
789 (krsf_segment_.EvalReferenceShapeFunctions(refcoords.row(0))).array();
790 Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic> segment_x1_eval =
791 (krsf_segment_.EvalReferenceShapeFunctions(refcoords.row(1))).array();
792
793 // Evaluate basis functions using the tensor product structure
794 for (int i = 0; i < 9; ++i) {
795 result.row(i) = (segment_x0_eval.row(ksegment_to_quad_mapping_(i, 0)) *
796 segment_x1_eval.row(ksegment_to_quad_mapping_(i, 1)))
797 .matrix();
798 }
799 return result;
800 }
801
811 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
813 const Eigen::MatrixXd& refcoords) const override {
814 LF_ASSERT_MSG(refcoords.rows() == 2,
815 "Reference coordinates must be 2-vectors");
816
817 const size_type n_pts(refcoords.cols());
818 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(9, 2 * n_pts);
819
820 // reshape into a 18xn_pts matrix.
821 Eigen::Map<Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic>,
822 Eigen::AutoAlign>
823 temp(&result(0, 0), 18, n_pts);
824
825 // evaluate "segment" shape functions (b^j and b^l)
826 Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic> segment_x0_eval =
827 (krsf_segment_.EvalReferenceShapeFunctions(refcoords.row(0))).array();
828 Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic> segment_x1_eval =
829 (krsf_segment_.EvalReferenceShapeFunctions(refcoords.row(1))).array();
830
831 // evaluate derivatives of "segment" shape functions (b^j and
832 // b^l)
833 Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic> segment_x0_grad =
834 (krsf_segment_.GradientsReferenceShapeFunctions(refcoords.row(0)))
835 .array();
836 Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic> segment_x1_grad =
837 (krsf_segment_.GradientsReferenceShapeFunctions(refcoords.row(1)))
838 .array();
839
840 // evaluate gradients using the product rule and
841 // the tensor product structure of the basis functions.
842 for (int i = 0; i < 9; ++i) {
843 // d/dx
844 temp.row(i) = (segment_x0_grad.row(ksegment_to_quad_mapping_(i, 0)) *
845 segment_x1_eval.row(ksegment_to_quad_mapping_(i, 1)))
846 .matrix();
847 // d/dy
848 temp.row(i + 9) = (segment_x1_grad.row(ksegment_to_quad_mapping_(i, 1)) *
849 segment_x0_eval.row(ksegment_to_quad_mapping_(i, 0)))
850 .matrix();
851 }
852 return result;
853 }
854
864 [[nodiscard]] Eigen::MatrixXd EvaluationNodes() const override {
865 Eigen::Matrix<double, 2, 9> nodes;
866
867 // clang-format off
868 nodes << 0.0, 1.0, 1.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.5,
869 0.0, 0.0, 1.0, 1.0, 0.0, 0.5, 1.0, 0.5, 0.5;
870 // clang-format on
871 return nodes;
872 }
873
878 [[nodiscard]] size_type NumEvaluationNodes() const override {
879 return NumRefShapeFunctions();
880 }
881
882 private:
887
905 const static Eigen::Matrix<int, 9, 2> ksegment_to_quad_mapping_;
906};
907
908template <typename SCALAR>
911
912template <typename SCALAR>
913const Eigen::Matrix<int, 9, 2>
915 // clang-format off
916 (Eigen::Matrix<int,9,2>() << 0, 0,
917 1, 0,
918 1, 1,
919 0, 1,
920 2, 0,
921 1, 2,
922 2, 1,
923 0, 2,
924 2, 2).finished();
925// clang-format on
926
944template <typename SCALAR>
946 : public lf::fe::ScalarReferenceFiniteElement<SCALAR> {
947 public:
949 FeLagrangeO3Tria(FeLagrangeO3Tria&&) noexcept = default;
950 FeLagrangeO3Tria& operator=(const FeLagrangeO3Tria&) = default;
951 FeLagrangeO3Tria& operator=(FeLagrangeO3Tria&&) noexcept = default;
952 FeLagrangeO3Tria() = default;
953 ~FeLagrangeO3Tria() override = default;
954
955 [[nodiscard]] base::RefEl RefEl() const override {
956 return base::RefEl::kTria();
957 }
958
962 [[nodiscard]] unsigned Degree() const override { return 3; }
963
969 [[nodiscard]] size_type NumRefShapeFunctions() const override { return 10; }
970
971 // clang-format off
978 // clang-format on
979 [[nodiscard]] size_type NumRefShapeFunctions(dim_t codim) const override {
980 LF_ASSERT_MSG(codim <= 2, "Illegal codim" << codim);
981 return (codim == 1) ? 2 : 1;
982 }
983
984 // clang-format off
991 // clang-format on
993 dim_t codim, sub_idx_t subidx) const override {
994 LF_ASSERT_MSG(codim <= 2, "Illegal codim" << codim);
995 LF_ASSERT_MSG((codim == 2 && subidx < 3) || (codim == 1 && subidx < 3) ||
996 (codim == 0 && subidx == 0),
997 "codim/subidx out of range.");
998 return (codim == 1) ? 2 : 1;
999 }
1000
1009 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
1010 EvalReferenceShapeFunctions(const Eigen::MatrixXd& refcoords) const override {
1011 LF_ASSERT_MSG(refcoords.rows() == 2,
1012 "Reference coordinates must be 2-vectors");
1013 // number of points for which evaluation is desired
1014 const size_type n_pts(refcoords.cols());
1015 // Matrix for returning the values of the reference shape functions at the
1016 // specified nodes. Each row corresponds to a reference shape function
1017 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(10, n_pts);
1018
1019 // evaluation of the barycentric coordinate functions
1020 const Eigen::Array<double, 1, Eigen::Dynamic> lambda0 =
1021 1 - refcoords.row(0).array() - refcoords.row(1).array();
1022 const Eigen::Array<double, 1, Eigen::Dynamic> lambda1 =
1023 refcoords.row(0).array();
1024 const Eigen::Array<double, 1, Eigen::Dynamic> lambda2 =
1025 refcoords.row(1).array();
1026
1027 // evaluation of the shape functions
1028 // The LSF associated with vertices
1029 result.row(0) = 4.5 * lambda0 * (lambda0 - 1 / 3.0) * (lambda0 - 2 / 3.0);
1030 result.row(1) = 4.5 * lambda1 * (lambda1 - 1 / 3.0) * (lambda1 - 2 / 3.0);
1031 result.row(2) = 4.5 * lambda2 * (lambda2 - 1 / 3.0) * (lambda2 - 2 / 3.0);
1032 // Six LSF are attached to edges
1033 result.row(3) = 13.5 * lambda0 * lambda1 * (lambda0 - 1 / 3.0);
1034 result.row(4) = 13.5 * lambda0 * lambda1 * (lambda1 - 1 / 3.0);
1035 result.row(5) = 13.5 * lambda1 * lambda2 * (lambda1 - 1 / 3.0);
1036 result.row(6) = 13.5 * lambda1 * lambda2 * (lambda2 - 1 / 3.0);
1037 result.row(7) = 13.5 * lambda2 * lambda0 * (lambda2 - 1 / 3.0);
1038 result.row(8) = 13.5 * lambda2 * lambda0 * (lambda0 - 1 / 3.0);
1039 // LSF associated with the cell: cubic bubble function
1040 result.row(9) = 27.0 * lambda0 * lambda1 * lambda2;
1041
1042 return result;
1043 }
1044
1052 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
1054 const Eigen::MatrixXd& refcoords) const override {
1055 LF_ASSERT_MSG(refcoords.rows() == 2,
1056 "Reference coordinates must be 2-vectors");
1057
1058 const size_type n_pts(refcoords.cols());
1059 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(10, 2 * n_pts);
1060
1061 // reshape into a 20xn_pts matrix.
1062 Eigen::Map<Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic>,
1063 Eigen::AutoAlign>
1064 temp(&result(0, 0), 20, n_pts);
1065
1066 // evaulate barycentric coordinate functions:
1067 const Eigen::Array<double, 1, Eigen::Dynamic> l0 =
1068 1 - refcoords.row(0).array() - refcoords.row(1).array();
1069 auto l1 = refcoords.row(0).array();
1070 auto l2 = refcoords.row(1).array();
1071
1072 // vertex-associated LSF: numbers 0,1,2
1073 temp.row(0) = -4.5 * ((l0 - 1 / 3.0) * (l0 - 2 / 3.0) +
1074 l0 * (l0 - 2 / 3.0) + l0 * (l0 - 1 / 3.0));
1075 temp.row(10) = -4.5 * ((l0 - 1 / 3.0) * (l0 - 2 / 3.0) +
1076 l0 * (l0 - 2 / 3.0) + l0 * (l0 - 1 / 3.0));
1077 temp.row(1) = 4.5 * ((l1 - 1 / 3.0) * (l1 - 2 / 3.0) + l1 * (l1 - 2 / 3.0) +
1078 l1 * (l1 - 1 / 3.0));
1079 temp.row(11) = 0.0;
1080 temp.row(2) = 0.0;
1081 temp.row(12) = 4.5 * ((l2 - 1 / 3.0) * (l2 - 2 / 3.0) +
1082 l2 * (l2 - 2 / 3.0) + l2 * (l2 - 1 / 3.0));
1083 // edge-associated basis functions
1084 temp.row(3) = 13.5 * (-l1 * (l0 - 1 / 3.0) + l0 * (l0 - 1 / 3.0) - l0 * l1);
1085 temp.row(13) = -13.5 * (l1 * (l0 - 1 / 3.0) + l0 * l1);
1086 temp.row(4) = 13.5 * (-l1 * (l1 - 1 / 3.0) + l0 * (l1 - 1 / 3.0) + l0 * l1);
1087 temp.row(14) = -13.5 * l1 * (l1 - 1 / 3.0);
1088 temp.row(5) = 13.5 * (l2 * (l1 - 1 / 3.0) + l1 * l2);
1089 temp.row(15) = 13.5 * (l1 * (l1 - 1 / 3.0));
1090 temp.row(6) = 13.5 * (l2 * (l2 - 1 / 3.0));
1091 temp.row(16) = 13.5 * (l1 * (l2 - 1 / 3.0) + l1 * l2);
1092 temp.row(7) = -13.5 * l2 * (l2 - 1 / 3.0);
1093 temp.row(17) = 13.5 * (l0 * (l2 - 1 / 3.0) - l2 * (l2 - 1 / 3.0) + l0 * l2);
1094 temp.row(8) = -13.5 * (l2 * (l0 - 1 / 3.0) + l2 * l0);
1095 temp.row(18) = 13.5 * (l0 * (l0 - 1 / 3.0) - l2 * (l0 - 1 / 3.0) - l2 * l0);
1096 // cell-associated LSF
1097 temp.row(9) = 27.0 * (-l1 * l2 + l0 * l2);
1098 temp.row(19) = 27.0 * (-l1 * l2 + l0 * l1);
1099
1100 return result;
1101 }
1102
1115 [[nodiscard]] Eigen::MatrixXd EvaluationNodes() const override {
1116 Eigen::Matrix<double, 2, 10> nodes;
1117 Eigen::Matrix<double, 2, 3> vertices;
1118 Eigen::Matrix<double, 2, 6> edges;
1119 Eigen::Matrix<double, 2, 1> center;
1120
1121 // clang-format off
1122 vertices << 0.0, 1.0, 0.0,
1123 0.0, 0.0, 1.0;
1124 edges << 1.0, 2.0, 2.0, 1.0, 0.0, 0.0,
1125 0.0, 0.0, 1.0, 2.0, 2.0, 1.0;
1126 center << 1.0 / 3.0,
1127 1.0 / 3.0;
1128 // clang-format on
1129 edges *= 1.0 / 3.0;
1130 nodes << vertices, edges, center;
1131
1132 return nodes;
1133 }
1137 [[nodiscard]] size_type NumEvaluationNodes() const override {
1138 return NumRefShapeFunctions();
1139 }
1140};
1141
1154template <typename SCALAR>
1156 : public lf::fe::ScalarReferenceFiniteElement<SCALAR> {
1157 public:
1160 FeLagrangeO3Segment& operator=(const FeLagrangeO3Segment&) = default;
1161 FeLagrangeO3Segment& operator=(FeLagrangeO3Segment&&) noexcept = default;
1163 ~FeLagrangeO3Segment() override = default;
1164
1165 [[nodiscard]] base::RefEl RefEl() const override {
1166 return base::RefEl::kSegment();
1167 }
1168
1169 [[nodiscard]] unsigned Degree() const override { return 3; }
1170
1171 [[nodiscard]] size_type NumRefShapeFunctions() const override { return 4; }
1172
1173 // clang-format off
1178 // clang-format on
1179 [[nodiscard]] size_type NumRefShapeFunctions(dim_t codim) const override {
1180 LF_ASSERT_MSG(codim <= 1, "Illegal codim " << codim);
1181 return (codim == 0) ? 2 : 1;
1182 }
1183
1184 // clang-format off
1189 // clang-format on
1191 dim_t codim, sub_idx_t subidx) const override {
1192 LF_ASSERT_MSG(codim <= 1, "Illegal codim " << codim);
1193 LF_ASSERT_MSG((codim == 1 && subidx < 2) || (codim == 0 && subidx == 0),
1194 "codim/subidx out of range.");
1195 return (codim == 0) ? 2 : 1;
1196 }
1197
1198 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
1199 EvalReferenceShapeFunctions(const Eigen::MatrixXd& refcoords) const override {
1200 LF_ASSERT_MSG(refcoords.rows() == 1,
1201 "Reference coordinates must be 1-vectors");
1202 const size_type n_pts(refcoords.cols());
1203
1204 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(4, n_pts);
1205
1206 auto x = refcoords.row(0).array();
1207
1208 // endpoints
1209 result.row(0) = 4.5 * (1.0 - x) * (1.0 / 3.0 - x) * (2.0 / 3.0 - x);
1210 result.row(1) = 4.5 * x * (x - 1.0 / 3.0) * (x - 2.0 / 3.0);
1211
1212 // interior
1213 result.row(2) = 13.5 * x * (1 - x) * (2.0 / 3.0 - x);
1214 result.row(3) = 13.5 * x * (1 - x) * (x - 1.0 / 3.0);
1215
1216 return result;
1217 }
1218
1219 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
1221 const Eigen::MatrixXd& refcoords) const override {
1222 LF_ASSERT_MSG(refcoords.rows() == 1,
1223 "Reference coordinates must be 1-vectors");
1224 const size_type n_pts(refcoords.cols());
1225
1226 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(4, n_pts);
1227
1228 auto x = refcoords.row(0).array();
1229
1230 // endpoints
1231 result.row(0) = -13.5 * x * x + 18.0 * x - 5.5;
1232 result.row(1) = 13.5 * x * x - 9.0 * x + 1.0;
1233 // interior
1234 result.row(2) = 40.5 * x * x - 45 * x + 9;
1235 result.row(3) = -40.5 * x * x + 36 * x - 4.5;
1236
1237 return result;
1238 }
1239
1240 [[nodiscard]] Eigen::MatrixXd EvaluationNodes() const override {
1241 Eigen::Matrix<double, 1, 4> nodes;
1242 nodes << 0.0, 1.0, 1.0 / 3.0, 2.0 / 3.0;
1243 return nodes;
1244 }
1245
1249 [[nodiscard]] size_type NumEvaluationNodes() const override {
1250 return NumRefShapeFunctions();
1251 }
1252};
1253
1278template <typename SCALAR>
1280 : public lf::fe::ScalarReferenceFiniteElement<SCALAR> {
1281 public:
1283 FeLagrangeO3Quad(FeLagrangeO3Quad&&) noexcept = default;
1284 FeLagrangeO3Quad& operator=(const FeLagrangeO3Quad&) = default;
1285 FeLagrangeO3Quad& operator=(FeLagrangeO3Quad&&) noexcept = default;
1286 FeLagrangeO3Quad() = default;
1287 ~FeLagrangeO3Quad() override = default;
1288
1289 [[nodiscard]] base::RefEl RefEl() const override {
1290 return base::RefEl::kQuad();
1291 }
1292
1293 [[nodiscard]] unsigned Degree() const override { return 3; }
1294
1295 [[nodiscard]] size_type NumRefShapeFunctions() const override { return 16; }
1296
1297 // clang-format off
1302 // clang-format on
1303 [[nodiscard]] size_type NumRefShapeFunctions(dim_t codim) const override {
1304 switch (codim) {
1305 case 0:
1306 return 4;
1307 case 1:
1308 return 2;
1309 case 2:
1310 return 1;
1311 default:
1312 LF_ASSERT_MSG(false, "Illegal codim" << codim);
1313 }
1314 return 0;
1315 }
1316
1317 // clang-format off
1322 // clang-format on
1324 dim_t codim, sub_idx_t subidx) const override {
1325 LF_ASSERT_MSG((codim == 2 && subidx < 4) || (codim == 1 && subidx < 4) ||
1326 (codim == 0 && subidx == 0),
1327 "codim/subidx out of range.");
1328 return NumRefShapeFunctions(codim);
1329 }
1330
1331 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
1332 EvalReferenceShapeFunctions(const Eigen::MatrixXd& refcoords) const override {
1333 LF_ASSERT_MSG(refcoords.rows() == 2,
1334 "Reference coordinates must be 2-vectors");
1335
1336 const size_type n_pts(refcoords.cols());
1337 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(
1338 16, refcoords.cols());
1339
1340 // evaluate "segment" shape functions (b^j and b^l)
1341 Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic> segment_x0_eval =
1342 (krsf_segment_.EvalReferenceShapeFunctions(refcoords.row(0))).array();
1343 Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic> segment_x1_eval =
1344 (krsf_segment_.EvalReferenceShapeFunctions(refcoords.row(1))).array();
1345
1346 // Evaluate basis functions using the tensor product structure
1347 for (int i = 0; i < 16; ++i) {
1348 result.row(i) = (segment_x0_eval.row(ksegment_to_quad_mapping_(i, 0)) *
1349 segment_x1_eval.row(ksegment_to_quad_mapping_(i, 1)))
1350 .matrix();
1351 }
1352 return result;
1353 }
1354
1355 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
1357 const Eigen::MatrixXd& refcoords) const override {
1358 LF_ASSERT_MSG(refcoords.rows() == 2,
1359 "Reference coordinates must be 2-vectors");
1360
1361 const size_type n_pts(refcoords.cols());
1362 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> result(16, 2 * n_pts);
1363
1364 // reshape into a 32xn_pts matrix.
1365 Eigen::Map<Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic>,
1366 Eigen::AutoAlign>
1367 temp(&result(0, 0), 32, n_pts);
1368
1369 // evaluate "segment" shape functions (b^j and b^l)
1370 Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic> segment_x0_eval =
1371 (krsf_segment_.EvalReferenceShapeFunctions(refcoords.row(0))).array();
1372 Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic> segment_x1_eval =
1373 (krsf_segment_.EvalReferenceShapeFunctions(refcoords.row(1))).array();
1374
1375 // evaluate derivatives of "segment" shape functions (b^j and
1376 // b^l)
1377 Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic> segment_x0_grad =
1378 (krsf_segment_.GradientsReferenceShapeFunctions(refcoords.row(0)))
1379 .array();
1380 Eigen::Array<SCALAR, Eigen::Dynamic, Eigen::Dynamic> segment_x1_grad =
1381 (krsf_segment_.GradientsReferenceShapeFunctions(refcoords.row(1)))
1382 .array();
1383
1384 // evaluate gradients using the product rule and
1385 // the tensor product structure of the basis functions.
1386 for (int i = 0; i < 16; ++i) {
1387 // d/dx
1388 temp.row(i) = (segment_x0_grad.row(ksegment_to_quad_mapping_(i, 0)) *
1389 segment_x1_eval.row(ksegment_to_quad_mapping_(i, 1)))
1390 .matrix();
1391 // d/dy
1392 temp.row(i + 16) = (segment_x1_grad.row(ksegment_to_quad_mapping_(i, 1)) *
1393 segment_x0_eval.row(ksegment_to_quad_mapping_(i, 0)))
1394 .matrix();
1395 }
1396
1397 return result;
1398 }
1399
1400 [[nodiscard]] Eigen::MatrixXd EvaluationNodes() const override {
1401 Eigen::Matrix<double, 2, 16> nodes;
1402
1403 Eigen::Matrix<double, 2, 4> vertices;
1404 Eigen::Matrix<double, 2, 8> midpoints;
1405 Eigen::Matrix<double, 2, 4> interior;
1406
1407 // clang-format off
1408 vertices << 0.0, 1.0, 1.0, 0.0,
1409 0.0, 0.0, 1.0, 1.0;
1410 midpoints << 1.0, 2.0, 3.0, 3.0, 2.0, 1.0, 0.0, 0.0,
1411 0.0, 0.0, 1.0, 2.0, 3.0, 3.0, 2.0, 1.0;
1412 interior << 1.0, 2.0, 2.0, 1.0,
1413 1.0, 1.0, 2.0, 2.0;
1414 // clang-format on
1415 midpoints *= 1.0 / 3.0;
1416 interior *= 1.0 / 3.0;
1417
1418 nodes << vertices, midpoints, interior;
1419
1420 return nodes;
1421 }
1422
1426 [[nodiscard]] size_type NumEvaluationNodes() const override {
1427 return NumRefShapeFunctions();
1428 }
1429
1430 private:
1435
1453 const static Eigen::Matrix<int, 16, 2> ksegment_to_quad_mapping_;
1454};
1455
1456template <typename SCALAR>
1459
1460template <typename SCALAR>
1461const Eigen::Matrix<int, 16, 2>
1463 // clang-format off
1464 (Eigen::Matrix<int,16,2>() << 0, 0,
1465 1, 0,
1466 1, 1,
1467 0, 1,
1468 2, 0,
1469 3, 0,
1470 1, 2,
1471 1, 3,
1472 3, 1,
1473 2, 1,
1474 0, 3,
1475 0, 2,
1476 2, 2,
1477 3, 2,
1478 3, 3,
1479 2, 3).finished();
1480// clang-format on
1481
1482} // namespace lf::uscalfe
1483
1484#endif
static constexpr RefEl kSegment()
Returns the (1-dimensional) reference segment.
Definition ref_el.h:153
const Eigen::MatrixXd & NodeCoords() const
Get the coordinates of the nodes of this reference element.
Definition ref_el.h:241
static constexpr RefEl kTria()
Returns the reference triangle.
Definition ref_el.h:161
constexpr size_type NumNodes() const
The number of nodes of this reference element.
Definition ref_el.h:213
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
Definition ref_el.h:169
Interface class for parametric scalar valued finite elements.
Linear Lagrange finite element on the quadrilateral reference element.
Definition lagr_fe.h:166
base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
Definition lagr_fe.h:175
FeLagrangeO1Quad(const FeLagrangeO1Quad &)=default
size_type NumRefShapeFunctions() const override
The local shape functions: four bilinear basis functions.
Definition lagr_fe.h:184
size_type NumRefShapeFunctions(dim_t codim) const override
Exactly one shape function attached to each node, none to other sub-entities.
Definition lagr_fe.h:192
size_type NumEvaluationNodes() const override
Tell the number of evaluation (interpolation) nodes.
Definition lagr_fe.h:260
FeLagrangeO1Quad(FeLagrangeO1Quad &&) noexcept=default
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t subidx) const override
Exactly one shape function attached to each node, none to other sub-entities.
Definition lagr_fe.h:202
Eigen::MatrixXd EvaluationNodes() const override
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom,...
Definition lagr_fe.h:256
unsigned Degree() const override
Request the maximal polynomial degree of the basis functions in this finite element.
Definition lagr_fe.h:179
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Computation of the gradients of all reference shape functions in a number of points.
Definition lagr_fe.h:232
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Evaluation of all reference shape functions in a number of points.
Definition lagr_fe.h:211
Linear Lagrange finite element on a line segment.
Definition lagr_fe.h:285
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t subidx) const override
All shape functions attached to nodes.
Definition lagr_fe.h:320
base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
Definition lagr_fe.h:294
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Evaluation of all reference shape functions in a number of points.
Definition lagr_fe.h:329
size_type NumEvaluationNodes() const override
Two evaluation nodes.
Definition lagr_fe.h:367
FeLagrangeO1Segment(const FeLagrangeO1Segment &)=default
size_type NumRefShapeFunctions(dim_t codim) const override
All shape functions attached to nodes.
Definition lagr_fe.h:310
unsigned Degree() const override
Request the maximal polynomial degree of the basis functions in this finite element.
Definition lagr_fe.h:298
Eigen::MatrixXd EvaluationNodes() const override
Evalutation nodes are just the vertices of the segment.
Definition lagr_fe.h:360
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Computation of the gradients of all reference shape functions in a number of points.
Definition lagr_fe.h:343
FeLagrangeO1Segment(FeLagrangeO1Segment &&) noexcept=default
size_type NumRefShapeFunctions() const override
The local shape functions: barycentric coordinate functions.
Definition lagr_fe.h:303
Linear Lagrange finite element on triangular reference element.
Definition lagr_fe.h:58
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t subidx) const override
Exactly one shape function attached to each node, none to other sub-entities.
Definition lagr_fe.h:94
size_type NumRefShapeFunctions(dim_t codim) const override
Exactly one shape function attached to each node, none to other sub-entities.
Definition lagr_fe.h:84
size_type NumEvaluationNodes() const override
Three evaluation nodes.
Definition lagr_fe.h:142
Eigen::MatrixXd EvaluationNodes() const override
Evalutation nodes are just the vertices of the triangle.
Definition lagr_fe.h:135
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Computation of the gradients of all reference shape functions in a number of points.
Definition lagr_fe.h:117
size_type NumRefShapeFunctions() const override
The local shape functions: barycentric coordinate functions.
Definition lagr_fe.h:76
unsigned Degree() const override
Request the maximal polynomial degree of the basis functions in this finite element.
Definition lagr_fe.h:71
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Evaluation of all reference shape functions in a number of points.
Definition lagr_fe.h:104
base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
Definition lagr_fe.h:67
FeLagrangeO1Tria(FeLagrangeO1Tria &&) noexcept=default
FeLagrangeO1Tria(const FeLagrangeO1Tria &)=default
Quadratic Lagrangian finite element on a quadrilateral reference element.
Definition lagr_fe.h:714
unsigned Degree() const override
Quadratic Lagrangian finite elements sport polynomials of degree 2.
Definition lagr_fe.h:730
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Point evaluation of reference local shape functions in unit square.
Definition lagr_fe.h:780
static const FeLagrangeO2Segment< SCALAR > krsf_segment_
Definition lagr_fe.h:886
FeLagrangeO2Quad(const FeLagrangeO2Quad &)=default
FeLagrangeO2Quad(FeLagrangeO2Quad &&) noexcept=default
size_type NumEvaluationNodes() const override
Nine evaluation nodes, same as the number of local shape functions.
Definition lagr_fe.h:878
Eigen::MatrixXd EvaluationNodes() const override
Evaluation nodes are the vertices, the midpoints of the edges and the center of the quadrilateral.
Definition lagr_fe.h:864
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t) const override
One shape function is attached to each node and each edge of the quadrilateral. One interior shape fu...
Definition lagr_fe.h:756
size_type NumRefShapeFunctions(dim_t codim) const override
One shape function is attached to each node and edge of the quadrilateral. One interior shape functio...
Definition lagr_fe.h:745
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Point evaluation of gradient of reference shape functions for quadratic Lagrangian finite elements on...
Definition lagr_fe.h:812
size_type NumRefShapeFunctions() const override
Nine local shape functions are associated with a quadrilateral cell.
Definition lagr_fe.h:738
base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
Definition lagr_fe.h:724
static const Eigen::Matrix< int, 9, 2 > ksegment_to_quad_mapping_
Definition lagr_fe.h:905
Quadratic Lagrangian finite element on a line segment.
Definition lagr_fe.h:564
size_type NumEvaluationNodes() const override
Three evaluation nodes, same number as local shape functions.
Definition lagr_fe.h:683
size_type NumRefShapeFunctions() const override
Three local shape functions are associated with an edge.
Definition lagr_fe.h:586
FeLagrangeO2Segment(FeLagrangeO2Segment &&) noexcept=default
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t) const override
One shape function attached to each node and one interior shape function.
Definition lagr_fe.h:606
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Evaluation of shape functions on reference segment (unit interval)
Definition lagr_fe.h:627
FeLagrangeO2Segment(const FeLagrangeO2Segment &)=default
base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
Definition lagr_fe.h:573
Eigen::MatrixXd EvaluationNodes() const override
Evaluation nodes are the endpoints of the segment and its midpoint.
Definition lagr_fe.h:673
unsigned Degree() const override
Quadratic Lagrangian finite element rely on polynomials of degree 2.
Definition lagr_fe.h:579
size_type NumRefShapeFunctions(dim_t codim) const override
One shape function attached to each node and one interior shape function.
Definition lagr_fe.h:596
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Evaluation of derivatives of local shape function on reference interval.
Definition lagr_fe.h:653
Quadratic Lagrangian finite element on a triangular reference element.
Definition lagr_fe.h:401
FeLagrangeO2Tria(const FeLagrangeO2Tria &)=default
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Point evaluation of reference shape functions.
Definition lagr_fe.h:459
size_type NumEvaluationNodes() const override
Six evaluation nodes, the same number as local shape functions.
Definition lagr_fe.h:543
Eigen::MatrixXd EvaluationNodes() const override
Evaluation nodes are the three vertices of the triangle and the three midpoints of its edges.
Definition lagr_fe.h:529
size_type NumRefShapeFunctions() const override
Six local shape functions are associated with a triangular cell in the case of quadratic Lagrangian f...
Definition lagr_fe.h:426
base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
Definition lagr_fe.h:412
FeLagrangeO2Tria(FeLagrangeO2Tria &&) noexcept=default
unsigned Degree() const override
Quadratic Lagrangian finite elements rely on polynomials of degree 2.
Definition lagr_fe.h:419
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t) const override
One shape function attached to each node and one to each edge of the triangle. There are no interior ...
Definition lagr_fe.h:445
size_type NumRefShapeFunctions(dim_t codim) const override
One shape function attached to each node and one to each edge of the triangle. There are no interior ...
Definition lagr_fe.h:434
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Point evaluations of gradients of reference shape functions.
Definition lagr_fe.h:486
Cubic Lagrangian finite element on a quadrilateral reference element.
Definition lagr_fe.h:1280
FeLagrangeO3Quad(const FeLagrangeO3Quad &)=default
Eigen::MatrixXd EvaluationNodes() const override
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom,...
Definition lagr_fe.h:1400
size_type NumEvaluationNodes() const override
Sixteen evaluation nodes.
Definition lagr_fe.h:1426
size_type NumRefShapeFunctions(dim_t codim) const override
One shape function attached to each node and two to each edge of the quadrilateral....
Definition lagr_fe.h:1303
static const FeLagrangeO3Segment< SCALAR > krsf_segment_
Definition lagr_fe.h:1434
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t subidx) const override
One shape function attached to each node and two to each edge of the quadrilateral....
Definition lagr_fe.h:1323
size_type NumRefShapeFunctions() const override
Total number of reference shape functions associated with the reference cell.
Definition lagr_fe.h:1295
unsigned Degree() const override
Request the maximal polynomial degree of the basis functions in this finite element.
Definition lagr_fe.h:1293
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Computation of the gradients of all reference shape functions in a number of points.
Definition lagr_fe.h:1356
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Evaluation of all reference shape functions in a number of points.
Definition lagr_fe.h:1332
base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
Definition lagr_fe.h:1289
static const Eigen::Matrix< int, 16, 2 > ksegment_to_quad_mapping_
Definition lagr_fe.h:1453
FeLagrangeO3Quad(FeLagrangeO3Quad &&) noexcept=default
Cubic Lagrangian finite element on a line segment.
Definition lagr_fe.h:1156
size_type NumRefShapeFunctions() const override
Total number of reference shape functions associated with the reference cell.
Definition lagr_fe.h:1171
unsigned Degree() const override
Request the maximal polynomial degree of the basis functions in this finite element.
Definition lagr_fe.h:1169
FeLagrangeO3Segment(const FeLagrangeO3Segment &)=default
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t subidx) const override
One shape function attached to each node and two interior shape functions.
Definition lagr_fe.h:1190
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Evaluation of all reference shape functions in a number of points.
Definition lagr_fe.h:1199
Eigen::MatrixXd EvaluationNodes() const override
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom,...
Definition lagr_fe.h:1240
FeLagrangeO3Segment(FeLagrangeO3Segment &&) noexcept=default
size_type NumEvaluationNodes() const override
Four evaluation nodes.
Definition lagr_fe.h:1249
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Computation of the gradients of all reference shape functions in a number of points.
Definition lagr_fe.h:1220
size_type NumRefShapeFunctions(dim_t codim) const override
One shape function attached to each node and two interior shape functions.
Definition lagr_fe.h:1179
base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
Definition lagr_fe.h:1165
Cubic Lagrangian finite elment on a triangular reference element.
Definition lagr_fe.h:946
FeLagrangeO3Tria(const FeLagrangeO3Tria &)=default
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t subidx) const override
One shape function attached to each node and two to each edge of the triangle. There is one interior ...
Definition lagr_fe.h:992
size_type NumRefShapeFunctions(dim_t codim) const override
One shape function attached to each node and two to each edge of the triangle. There is one interior ...
Definition lagr_fe.h:979
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Point evaluations of gradients of reference shape functions.
Definition lagr_fe.h:1053
FeLagrangeO3Tria(FeLagrangeO3Tria &&) noexcept=default
unsigned Degree() const override
Quadratic Lagrangian finite elements rely on polynomials of degree 3.
Definition lagr_fe.h:962
size_type NumEvaluationNodes() const override
Ten evaluation nodes.
Definition lagr_fe.h:1137
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Point evaluation of reference shape functions.
Definition lagr_fe.h:1010
base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
Definition lagr_fe.h:955
size_type NumRefShapeFunctions() const override
Ten local shape functions are associated with a triangular cell in the case of cubic Lagrangian finit...
Definition lagr_fe.h:969
Eigen::MatrixXd EvaluationNodes() const override
Evaluation nodes are the three vertices of the triangle, two equispaced interio nodes for each edge,...
Definition lagr_fe.h:1115
unsigned int sub_idx_t
type for local indices of sub-entities
Definition types.h:28
lf::base::glb_idx_t glb_idx_t
lf::base::size_type size_type
Eigen::Index ldof_idx_t
lf::base::dim_t dim_t
Eigen::Index gdof_idx_t
Collects data structures and algorithms designed for scalar finite element methods primarily meant fo...
lf::assemble::gdof_idx_t gdof_idx_t
Definition lagr_fe.h:26
lf::assemble::size_type size_type
Definition lagr_fe.h:30
lf::assemble::glb_idx_t glb_idx_t
Definition lagr_fe.h:34
lf::assemble::dim_t dim_t
Definition lagr_fe.h:32
lf::base::sub_idx_t sub_idx_t
Definition lagr_fe.h:36
lf::assemble::ldof_idx_t ldof_idx_t
Definition lagr_fe.h:28