5#ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_PYRAMID_LOCALBASIS_HH
6#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_PYRAMID_LOCALBASIS_HH
11#include <dune/common/fmatrix.hh>
26 template<
class D,
class R>
45 offset_[0][0] = {0.0, 0.0, -1.0};
47 offset_[0][1] = {0.0, 0.0, -1.0};
50 offset_[1][0] = {-2.0, -2.0, 0.0};
52 offset_[1][1] = {0.0, 0.0, 0.0};
55 offset_[2][0] = {0.0, 0.0, 0.0};
57 offset_[2][1] = {0.0, 0.0, 0.0};
60 offset_[3][0] = {0.0, 0.0, 0.0};
62 offset_[3][1] = {-2.0, -2.0, 0.0};
65 offset_[4][0] = {0.0, 0.0, 0.0};
67 offset_[4][1] = {0.0, 0.0, 0.0};
73 offset_[5][0] = {0.0, -2.0, 0.0};
75 offset_[5][1] = {2.0, 0.0, 0.0};
78 for (
size_t i=0; i<5; i++)
79 sign_[i] = s[i] ? -1.0 : 1.0;
98 std::vector<typename Traits::RangeType>& out)
const
101 bool compositeElement = x[0] > x[1];
102 for (std::size_t i=0; i<
size(); i++)
105 out[i] *= factor_[i][compositeElement];
106 out[i] += offset_[i][compositeElement];
118 std::vector<typename Traits::JacobianType>& out)
const
121 bool compositeElement = x[0] > x[1];
122 for (std::size_t i=0; i<
size(); i++)
123 for(std::size_t j=0; j<3; j++)
124 for(std::size_t k=0; k<3; k++)
125 out[i][j][k] = (j==k) * factor_[i][compositeElement] * sign_[i];
131 std::vector<typename Traits::RangeType>& out)
const
133 auto totalOrder = std::accumulate(
order.begin(),
order.end(), 0);
134 if (totalOrder == 0) {
137 DUNE_THROW(NotImplemented,
"Desired derivative order is not implemented");
148 std::array<R,6> sign_;
149 std::array<std::array<typename Traits::RangeType, 2>, 6> offset_;
150 std::array<std::array<R, 2>, 6> factor_;
Definition bdfmcube.hh:18
Type traits for LocalBasisVirtualInterface.
Definition common/localbasis.hh:35
Dune::LocalBasisTraits< D, 3, Dune::FieldVector< D, 3 >, R, 3, Dune::FieldVector< R, 3 >, Dune::FieldMatrix< R, 3, 3 > >::DomainType D DomainType
Definition common/localbasis.hh:43
RT0PyramidLocalBasis(std::bitset< 5 > s=0)
Make set number s, where 0 <= s < 32.
Definition raviartthomas0pyramidlocalbasis.hh:39
void evaluateFunction(const typename Traits::DomainType &x, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition raviartthomas0pyramidlocalbasis.hh:97
LocalBasisTraits< D, 3, Dune::FieldVector< D, 3 >, R, 3, Dune::FieldVector< R, 3 >, Dune::FieldMatrix< R, 3, 3 > > Traits
Definition raviartthomas0pyramidlocalbasis.hh:32
unsigned int order() const
Polynomial order of the shape functions.
Definition raviartthomas0pyramidlocalbasis.hh:142
unsigned int size() const
number of shape functions
Definition raviartthomas0pyramidlocalbasis.hh:86
void evaluateJacobian(const typename Traits::DomainType &x, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition raviartthomas0pyramidlocalbasis.hh:117
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 raviartthomas0pyramidlocalbasis.hh:129