3 #ifndef DUNE_P1_LOCALBASIS_HH 4 #define DUNE_P1_LOCALBASIS_HH 9 #include <dune/common/fmatrix.hh> 26 template<
class D,
class R,
int dim>
42 std::vector<typename Traits::RangeType>& out)
const 46 for (
size_t i=0; i<dim; i++) {
55 std::vector<typename Traits::JacobianType>& out)
const 59 for (
int i=0; i<dim; i++)
62 for (
int i=0; i<dim; i++)
63 for (
int j=0; j<dim; j++)
64 out[i+1][0][j] = (i==j);
75 std::vector<typename Traits::RangeType>& out)
const 77 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
81 else if (totalOrder==1)
83 auto direction = std::find(order.begin(), order.end(), 1);
87 for (
int i=0; i<dim; i++)
88 out[i+1] = (i==(direction-order.begin()));
94 for (
int i=0; i<dim+1; i++)
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:31
D DomainType
domain type
Definition: localbasis.hh:43
Linear Lagrange shape functions on the simplex.
Definition: p1localbasis.hh:27
LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim > > Traits
export type traits for function signature
Definition: p1localbasis.hh:32
unsigned int size() const
number of shape functions
Definition: p1localbasis.hh:35
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: p1localbasis.hh:54
void partial(const std::array< unsigned int, dim > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of any order of all shape functions.
Definition: p1localbasis.hh:73
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:15
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: p1localbasis.hh:41
unsigned int order() const
Polynomial order of the shape functions.
Definition: p1localbasis.hh:100