3 #ifndef DUNE_P1_LOCALBASIS_HH 4 #define DUNE_P1_LOCALBASIS_HH 8 #include <dune/common/fmatrix.hh> 25 template<
class D,
class R,
int dim>
31 Dune::FieldMatrix<R,1,dim>, 2>
Traits;
41 std::vector<typename Traits::RangeType>& out)
const 45 for (
size_t i=0; i<dim; i++) {
54 std::vector<typename Traits::JacobianType>& out)
const 58 for (
int i=0; i<dim; i++)
61 for (
int i=0; i<dim; i++)
62 for (
int j=0; j<dim; j++)
63 out[i+1][0][j] = (i==j);
68 template<
unsigned int k>
69 inline void evaluate (
const typename std::array<int,k>& directions,
71 std::vector<typename Traits::RangeType>& out)
const 80 for (
int i=0; i<dim; i++)
81 out[i+1] = (i==directions[0]);
87 for (
int i=0; i<dim+1; i++)
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
Linear Lagrange shape functions on the simplex.
Definition: p1localbasis.hh:26
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: p1localbasis.hh:40
unsigned int order() const
Polynomial order of the shape functions.
Definition: p1localbasis.hh:93
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:14
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: p1localbasis.hh:53
unsigned int size() const
number of shape functions
Definition: p1localbasis.hh:34
D DomainType
domain type
Definition: localbasis.hh:49
LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim >, 2 > Traits
export type traits for function signature
Definition: p1localbasis.hh:31
void evaluate(const typename std::array< int, k > &directions, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: p1localbasis.hh:69