3#ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMAS03D_RAVIARTTHOMAS03DLOCALBASIS_HH
4#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMAS03D_RAVIARTTHOMAS03DLOCALBASIS_HH
10#include <dune/localfunctions/common/localbasis.hh>
22 template<
class D,
class R>
32 for (
int i=0; i<4; i++)
33 sign_[i] = s[i] ? -1.0 : 1.0;
44 std::vector<typename Traits::RangeType>& out)
const
47 auto c = std::sqrt(2.0);
48 out[0] = {sign_[0]*c* in[0], sign_[0]*c* in[1], sign_[0]*c*(in[2]-D(1))};
49 out[1] = {sign_[1]*c* in[0], sign_[1]*c*(in[1]-D(1)), sign_[1]*c* in[2] };
50 out[2] = {sign_[2]*c*(in[0]-D(1)), sign_[2]*c* in[1], sign_[2]*c* in[2] };
51 out[3] = {sign_[3]*c* in[0], sign_[3]*c* in[1], sign_[3]*c* in[2] };
57 std::vector<typename Traits::JacobianType>& out)
const
60 for (
int i=0; i<4; i++)
62 auto c = std::sqrt(2.0);
63 out[i][0] = {c*sign_[i], 0, 0};
64 out[i][1] = { 0,c*sign_[i], 0};
65 out[i][2] = { 0, 0,c*sign_[i]};
72 std::vector<typename Traits::RangeType>& out)
const
75 if (totalOrder == 0) {
77 }
else if (totalOrder == 1) {
78 auto const direction = std::distance(
order.begin(), std::find(
order.begin(),
order.end(), 1));
81 for (
int i=0; i<
size(); i++)
83 out[i][direction] = sign_[i]* std::sqrt(2.0) ;
84 out[i][(direction+1)%3] = 0;
85 out[i][(direction+2)%3] = 0;
89 for (std::size_t i = 0; i <
size(); ++i)
90 for (std::size_t j = 0; j < 3; ++j)
105 std::array<R,4> sign_;
A dense n x m matrix.
Definition: fmatrix.hh:69
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Lowest order Raviart-Thomas shape functions on the reference tetrahedron.
Definition: raviartthomas03dlocalbasis.hh:24
RT03DLocalBasis(std::bitset< 4 > s=0)
Make set number s, where 0 <= s < 16.
Definition: raviartthomas03dlocalbasis.hh:30
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: raviartthomas03dlocalbasis.hh:70
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas03dlocalbasis.hh:97
unsigned int size() const
number of shape functions
Definition: raviartthomas03dlocalbasis.hh:37
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas03dlocalbasis.hh:56
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas03dlocalbasis.hh:43
Implements a matrix constructed from a given type representing a field and compile-time given number ...
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:289
Dune namespace.
Definition: alignedallocator.hh:11
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:32
D DomainType
domain type
Definition: localbasis.hh:43