dune-localfunctions  2.3.1
brezzidouglasmarini1cube2dlocalbasis.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH
4 #define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH
5 
6 #include <vector>
7 #include <bitset>
8 
9 #include <dune/common/fmatrix.hh>
10 
11 #include "../../common/localbasis.hh"
12 
13 namespace Dune
14 {
24  template<class D, class R>
26  {
27 
28  public:
29  typedef LocalBasisTraits<D,2,Dune::FieldVector<D,2>,R,2,Dune::FieldVector<R,2>,
30  Dune::FieldMatrix<R,2,2> > Traits;
31 
34  {
35  for (size_t i=0; i<4; i++)
36  sign_[i] = 1.0;
37  }
38 
44  BDM1Cube2DLocalBasis (std::bitset<4> s)
45  {
46  for (size_t i=0; i<4; i++)
47  sign_[i] = s[i] ? -1.0 : 1.0;
48  }
49 
51  unsigned int size () const
52  {
53  return 8;
54  }
55 
62  inline void evaluateFunction (const typename Traits::DomainType& in,
63  std::vector<typename Traits::RangeType>& out) const
64  {
65  out.resize(8);
66 
67  out[0][0] = sign_[0]*(in[0] - 1.0);
68  out[0][1] = 0.0;
69  out[1][0] = 6.0*in[0]*in[1] - 3.0*in[0]-6*in[1] + 3.0;
70  out[1][1] = -3.0*in[1]*in[1] + 3.0*in[1];
71  out[2][0] = sign_[1]*(in[0]);
72  out[2][1] = 0.0;
73  out[3][0] = -6.0*in[0]*in[1] + 3.0*in[0];
74  out[3][1] = 3.0*in[1]*in[1] - 3.0*in[1];
75  out[4][0] = 0.0;
76  out[4][1] = sign_[2]*(in[1] - 1.0);
77  out[5][0] = 3.0*in[0]*in[0] - 3.0*in[0];
78  out[5][1] = -6.0*in[0]*in[1] + 6.0*in[0] + 3.0*in[1] - 3.0;
79  out[6][0] = 0.0;
80  out[6][1] = sign_[3]*(in[1]);
81  out[7][0] = -3.0*in[0]*in[0] + 3.0*in[0];
82  out[7][1] = 6.0*in[0]*in[1] - 3.0*in[1];
83  }
84 
91  inline void evaluateJacobian (const typename Traits::DomainType& in,
92  std::vector<typename Traits::JacobianType>& out) const
93  {
94  out.resize(8);
95 
96  out[0][0][0] = sign_[0];
97  out[0][0][1] = 0.0;
98  out[0][1][0] = 0.0;
99  out[0][1][1] = 0.0;
100 
101  out[1][0][0] = 6.0*in[1] - 3.0;
102  out[1][0][1] = 6.0*in[0] - 6.0;
103  out[1][1][0] = 0.0;
104  out[1][1][1] = -6.0*in[1] + 3.0;
105 
106  out[2][0][0] = sign_[1];
107  out[2][0][1] = 0.0;
108  out[2][1][0] = 0.0;
109  out[2][1][1] = 0.0;
110 
111  out[3][0][0] = -6.0*in[1] + 3.0;
112  out[3][0][1] = -6.0*in[0];
113  out[3][1][0] = 0.0;
114  out[3][1][1] = 6.0*in[1] - 3.0;
115 
116  out[4][0][0] = 0.0;
117  out[4][0][1] = 0.0;
118  out[4][1][0] = 0.0;
119  out[4][1][1] = sign_[2];
120 
121  out[5][0][0] = 6.0*in[0] - 3.0;
122  out[5][0][1] = 0.0;
123  out[5][1][0] = -6.0*in[1] + 6.0;
124  out[5][1][1] = -6.0*in[0] + 3.0;
125 
126  out[6][0][0] = 0.0;
127  out[6][0][1] = 0.0;
128  out[6][1][0] = 0.0;
129  out[6][1][1] = sign_[3];
130 
131  out[7][0][0] = -6.0*in[0] + 3.0;
132  out[7][0][1] = 0.0;
133  out[7][1][0] = 6.0*in[1];
134  out[7][1][1] = 6.0*in[0] - 3.0;
135  }
136 
138  unsigned int order () const
139  {
140  return 2;
141  }
142 
143  private:
144  array<R,4> sign_;
145  };
146 }
147 #endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:39
unsigned int order() const
Polynomial order of the shape functions.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:138
BDM1Cube2DLocalBasis(std::bitset< 4 > s)
Make set number s, where 0 <= s < 16.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:44
unsigned int size() const
number of shape functions
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:51
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:30
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:62
BDM1Cube2DLocalBasis()
Standard constructor.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:33
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:91
First order Brezzi-Douglas-Marini shape functions on the reference quadrilateral. ...
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:25
D DomainType
domain type
Definition: localbasis.hh:51