1#ifndef PROJECTS_DPG_DISCONTINUOUS_SCALAR_REFERENCE_FINITE_ELEMENT
2#define PROJECTS_DPG_DISCONTINUOUS_SCALAR_REFERENCE_FINITE_ELEMENT
14#include <lf/uscalfe/uscalfe.h>
41template <
typename SCALAR>
76 [[nodiscard]]
unsigned Degree()
const override {
78 return cfe_->Degree();
83 return cfe_->NumRefShapeFunctions();
90 return (codim == 0) ?
cfe_->NumRefShapeFunctions() : 0;
98 return (codim == 0) ?
cfe_->NumRefShapeFunctions() : 0;
101 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
104 return cfe_->EvalReferenceShapeFunctions(local);
107 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
109 const Eigen::MatrixXd& local)
const override {
111 return cfe_->GradientsReferenceShapeFunctions(local);
116 return cfe_->EvaluationNodes();
121 return cfe_->NumEvaluationNodes();
125 const Eigen::Matrix<SCALAR, 1, Eigen::Dynamic>& nodvals)
const override {
127 return cfe_->NodalValuesToDofs(nodvals);
135 std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR>>
cfe_;
Represents a reference element with all its properties.
Interface class for parametric scalar valued finite elements.
ScalarReferenceFiniteElement()=default
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()=default
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t) const override
DiscontinuousScalarReferenceFiniteElement(DiscontinuousScalarReferenceFiniteElement &&) noexcept=default
~DiscontinuousScalarReferenceFiniteElement() override=default
size_type NumRefShapeFunctions(dim_t codim) const override
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...
Contains functionality for the implementation of DPG methods.
lf::uscalfe::sub_idx_t sub_idx_t
Type for indexing of sub-entities.
lf::uscalfe::dim_t dim_t
Tpe for (co)-dimensions.
lf::uscalfe::size_type size_type
Type for vector length/matrix sizes.