1#ifndef LF_FE_SCALAR_REFERENCE_FINITE_ELEMENT_H_
2#define LF_FE_SCALAR_REFERENCE_FINITE_ELEMENT_H_
12#include <lf/assemble/assemble.h>
81template <
typename SCALAR>
114 [[nodiscard]] virtual
unsigned int Degree() const = 0;
153 LF_VERIFY_MSG(
false,
"Illegal call for non-uniform sub-entity dof numbers");
201 [[nodiscard]]
virtual Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
233 [[nodiscard]]
virtual Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
310 [[nodiscard]]
virtual Eigen::Matrix<SCALAR, 1, Eigen::Dynamic>
312 const Eigen::Matrix<SCALAR, 1, Eigen::Dynamic>& nodvals)
const {
336template <
class SCALAR>
339 unsigned int ctrl = 0) {
340 o <<
typeid(srfe).name() <<
", degree = " << srfe.
Degree()
352template <
typename SCALAR>
355 PrintInfo(o, fe_desc, 0);
366template <
class SCALAR>
368 : ostream_formatter {};
constexpr size_type NumSubEntities(dim_t sub_codim) const
Get the number of sub-entities of this RefEl with the given codimension.
constexpr dim_t Dimension() const
Return the dimension of this reference element.
Interface class for parametric scalar valued finite elements.
ScalarReferenceFiniteElement(const ScalarReferenceFiniteElement &)=default
virtual Eigen::MatrixXd EvaluationNodes() const =0
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom,...
ScalarReferenceFiniteElement(ScalarReferenceFiniteElement &&) noexcept=default
virtual Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const =0
Evaluation of all reference shape functions in a number of points.
virtual size_type NumEvaluationNodes() const =0
Tell the number of evaluation (interpolation) nodes.
dim_t Dimension() const
Returns the spatial dimension of the reference cell.
virtual Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > NodalValuesToDofs(const Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > &nodvals) const
Computes the linear combination of reference shape functions matching function values at evaluation n...
virtual size_type NumRefShapeFunctions(dim_t codim) const
The number of interior reference shape functions for sub-entities of a particular co-dimension.
virtual unsigned int Degree() const =0
Request the maximal polynomial degree of the basis functions in this finite element.
virtual Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const =0
Computation of the gradients of all reference shape functions in a number of points.
std::ostream & operator<<(std::ostream &o, const ScalarReferenceFiniteElement< SCALAR > &fe_desc)
Stream output operator: just calls the ScalarReferenceFiniteElement::print() method.
virtual size_type NumRefShapeFunctions(dim_t codim, sub_idx_t subidx) const =0
The number of interior reference shape functions for every sub-entity.
SCALAR Scalar
The underlying scalar type.
void PrintInfo(std::ostream &o, const ScalarReferenceFiniteElement< SCALAR > &srfe, unsigned int ctrl=0)
Print information about a ScalarReferenceFiniteElement to the given stream.
virtual size_type NumRefShapeFunctions() const
Total number of reference shape functions associated with the reference cell.
ScalarReferenceFiniteElement()=default
virtual base::RefEl RefEl() const =0
Tells the type of reference cell underlying the parametric finite element.
unsigned int sub_idx_t
type for local indices of sub-entities
lf::base::glb_idx_t glb_idx_t
lf::base::size_type size_type
Collects data structures and algorithms designed for scalar finite element methods primarily meant fo...
lf::assemble::ldof_idx_t ldof_idx_t
lf::assemble::gdof_idx_t gdof_idx_t
lf::assemble::size_type size_type
lf::assemble::dim_t dim_t
lf::base::sub_idx_t sub_idx_t