LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
discontinuous_scalar_reference_finite_element.h
1#ifndef PROJECTS_DPG_DISCONTINUOUS_SCALAR_REFERENCE_FINITE_ELEMENT
2#define PROJECTS_DPG_DISCONTINUOUS_SCALAR_REFERENCE_FINITE_ELEMENT
3
13#include <lf/fe/fe.h>
14#include <lf/uscalfe/uscalfe.h>
15
16#include <iostream>
17#include <typeinfo>
18
19#include "dpg.h"
20
21namespace projects::dpg {
22
41template <typename SCALAR>
44 public:
50
59
61 std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR>> cfe)
62 : lf::fe::ScalarReferenceFiniteElement<SCALAR>(), cfe_(std::move(cfe)) {}
63
69 [[nodiscard]] inline bool isInitialized() const { return (cfe_ != nullptr); }
70
71 [[nodiscard]] lf::base::RefEl RefEl() const override {
72 LF_ASSERT_MSG(isInitialized(), "Not initialized!");
73 return cfe_->RefEl();
74 }
75
76 [[nodiscard]] unsigned Degree() const override {
77 LF_ASSERT_MSG(isInitialized(), "Not Initialized");
78 return cfe_->Degree();
79 }
80
81 [[nodiscard]] size_type NumRefShapeFunctions() const override {
82 LF_ASSERT_MSG(isInitialized(), "Not initialized");
83 return cfe_->NumRefShapeFunctions();
84 }
85
88 [[nodiscard]] size_type NumRefShapeFunctions(dim_t codim) const override {
89 LF_ASSERT_MSG(isInitialized(), "Not initialized");
90 return (codim == 0) ? cfe_->NumRefShapeFunctions() : 0;
91 }
92
96 dim_t codim, sub_idx_t /*subidx*/) const override {
97 LF_ASSERT_MSG(isInitialized(), "Not initialized");
98 return (codim == 0) ? cfe_->NumRefShapeFunctions() : 0;
99 }
100
101 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
102 EvalReferenceShapeFunctions(const Eigen::MatrixXd& local) const override {
103 LF_ASSERT_MSG(isInitialized(), "Not initialized");
104 return cfe_->EvalReferenceShapeFunctions(local);
105 }
106
107 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
109 const Eigen::MatrixXd& local) const override {
110 LF_ASSERT_MSG(isInitialized(), "Not initialized");
111 return cfe_->GradientsReferenceShapeFunctions(local);
112 }
113
114 [[nodiscard]] Eigen::MatrixXd EvaluationNodes() const override {
115 LF_ASSERT_MSG(isInitialized(), "Not initialized");
116 return cfe_->EvaluationNodes();
117 }
118
119 [[nodiscard]] size_type NumEvaluationNodes() const override {
120 LF_ASSERT_MSG(isInitialized(), "Not initialized");
121 return cfe_->NumEvaluationNodes();
122 }
123
124 [[nodiscard]] Eigen::Matrix<SCALAR, 1, Eigen::Dynamic> NodalValuesToDofs(
125 const Eigen::Matrix<SCALAR, 1, Eigen::Dynamic>& nodvals) const override {
126 LF_ASSERT_MSG(isInitialized(), "Not initialized");
127 return cfe_->NodalValuesToDofs(nodvals);
128 }
129
132
133 private:
135 std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR>> cfe_;
136};
137
138} // namespace projects::dpg
139
140#endif // PROJECTS_DPG_DISCONTINUOUS_SCALAR_REFERENCE_FINITE_ELEMENT
Represents a reference element with all its properties.
Definition: ref_el.h:106
Interface class for parametric scalar valued finite elements.
Decorator class around a ScalarReferenceFiniteElement to represent discontinuous shape functions.
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions(const Eigen::MatrixXd &local) const override
Evaluation of all reference shape functions in a number of points.
lf::base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > cfe_
DiscontinuousScalarReferenceFiniteElement(const DiscontinuousScalarReferenceFiniteElement &)=delete
DiscontinuousScalarReferenceFiniteElement(DiscontinuousScalarReferenceFiniteElement &&) noexcept=default
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions(const Eigen::MatrixXd &local) const override
Computation of the gradients of all reference shape functions in a number of points.
bool isInitialized() const
Reports initialization status of the object.
Eigen::MatrixXd EvaluationNodes() const override
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom,...
unsigned Degree() const override
Request the maximal polynomial degree of the basis functions in this finite element.
size_type NumEvaluationNodes() const override
Tell the number of evaluation (interpolation) nodes.
size_type NumRefShapeFunctions() const override
Total number of reference shape functions associated with the reference cell.
Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > NodalValuesToDofs(const Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > &nodvals) const override
Computes the linear combination of reference shape functions matching function values at evaluation n...
Definition: assemble.h:30
Contains functionality for the implementation of DPG methods.
Definition: primal_dpg.h:33
lf::uscalfe::sub_idx_t sub_idx_t
Type for indexing of sub-entities.
Definition: dpg.h:26
lf::uscalfe::dim_t dim_t
Tpe for (co)-dimensions.
Definition: dpg.h:20
lf::uscalfe::size_type size_type
Type for vector length/matrix sizes.
Definition: dpg.h:17