LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
scalar_reference_finite_element.h
1#ifndef LF_FE_SCALAR_REFERENCE_FINITE_ELEMENT_H_
2#define LF_FE_SCALAR_REFERENCE_FINITE_ELEMENT_H_
3
12#include <lf/assemble/assemble.h>
13
14namespace lf::fe {
20using size_type = lf::assemble::size_type;
22using dim_t = lf::assemble::dim_t;
24using glb_idx_t = lf::assemble::glb_idx_t;
26using sub_idx_t = lf::base::sub_idx_t;
27
81template <typename SCALAR>
83 protected:
86 // NOLINTNEXTLINE
88 default;
89 // NOLINTNEXTLINE
91 default;
92 // NOLINTNEXTLINE
94 ScalarReferenceFiniteElement&&) noexcept = default;
95
96 public:
98 using Scalar = SCALAR;
99
101 virtual ~ScalarReferenceFiniteElement() = default;
102
107 [[nodiscard]] virtual base::RefEl RefEl() const = 0;
108
114 [[nodiscard]] virtual unsigned int Degree() const = 0;
115
119 [[nodiscard]] dim_t Dimension() const { return RefEl().Dimension(); }
120
128 [[nodiscard]] virtual size_type NumRefShapeFunctions() const {
129 size_type cnt = 0;
130 for (dim_t codim = 0; codim <= Dimension(); ++codim) {
131 for (sub_idx_t subidx = 0; subidx < RefEl().NumSubEntities(codim);
132 ++subidx) {
133 cnt += NumRefShapeFunctions(codim, subidx);
134 }
135 }
136 return cnt;
137 }
138
151 // NOLINTNEXTLINE(misc-unused-parameters)
152 [[nodiscard]] virtual size_type NumRefShapeFunctions(dim_t codim) const {
153 LF_VERIFY_MSG(false, "Illegal call for non-uniform sub-entity dof numbers");
154 return 0;
155 }
156
166 [[nodiscard]] virtual size_type NumRefShapeFunctions(
167 dim_t codim, sub_idx_t subidx) const = 0;
168
201 [[nodiscard]] virtual Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
202 EvalReferenceShapeFunctions(const Eigen::MatrixXd& refcoords) const = 0;
203
233 [[nodiscard]] virtual Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
234 GradientsReferenceShapeFunctions(const Eigen::MatrixXd& refcoords) const = 0;
235
278 [[nodiscard]] virtual Eigen::MatrixXd EvaluationNodes() const = 0;
279
283 [[nodiscard]] virtual size_type NumEvaluationNodes() const = 0;
284
310 [[nodiscard]] virtual Eigen::Matrix<SCALAR, 1, Eigen::Dynamic>
312 const Eigen::Matrix<SCALAR, 1, Eigen::Dynamic>& nodvals) const {
313 LF_ASSERT_MSG(nodvals.cols() == NumEvaluationNodes(),
314 "nodvals = " << nodvals << " <-> " << NumEvaluationNodes());
315 return nodvals;
316 }
317};
318
336template <class SCALAR>
337void PrintInfo(std::ostream& o,
339 unsigned int ctrl = 0) {
340 o << typeid(srfe).name() << ", degree = " << srfe.Degree()
341 << ", n_rsf = " << srfe.NumRefShapeFunctions()
342 << ", n_evln = " << srfe.NumEvaluationNodes();
343 if (ctrl > 0) {
344 o << "\n evl nodes = " << srfe.EvaluationNodes();
345 }
346}
347
352template <typename SCALAR>
353std::ostream& operator<<(std::ostream& o,
355 PrintInfo(o, fe_desc, 0);
356 return o;
357}
358
359} // end namespace lf::fe
360
362
366template <class SCALAR>
367struct fmt::formatter<lf::fe::ScalarReferenceFiniteElement<SCALAR>>
368 : ostream_formatter {};
370
371#endif // LF_FE_SCALAR_REFERENCE_FINITE_ELEMENT_H_
constexpr size_type NumSubEntities(dim_t sub_codim) const
Get the number of sub-entities of this RefEl with the given codimension.
Definition ref_el.h:308
constexpr dim_t Dimension() const
Return the dimension of this reference element.
Definition ref_el.h:204
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.
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.
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
Definition types.h:28
lf::base::glb_idx_t glb_idx_t
lf::base::size_type size_type
Eigen::Index ldof_idx_t
lf::base::dim_t dim_t
Eigen::Index gdof_idx_t
Collects data structures and algorithms designed for scalar finite element methods primarily meant fo...
Definition fe.h:47
lf::assemble::ldof_idx_t ldof_idx_t
Definition fe.h:52
lf::assemble::gdof_idx_t gdof_idx_t
Definition fe.h:50
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
Definition assemble.h:31