21 #ifndef fiber_constant_scalar_shapeset_hpp
22 #define fiber_constant_scalar_shapeset_hpp
24 #include "../common/common.hpp"
27 #include "basis_data.hpp"
35 template <
typename ValueType>
39 typedef typename Basis<ValueType>::CoordinateType CoordinateType;
50 const arma::Mat<CoordinateType>& points,
51 LocalDofIndex localDofIndex,
53 if (localDofIndex != ALL_DOFS && localDofIndex != 0)
54 throw std::invalid_argument(
"ConstantScalarShapeset::evaluate(): "
55 "Invalid localDofIndex");
59 const int componentCount = 1;
60 const int functionCount = 1;
61 const int pointCount = points.n_cols;
64 data.
values.set_size(componentCount, functionCount, pointCount);
67 if (what & DERIVATIVES)
69 const int coordCount = points.n_rows;
70 data.
derivatives.set_size(componentCount, coordCount,
71 functionCount, pointCount);
_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: constant_scalar_shapeset.hpp:45
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: constant_scalar_shapeset.hpp:49
virtual int size() const
Return the number of shape functions.
Definition: constant_scalar_shapeset.hpp:41
_3dArray< ValueType > values
Values of shape functions.
Definition: basis_data.hpp:51
A shapeset containing only one function: the constant function.
Definition: constant_scalar_shapeset.hpp:36
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