1#ifndef LF_FE_FE_POINT_H_
2#define LF_FE_FE_POINT_H_
12#include "scalar_reference_finite_element.h"
24template <
class SCALAR>
40 LF_ASSERT_MSG(codim == 0,
"Codim out of bounds");
41 LF_ASSERT_MSG(subidx == 0,
"subidx out of bounds.");
44 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
46 LF_ASSERT_MSG(refcoords.rows() == 0,
"refcoords has too many rows.");
47 return Eigen::Matrix<SCALAR, 1, Eigen::Dynamic>::Ones(1, refcoords.cols());
49 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
51 const Eigen::MatrixXd & )
const override {
52 LF_VERIFY_MSG(
false,
"gradients not defined in points of mesh.");
Represents a reference element with all its properties.
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
Finite element on a point.
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t subidx) const override
The number of interior reference shape functions for every sub-entity.
unsigned Degree() const override
Request the maximal polynomial degree of the basis functions in this finite element.
base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions(const Eigen::MatrixXd &) const override
Computation of the gradients of all reference shape functions in a number of points.
Eigen::MatrixXd EvaluationNodes() const override
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom,...
size_type NumEvaluationNodes() const override
Tell the number of evaluation (interpolation) nodes.
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.
FePoint(unsigned degree)
Create a new FePoint by specifying the degree of the shape functions.
Interface class for parametric scalar valued finite elements.
Collects data structures and algorithms designed for scalar finite element methods primarily meant fo...
lf::assemble::size_type size_type
lf::assemble::dim_t dim_t
lf::base::sub_idx_t sub_idx_t