21 #ifndef fiber_lagrange_scalar_shapeset_hpp
22 #define fiber_lagrange_scalar_shapeset_hpp
24 #include "../common/common.hpp"
28 #include "basis_data.hpp"
29 #include "dune_basis_helper.hpp"
31 #include <dune/localfunctions/lagrange/pk2d/pk2dlocalbasis.hh>
39 template <
int elementVertexCount,
typename CoordinateType,
typename ValueType,
46 template <
typename CoordinateType,
typename ValueType,
int polynomialOrder>
50 typedef Dune::Pk2DLocalBasis<CoordinateType, ValueType, polynomialOrder>
52 typedef DuneShapeset DuneBasis;
56 template <
int elementVertexCount,
typename ValueType,
int polynomialOrder>
60 typedef typename Basis<ValueType>::CoordinateType CoordinateType;
64 elementVertexCount, CoordinateType, ValueType, polynomialOrder>::DuneShapeset
66 DuneBasis m_duneBasis;
70 return m_duneBasis.size();
74 return polynomialOrder;
78 const arma::Mat<CoordinateType>& points,
79 LocalDofIndex localDofIndex,
81 if (localDofIndex != ALL_DOFS &&
82 (localDofIndex < 0 ||
size() <= localDofIndex))
83 throw std::invalid_argument(
"LagrangeScalarShapeset::"
84 "evaluate(): Invalid localDofIndex");
87 evaluateShapeFunctionsWithDune<CoordinateType, ValueType, DuneBasis>(
88 points, localDofIndex, data.
values, m_duneBasis);
89 if (what & DERIVATIVES)
90 evaluateShapeFunctionDerivativesWithDune<
91 CoordinateType, ValueType, DuneBasis>(
92 points, localDofIndex, data.
derivatives, m_duneBasis);
_4dArray< ValueType > derivatives
Derivatives of shape functions.
Definition: basis_data.hpp:59
virtual int order() const
Return the maximum polynomial order of shape functions.
Definition: lagrange_scalar_shapeset.hpp:73
Definition: lagrange_scalar_shapeset.hpp:41
_3dArray< ValueType > values
Values of shape functions.
Definition: basis_data.hpp:51
virtual int size() const
Return the number of shape functions.
Definition: lagrange_scalar_shapeset.hpp:69
Collection of shape functions defined on a reference element.
Definition: basis.hpp:36
Storage of values and/or derivatives of shape functions.
Definition: basis_data.hpp:44
Shapeset composed of the Lagrange polynomials up to a specified order.
Definition: lagrange_scalar_shapeset.hpp:57
virtual void evaluate(size_t what, const arma::Mat< CoordinateType > &points, LocalDofIndex localDofIndex, BasisData< ValueType > &data) const
Evaluate the shape functions making up this shapeset and/or their derivatives at specified points...
Definition: lagrange_scalar_shapeset.hpp:77