3 #ifndef DUNE_PRISM_P1_LOCALBASIS_HH 4 #define DUNE_PRISM_P1_LOCALBASIS_HH 8 #include <dune/common/fmatrix.hh> 24 template<
class D,
class R>
40 std::vector<typename Traits::RangeType>& out)
const 43 out[0] = (1.0-in[0]-in[1])*(1.0-in[2]);
44 out[1] = in[0]*(1-in[2]);
45 out[2] = in[1]*(1-in[2]);
46 out[3] = in[2]*(1.0-in[0]-in[1]);
54 std::vector<typename Traits::JacobianType>& out)
const 57 out[0][0][0] = in[2]-1; out[0][0][1] = in[2]-1; out[0][0][2] = in[0]+in[1]-1;
58 out[1][0][0] = 1-in[2]; out[1][0][1] = 0; out[1][0][2] = -in[0];
59 out[2][0][0] = 0; out[2][0][1] = 1-in[2]; out[2][0][2] = -in[1];
60 out[3][0][0] = -in[2]; out[3][0][1] = -in[2]; out[3][0][2] = 1-in[0]-in[1];
61 out[4][0][0] = in[2]; out[4][0][1] = 0; out[4][0][2] = in[0];
62 out[5][0][0] = 0; out[5][0][1] = in[2]; out[5][0][2] = in[1];
68 std::vector<typename Traits::RangeType>& out)
const 70 auto totalOrder = std::accumulate(
order.begin(),
order.end(), 0);
71 if (totalOrder == 0) {
73 }
else if (totalOrder == 1) {
74 auto direction = std::distance(
order.begin(), std::find(
order.begin(),
order.end(), 1));
95 out[0] = in[0]+in[1]-1;
98 out[3] = 1-in[0]-in[1];
103 DUNE_THROW(RangeError,
"Component out of range.");
105 }
else if (totalOrder == 2) {
122 for (std::size_t i = 0; i <
size(); ++i)
127 for (std::size_t i = 0; i <
size(); ++i)
unsigned int order() const
Polynomial order of the shape functions.
Definition: prismp1localbasis.hh:133
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: prismp1localbasis.hh:66
unsigned int size() const
number of shape functions
Definition: prismp1localbasis.hh:33
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: prismp1localbasis.hh:39
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:31
Linear Lagrange shape functions on the prism.
Definition: prismp1localbasis.hh:25
D DomainType
domain type
Definition: localbasis.hh:43
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: prismp1localbasis.hh:30
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:15
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: prismp1localbasis.hh:53