3#ifndef DUNE_DUAL_P1_LOCALBASIS_HH
4#define DUNE_DUAL_P1_LOCALBASIS_HH
8#include <dune/common/fvector.hh>
9#include <dune/common/fmatrix.hh>
30 template<
class D,
class R,
int dim,
bool faceDualT=false>
48 std::vector<typename Traits::RangeType>& out)
const
51 std::vector<typename Traits::RangeType> p1Values(
size());
55 for (
int i=0; i<dim; i++) {
57 p1Values[i+1] = in[i];
63 for (
int i=0; i<=dim; i++) {
64 out[i] = (dim+!
faceDual)*p1Values[i];
65 for (
int j=0; j<i; j++)
66 out[i] -= p1Values[j];
68 for (
int j=i+1; j<=dim; j++)
69 out[i] -= p1Values[j];
76 std::vector<typename Traits::JacobianType>& out)
const
79 std::vector<typename Traits::JacobianType> p1Jacs(
size());
81 for (
int i=0; i<dim; i++)
84 for (
int i=0; i<dim; i++)
85 for (
int j=0; j<dim; j++)
86 p1Jacs[i+1][0][j] = (i==j);
91 for (
size_t i=0; i<=dim; i++) {
93 out[i][0].axpy(dim+!
faceDual,p1Jacs[i][0]);
95 for (
size_t j=0; j<i; j++)
96 out[i][0] -= p1Jacs[j][0];
98 for (
int j=i+1; j<=dim; j++)
99 out[i][0] -= p1Jacs[j][0];
106 std::vector<typename Traits::RangeType>& out)
const
108 auto totalOrder = std::accumulate(
order.begin(),
order.end(), 0);
109 if (totalOrder == 0) {
112 DUNE_THROW(NotImplemented,
"Desired derivative order is not implemented");
Definition: bdfmcube.hh:16
Type traits for LocalBasisVirtualInterface.
Definition: common/localbasis.hh:32
D DomainType
domain type
Definition: common/localbasis.hh:43
Dual Lagrange shape functions on the simplex.
Definition: dualp1localbasis.hh:32
unsigned int order() const
Polynomial order of the shape functions.
Definition: dualp1localbasis.hh:117
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: dualp1localbasis.hh:47
static const bool faceDual
Determines if the basis is only biorthogonal on adjacent faces.
Definition: dualp1localbasis.hh:35
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 all shape functions.
Definition: dualp1localbasis.hh:104
unsigned int size() const
number of shape functions
Definition: dualp1localbasis.hh:41
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: dualp1localbasis.hh:38
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: dualp1localbasis.hh:75