22 #ifndef linear_scalar_shapeset_barycentric_hpp
23 #define linear_scalar_shapeset_barycentric_hpp
26 #include "basis_data.hpp"
27 #include "piecewise_linear_continuous_scalar_basis.hpp"
32 template <
typename ValueType>
36 typedef typename Basis<ValueType>::CoordinateType CoordinateType;
37 enum BasisType {TYPE1,TYPE2};
55 const arma::Mat<CoordinateType>& points,
56 LocalDofIndex localDofIndex,
61 linearBasis.
evaluate(what,points,ALL_DOFS,temp);
63 ValueType coeffs[2][3][3] =
65 {{1,1./3,1./2},{0,1./3,0},{0,1./3,1./2}},
66 {{1,1./2,1./3},{0,1./2,1./3},{0,0,1./3}}
75 if (localDofIndex != ALL_DOFS){
79 for (
int i=0;i<data.
values.extent(2);++i){
82 data.
values(0,0,i) +=coeffs[type_index][localDofIndex][j]*
86 if (what & DERIVATIVES){
93 coeffs[type_index][localDofIndex][k]*
102 for (
int dofIndex=0;dofIndex<3;++dofIndex){
103 for (
int i=0;i<data.
values.extent(2);++i){
104 data.
values(0,dofIndex,i) = 0;
105 for (
int j=0;j<3;++j)
106 data.
values(0,dofIndex,i) +=coeffs[type_index][dofIndex][j]*
113 if (what & DERIVATIVES){
115 for (
int dofIndex=0;dofIndex<3;++dofIndex)
119 for (
int k=0;k<3;++k)
121 coeffs[type_index][dofIndex][k]*
133 virtual std::pair<const char*,int>
clCodeString (
bool isTestBasis)
const {
134 throw std::runtime_error(
"PiecewiseLinearContinuousScalarBasisBarycentric::clCodeString():"
135 "OpenCL not supported for this basis type.");
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: linear_scalar_shapeset_barycentric.hpp:54
Definition: linear_scalar_shapeset_barycentric.hpp:33
_4dArray< ValueType > derivatives
Derivatives of shape functions.
Definition: basis_data.hpp:59
virtual int size() const
Return the number of shape functions.
Definition: linear_scalar_shapeset_barycentric.hpp:46
virtual int order() const
Return the maximum polynomial order of shape functions.
Definition: linear_scalar_shapeset_barycentric.hpp:50
_3dArray< ValueType > values
Values of shape functions.
Definition: basis_data.hpp:51
virtual std::pair< const char *, int > clCodeString(bool isTestBasis) const
Returns an OpenCL code snippet for shape function evaluation.
Definition: linear_scalar_shapeset_barycentric.hpp:133
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
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: linear_scalar_shapeset.hpp:88