3#ifndef DUNE_PYRAMID_P1_LOCALBASIS_HH
4#define DUNE_PYRAMID_P1_LOCALBASIS_HH
10#include <dune/localfunctions/common/localbasis.hh>
25 template<
class D,
class R>
41 std::vector<typename Traits::RangeType>& out)
const
47 out[0] = (1-in[0])*(1-in[1])-in[2]*(1-in[1]);
48 out[1] = in[0]*(1-in[1])-in[2]*in[1];
49 out[2] = (1-in[0])*in[1]-in[2]*in[1];
50 out[3] = in[0]*in[1]+in[2]*in[1];
54 out[0] = (1-in[0])*(1-in[1])-in[2]*(1-in[0]);
55 out[1] = in[0]*(1-in[1])-in[2]*in[0];
56 out[2] = (1-in[0])*in[1]-in[2]*in[0];
57 out[3] = in[0]*in[1]+in[2]*in[0];
69 std::vector<typename Traits::JacobianType>& out)
const
75 out[0][0][0] = -1 + in[1]; out[0][0][1] = -1 + in[0] + in[2]; out[0][0][2] = -1 + in[1];
76 out[1][0][0] = 1 - in[1]; out[1][0][1] = -in[0] - in[2]; out[1][0][2] = -in[1];
77 out[2][0][0] = -in[1]; out[2][0][1] = 1 - in[0] - in[2]; out[2][0][2] = -in[1];
78 out[3][0][0] = in[1]; out[3][0][1] = in[0]+in[2]; out[3][0][2] = in[1];
82 out[0][0][0] = -1 + in[1] + in[2]; out[0][0][1] = -1 + in[0]; out[0][0][2] = -1 + in[0];
83 out[1][0][0] = 1 - in[1] - in[2]; out[1][0][1] = -in[0]; out[1][0][2] = -in[0];
84 out[2][0][0] = -in[1] - in[2]; out[2][0][1] = 1 - in[0]; out[2][0][2] = -in[0];
85 out[3][0][0] = in[1] + in[2]; out[3][0][1] = in[0]; out[3][0][2] = in[0];
89 out[4][0][0] = 0; out[4][0][1] = 0; out[4][0][2] = 1;
95 std::vector<typename Traits::RangeType>& out)
const
98 if (totalOrder == 0) {
100 }
else if (totalOrder == 1) {
103 auto const direction = std::distance(
order.begin(), std::find(
order.begin(),
order.end(), 1));
115 out[0] = -1 + in[0] + in[2];
116 out[1] = -in[0] - in[2];
117 out[2] = 1 - in[0] - in[2];
118 out[3] = in[0]+in[2];
136 out[0] = -1 + in[1] + in[2];
137 out[1] = 1 - in[1] - in[2];
138 out[2] = -in[1] - in[2];
139 out[3] = in[1] + in[2];
160 }
else if (totalOrder == 2) {
164 (
order[1] == 1 &&
order[2] == 1 && in[0] > in[1]) ||
165 (
order[0] == 1 &&
order[2] == 1 && in[0] <=in[1])) {
172 for (std::size_t i = 0; i <
size(); ++i)
178 for (std::size_t i = 0; i <
size(); ++i)
A dense n x m matrix.
Definition: fmatrix.hh:68
vector space out of a tensor product of fields.
Definition: fvector.hh:93
Linear Lagrange shape functions on the pyramid.
Definition: pyramidp1localbasis.hh:27
void partial(const std::array< unsigned int, 3 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: pyramidp1localbasis.hh:93
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: pyramidp1localbasis.hh:68
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: pyramidp1localbasis.hh:40
LocalBasisTraits< D, 3, Dune::FieldVector< D, 3 >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, 3 > > Traits
export type traits for function signature
Definition: pyramidp1localbasis.hh:31
unsigned int order() const
Polynomial order of the shape functions.
Definition: pyramidp1localbasis.hh:184
unsigned int size() const
number of shape functions
Definition: pyramidp1localbasis.hh:34
Default exception class for range errors.
Definition: exceptions.hh:252
Implements a matrix constructed from a given type representing a field and compile-time given number ...
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:331
Dune namespace.
Definition: alignedallocator.hh:10
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:32
D DomainType
domain type
Definition: localbasis.hh:43