3 #ifndef DUNE_P1_LOCALBASIS_HH
4 #define DUNE_P1_LOCALBASIS_HH
6 #include <dune/common/fmatrix.hh>
23 template<
class D,
class R,
int dim>
29 Dune::FieldMatrix<R,1,dim>, 2>
Traits;
39 std::vector<typename Traits::RangeType>& out)
const
43 for (
size_t i=0; i<dim; i++) {
52 std::vector<typename Traits::JacobianType>& out)
const
56 for (
int i=0; i<dim; i++)
59 for (
int i=0; i<dim; i++)
60 for (
int j=0; j<dim; j++)
61 out[i+1][0][j] = (i==j);
66 template<
unsigned int k>
67 inline void evaluate (
const typename Dune::array<int,k>& directions,
69 std::vector<typename Traits::RangeType>& out)
const
78 for (
int i=0; i<dim; i++)
79 out[i+1] = (i==directions[0]);
85 for (
int i=0; i<dim+1; i++)
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:39
unsigned int order() const
Polynomial order of the shape functions.
Definition: p1localbasis.hh:91
void evaluate(const typename Dune::array< int, k > &directions, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: p1localbasis.hh:67
Linear Lagrange shape functions on the simplex.
Definition: p1localbasis.hh:24
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: p1localbasis.hh:51
unsigned int size() const
number of shape functions
Definition: p1localbasis.hh:32
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: p1localbasis.hh:38
D DomainType
domain type
Definition: localbasis.hh:51
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:29