5#ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH 
    6#define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH 
   15#include "../../common/localbasis.hh" 
   28  template<
class D, 
class R>
 
   39      for (
size_t i=0; i<4; i++)
 
   50      for (
size_t i=0; i<4; i++)
 
   51        sign_[i] = s[i] ? -1.0 : 1.0;
 
   67                                  std::vector<typename Traits::RangeType>& out)
 const 
   71      out[0][0] = sign_[0]*(in[0] - 1.0);
 
   73      out[1][0] = 6.0*in[0]*in[1] - 3.0*in[0]-6*in[1] + 3.0;
 
   74      out[1][1] = -3.0*in[1]*in[1] + 3.0*in[1];
 
   75      out[2][0] = sign_[1]*(in[0]);
 
   77      out[3][0] = -6.0*in[0]*in[1] + 3.0*in[0];
 
   78      out[3][1] = 3.0*in[1]*in[1] - 3.0*in[1];
 
   80      out[4][1] = sign_[2]*(in[1] - 1.0);
 
   81      out[5][0] = 3.0*in[0]*in[0] - 3.0*in[0];
 
   82      out[5][1] = -6.0*in[0]*in[1] + 6.0*in[0] + 3.0*in[1] - 3.0;
 
   84      out[6][1] = sign_[3]*(in[1]);
 
   85      out[7][0] = -3.0*in[0]*in[0] + 3.0*in[0];
 
   86      out[7][1] = 6.0*in[0]*in[1] - 3.0*in[1];
 
   96                                  std::vector<typename Traits::JacobianType>& out)
 const 
  100      out[0][0][0] = sign_[0];
 
  105      out[1][0][0] = 6.0*in[1] - 3.0;
 
  106      out[1][0][1] = 6.0*in[0] - 6.0;
 
  108      out[1][1][1] = -6.0*in[1] + 3.0;
 
  110      out[2][0][0] = sign_[1];
 
  115      out[3][0][0] = -6.0*in[1] + 3.0;
 
  116      out[3][0][1] = -6.0*in[0];
 
  118      out[3][1][1] = 6.0*in[1] - 3.0;
 
  123      out[4][1][1] = sign_[2];
 
  125      out[5][0][0] = 6.0*in[0] - 3.0;
 
  127      out[5][1][0] = -6.0*in[1] + 6.0;
 
  128      out[5][1][1] = -6.0*in[0] + 3.0;
 
  133      out[6][1][1] = sign_[3];
 
  135      out[7][0][0] = -6.0*in[0] + 3.0;
 
  137      out[7][1][0] = 6.0*in[1];
 
  138      out[7][1][1] = 6.0*in[0] - 3.0;
 
  144                  std::vector<typename Traits::RangeType>& out) 
const       
  147      if (totalOrder == 0) {
 
  149      } 
else if (totalOrder == 1) {
 
  151        auto const direction = std::distance(
order.begin(), std::find(
order.begin(), 
order.end(), 1));
 
  155          out[0][0] = sign_[0];
 
  158          out[1][0] = 6.0*in[1] - 3.0;
 
  161          out[2][0] = sign_[1];
 
  164          out[3][0] = -6.0*in[1] + 3.0;
 
  170          out[5][0] = 6.0*in[0] - 3.0;
 
  171          out[5][1] = -6.0*in[1] + 6.0;
 
  176          out[7][0] = -6.0*in[0] + 3.0;
 
  177          out[7][1] = 6.0*in[1];
 
  183          out[1][0] = 6.0*in[0] - 6.0;
 
  184          out[1][1] = -6.0*in[1] + 3.0;
 
  189          out[3][0] = -6.0*in[0];
 
  190          out[3][1] = 6.0*in[1] - 3.0;
 
  193          out[4][1] = sign_[2];
 
  196          out[5][1] = -6.0*in[0] + 3.0;
 
  199          out[6][1] = sign_[3];
 
  202          out[7][1] = 6.0*in[0] - 3.0;
 
  219    std::array<R,4> sign_;
 
First order Brezzi-Douglas-Marini shape functions on the reference quadrilateral.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:30
 
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:95
 
void partial(const std::array< unsigned int, 2 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:142
 
BDM1Cube2DLocalBasis()
Standard constructor.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:37
 
unsigned int size() const
number of shape functions
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:55
 
BDM1Cube2DLocalBasis(std::bitset< 4 > s)
Make set number s, where 0 <= s < 16.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:48
 
unsigned int order() const
Polynomial order of the shape functions.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:213
 
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:66
 
A dense n x m matrix.
Definition: fmatrix.hh:117
 
Default exception for dummy implementations.
Definition: exceptions.hh:357
 
Default exception class for range errors.
Definition: exceptions.hh:348
 
Implements a matrix constructed from a given type representing a field and compile-time given number ...
 
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
 
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:280
 
Dune namespace.
Definition: alignedallocator.hh:13
 
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:35
 
D DomainType
domain type
Definition: localbasis.hh:43