LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
precomputed_scalar_reference_finite_element.h
1
10#ifndef INCGe281cd0ab7fb476e9315a3dda7f45ffe
11#define INCGe281cd0ab7fb476e9315a3dda7f45ffe
12
13#include <Eigen/src/Core/util/ForwardDeclarations.h>
14#include <lf/quad/quad.h>
15
16#include <memory>
17
18#include "uniform_scalar_fe_space.h"
19
20namespace lf::uscalfe {
21
43template <class SCALAR>
46 public:
52
55
57 PrecomputedScalarReferenceFiniteElement&&) noexcept = default;
58
61
63 PrecomputedScalarReferenceFiniteElement&&) noexcept = default;
73 lf::fe::ScalarReferenceFiniteElement<SCALAR> const* fe, quad::QuadRule qr)
74 : lf::fe::ScalarReferenceFiniteElement<SCALAR>(),
75 fe_(fe),
76 qr_(std::move(qr)),
79
85 [[nodiscard]] bool isInitialized() const { return (fe_ != nullptr); }
86
87 [[nodiscard]] base::RefEl RefEl() const override {
88 LF_ASSERT_MSG(fe_ != nullptr, "Not initialized.");
89 return fe_->RefEl();
90 }
91
92 [[nodiscard]] unsigned Degree() const override {
93 LF_ASSERT_MSG(fe_ != nullptr, "Not initialized.");
94 return fe_->Degree();
95 }
96
97 [[nodiscard]] size_type NumRefShapeFunctions() const override {
98 LF_ASSERT_MSG(fe_ != nullptr, "Not initialized.");
99 return fe_->NumRefShapeFunctions();
100 }
101
102 [[nodiscard]] size_type NumRefShapeFunctions(dim_t codim) const override {
103 LF_ASSERT_MSG(fe_ != nullptr, "Not initialized.");
104 return fe_->NumRefShapeFunctions(codim);
105 }
106
108 dim_t codim, sub_idx_t subidx) const override {
109 LF_ASSERT_MSG(fe_ != nullptr, "Not initialized.");
110 return fe_->NumRefShapeFunctions(codim, subidx);
111 }
112
113 [[nodiscard]] Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>
114 EvalReferenceShapeFunctions(const Eigen::MatrixXd& local) const override {
115 LF_ASSERT_MSG(fe_ != nullptr, "Not initialized.");
116 return fe_->EvalReferenceShapeFunctions(local);
117 }
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);
123 }
124 [[nodiscard]] Eigen::MatrixXd EvaluationNodes() const override {
125 LF_ASSERT_MSG(fe_ != nullptr, "Not initialized.");
126 return fe_->EvaluationNodes();
127 }
128 [[nodiscard]] size_type NumEvaluationNodes() const override {
129 LF_ASSERT_MSG(fe_ != nullptr, "Not initialized.");
130 return fe_->NumEvaluationNodes();
131 }
132
133 [[nodiscard]] Eigen::Matrix<SCALAR, 1, Eigen::Dynamic> NodalValuesToDofs(
134 const Eigen::Matrix<SCALAR, 1, Eigen::Dynamic>& nodvals) const override {
135 LF_ASSERT_MSG(fe_ != nullptr, "Not initialized.");
136 return fe_->NodalValuesToDofs(nodvals);
137 }
138
143 [[nodiscard]] const quad::QuadRule& Qr() const {
144 LF_ASSERT_MSG(fe_ != nullptr, "Not initialized.");
145 return qr_;
146 }
147
151 [[nodiscard]] const Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>&
153 LF_ASSERT_MSG(fe_ != nullptr, "Not initialized.");
154 return shap_fun_;
155 }
156
157 // clang-format off
164 // clang-format on
165 [[nodiscard]] const Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic>&
167 LF_ASSERT_MSG(fe_ != nullptr, "Not initialized.");
168 return grad_shape_fun_;
169 }
170
173
174 private:
181 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> shap_fun_;
183 Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> grad_shape_fun_;
184};
185
186} // namespace lf::uscalfe
187
188#endif // INCGe281cd0ab7fb476e9315a3dda7f45ffe
Represents a reference element with all its properties.
Definition ref_el.h:109
Interface class for parametric scalar valued finite elements.
Represents a Quadrature Rule over one of the Reference Elements.
Definition quad_rule.h:58
Helper class which stores a ScalarReferenceFiniteElement with precomputed values at the nodes of a qu...
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
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...
lf::assemble::size_type size_type
Definition lagr_fe.h:30
lf::assemble::dim_t dim_t
Definition lagr_fe.h:32
lf::base::sub_idx_t sub_idx_t
Definition lagr_fe.h:36
Definition assemble.h:31