10#ifndef INCGb9e80d63ee55424493f538a7971df592
11#define INCGb9e80d63ee55424493f538a7971df592
14#include <lf/uscalfe/uscalfe.h>
22template <
class SCALAR>
31 [[nodiscard]]
unsigned Degree()
const override {
return inner_->Degree(); }
34 return inner_->NumRefShapeFunctions(codim, subidx);
36 [[nodiscard]] Eigen::Matrix<std::complex<SCALAR>, Eigen::Dynamic,
39 return std::complex<SCALAR>(0, 1) *
40 inner_->EvalReferenceShapeFunctions(refcoords);
42 [[nodiscard]] Eigen::Matrix<std::complex<SCALAR>, Eigen::Dynamic,
45 const Eigen::MatrixXd &refcoords)
const override {
46 return std::complex<SCALAR>(0, 1) *
47 inner_->GradientsReferenceShapeFunctions(refcoords);
50 return inner_->EvaluationNodes();
53 return inner_->NumEvaluationNodes();
56 [[nodiscard]] Eigen::Matrix<std::complex<SCALAR>, 1, Eigen::Dynamic>
58 &nodvals)
const override {
59 return inner_->NodalValuesToDofs(
60 (nodvals / std::complex<SCALAR>(0, 1)).real());
64 return inner_->NumRefShapeFunctions(dim);
68 std::unique_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR>>
inner_;
75inline std::shared_ptr<uscalfe::UniformScalarFESpace<std::complex<double>>>
77 return std::make_shared<uscalfe::UniformScalarFESpace<std::complex<double>>>(
79 std::make_shared<ComplexScalarReferenceFiniteElement<double>>(
80 std::make_unique<uscalfe::FeLagrangeO1Tria<double>>()),
81 std::make_shared<ComplexScalarReferenceFiniteElement<double>>(
82 std::make_unique<uscalfe::FeLagrangeO1Quad<double>>()),
83 std::make_shared<ComplexScalarReferenceFiniteElement<double>>(
84 std::make_unique<uscalfe::FeLagrangeO1Segment<double>>()),
Represents a reference element with all its properties.
constexpr RefEl(RefElType type) noexcept
Create a RefEl from a lf::base::RefElType enum.
Interface class for parametric scalar valued finite elements.
Wraps another ScalarReferenceFiniteElement and multiplies the shape functions with the imaginary unit...
ComplexScalarReferenceFiniteElement(std::unique_ptr< lf::fe::ScalarReferenceFiniteElement< SCALAR > > fe)
size_type NumRefShapeFunctions(dim_t codim, sub_idx_t subidx) const override
The number of interior reference shape functions for every sub-entity.
size_type NumEvaluationNodes() const override
Tell the number of evaluation (interpolation) nodes.
std::unique_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > inner_
unsigned Degree() const override
Request the maximal polynomial degree of the basis functions in this finite element.
Eigen::Matrix< std::complex< SCALAR >, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Computation of the gradients of all reference shape functions in a number of points.
base::RefEl RefEl() const override
Tells the type of reference cell underlying the parametric finite element.
size_type NumRefShapeFunctions(dim_t dim) const override
The number of interior reference shape functions for sub-entities of a particular co-dimension.
Eigen::MatrixXd EvaluationNodes() const override
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom,...
Eigen::Matrix< std::complex< SCALAR >, 1, Eigen::Dynamic > NodalValuesToDofs(const Eigen::Matrix< std::complex< SCALAR >, 1, Eigen::Dynamic > &nodvals) const override
Eigen::Matrix< std::complex< SCALAR >, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Evaluation of all reference shape functions in a number of points.
Includes utilities to test classes from lf::fe.
std::shared_ptr< uscalfe::UniformScalarFESpace< std::complex< double > > > MakeComplexLagrangeO1FeSpace(std::shared_ptr< const mesh::Mesh > mesh_p)
Returns a UniformScalarFESpace that is made up of "complexified" (via ComplexScalarReferenceFiniteEle...
lf::assemble::size_type size_type
lf::assemble::dim_t dim_t
lf::base::sub_idx_t sub_idx_t