LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
fe_point.h
1#ifndef LF_FE_FE_POINT_H_
2#define LF_FE_FE_POINT_H_
3
12#include "scalar_reference_finite_element.h"
13
14namespace lf::fe {
15
24template <class SCALAR>
25class FePoint : public ScalarReferenceFiniteElement<SCALAR> {
26 public:
32 explicit FePoint(unsigned degree) : degree_(degree) {}
33
34 [[nodiscard]] base::RefEl RefEl() const override {
35 return base::RefEl::kPoint();
36 }
37 [[nodiscard]] unsigned Degree() const override { return degree_; }
39 dim_t codim, sub_idx_t subidx) const override {
40 LF_ASSERT_MSG(codim == 0, "Codim out of bounds");
41 LF_ASSERT_MSG(subidx == 0, "subidx out of bounds.");
42 return 1;
43 }
44 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
45 EvalReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override {
46 LF_ASSERT_MSG(refcoords.rows() == 0, "refcoords has too many rows.");
47 return Eigen::Matrix<SCALAR, 1, Eigen::Dynamic>::Ones(1, refcoords.cols());
48 }
49 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
51 const Eigen::MatrixXd & /*refcoords*/) const override {
52 LF_VERIFY_MSG(false, "gradients not defined in points of mesh.");
53 }
54 [[nodiscard]] Eigen::MatrixXd EvaluationNodes() const override {
55 return {0, 1};
56 }
57 [[nodiscard]] size_type NumEvaluationNodes() const override { return 1; }
58
59 private:
60 unsigned degree_;
61};
62
63} // end namespace lf::fe
64
65#endif // LF_FE_FE_POINT_H_
Represents a reference element with all its properties.
Definition ref_el.h:109
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
Definition ref_el.h:144
Finite element on a point.
Definition fe_point.h:25
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t subidx) const override
The number of interior reference shape functions for every sub-entity.
Definition fe_point.h:38
unsigned Degree() const override
Request the maximal polynomial degree of the basis functions in this finite element.
Definition fe_point.h:37
unsigned degree_
Definition fe_point.h:60
base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
Definition fe_point.h:34
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.
Definition fe_point.h:50
Eigen::MatrixXd EvaluationNodes() const override
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom,...
Definition fe_point.h:54
size_type NumEvaluationNodes() const override
Tell the number of evaluation (interpolation) nodes.
Definition fe_point.h:57
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 fe_point.h:45
FePoint(unsigned degree)
Create a new FePoint by specifying the degree of the shape functions.
Definition fe_point.h:32
Interface class for parametric scalar valued finite elements.
Collects data structures and algorithms designed for scalar finite element methods primarily meant fo...
Definition fe.h:47
lf::assemble::size_type size_type
Definition fe.h:54
lf::assemble::dim_t dim_t
Definition fe.h:56
lf::base::sub_idx_t sub_idx_t
Definition fe.h:60