10#ifndef INCGe281cd0ab7fb476e9315a3dda7f45ffe
11#define INCGe281cd0ab7fb476e9315a3dda7f45ffe
13#include <Eigen/src/Core/util/ForwardDeclarations.h>
14#include <lf/quad/quad.h>
18#include "uniform_scalar_fe_space.h"
43template <
class SCALAR>
88 LF_ASSERT_MSG(
fe_ !=
nullptr,
"Not initialized.");
92 [[nodiscard]]
unsigned Degree()
const override {
93 LF_ASSERT_MSG(
fe_ !=
nullptr,
"Not initialized.");
98 LF_ASSERT_MSG(
fe_ !=
nullptr,
"Not initialized.");
99 return fe_->NumRefShapeFunctions();
103 LF_ASSERT_MSG(
fe_ !=
nullptr,
"Not initialized.");
104 return fe_->NumRefShapeFunctions(codim);
109 LF_ASSERT_MSG(
fe_ !=
nullptr,
"Not initialized.");
110 return fe_->NumRefShapeFunctions(codim, subidx);
113 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
115 LF_ASSERT_MSG(
fe_ !=
nullptr,
"Not initialized.");
116 return fe_->EvalReferenceShapeFunctions(local);
118 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
120 const Eigen::MatrixXd& local)
const override {
121 LF_ASSERT_MSG(
fe_ !=
nullptr,
"Not initialized.");
122 return fe_->GradientsReferenceShapeFunctions(local);
125 LF_ASSERT_MSG(
fe_ !=
nullptr,
"Not initialized.");
126 return fe_->EvaluationNodes();
129 LF_ASSERT_MSG(
fe_ !=
nullptr,
"Not initialized.");
130 return fe_->NumEvaluationNodes();
134 const Eigen::Matrix<SCALAR, 1, Eigen::Dynamic>& nodvals)
const override {
135 LF_ASSERT_MSG(
fe_ !=
nullptr,
"Not initialized.");
136 return fe_->NodalValuesToDofs(nodvals);
144 LF_ASSERT_MSG(
fe_ !=
nullptr,
"Not initialized.");
151 [[nodiscard]]
const Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>&
153 LF_ASSERT_MSG(
fe_ !=
nullptr,
"Not initialized.");
165 [[nodiscard]]
const Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>&
167 LF_ASSERT_MSG(
fe_ !=
nullptr,
"Not initialized.");
181 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
shap_fun_;
Represents a reference element with all its properties.
Interface class for parametric scalar valued finite elements.
ScalarReferenceFiniteElement()=default
Represents a Quadrature Rule over one of the Reference Elements.
bool isInitialized() const
Tells initialization status of object.
PrecomputedScalarReferenceFiniteElement()=default
Default constructor which does not initialize this class at all (invalid state). If any method is cal...
PrecomputedScalarReferenceFiniteElement(PrecomputedScalarReferenceFiniteElement &&) noexcept=default
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.
PrecomputedScalarReferenceFiniteElement(const PrecomputedScalarReferenceFiniteElement &)=delete
fe::ScalarReferenceFiniteElement< SCALAR > const * fe_
~PrecomputedScalarReferenceFiniteElement() override=default
base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
size_type NumRefShapeFunctions(dim_t codim) const override
The number of interior reference shape functions for sub-entities of a particular co-dimension.
const Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > & PrecompReferenceShapeFunctions() const
Value of EvalReferenceShapeFunctions(Qr().Points())
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.
const quad::QuadRule & Qr() const
Return the Quadrature rule at which the shape functions (and their gradients) have been precomputed.
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...
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > grad_shape_fun_
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.
const Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > & PrecompGradientsReferenceShapeFunctions() const
Value of EvalGradientsReferenceShapeFunctions(Qr().Points())
Eigen::MatrixXd EvaluationNodes() const override
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom,...
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > shap_fun_
unsigned Degree() const override
Request the maximal polynomial degree of the basis functions in this finite element.
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t subidx) const override
The number of interior reference shape functions for every sub-entity.
Collects data structures and algorithms designed for scalar finite element methods primarily meant fo...
Rules for numerical quadrature on reference entity shapes.
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