3#ifndef DUNE_P0LOCALBASIS_HH
4#define DUNE_P0LOCALBASIS_HH
8#include <dune/common/fmatrix.hh>
26 template<
class D,
class R,
int d>
42 std::vector<typename Traits::RangeType>& out)
const
51 std::vector<typename Traits::JacobianType>& out)
const
54 for (
int i=0; i<d; i++)
65 std::vector<typename Traits::RangeType>& out)
const
67 auto totalOrder = std::accumulate(
order.begin(),
order.end(), 0);
68 if (totalOrder == 0) {
Definition brezzidouglasmarini1cube2d.hh:14
Type traits for LocalBasisVirtualInterface.
Definition localbasis.hh:38
D DomainType
domain type
Definition localbasis.hh:49
Constant shape function.
Definition p0localbasis.hh:28
LocalBasisTraits< D, d, Dune::FieldVector< D, d >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, d >, 0 > Traits
export type traits for function signature
Definition p0localbasis.hh:32
unsigned int order() const
Polynomial order of the shape functions.
Definition p0localbasis.hh:77
void partial(const std::array< unsigned int, d > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of any order of all shape functions.
Definition p0localbasis.hh:63
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition p0localbasis.hh:50
unsigned int size() const
number of shape functions
Definition p0localbasis.hh:35
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition p0localbasis.hh:41