dune-localfunctions  2.8.0
raviartthomas1cube2dlocalbasis.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_RAVIARTTHOMAS1_CUBE2D_LOCALBASIS_HH
4 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALBASIS_HH
5 
6 #include <numeric>
7 #include <vector>
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 
37  RT1Cube2DLocalBasis (std::bitset<4> s = 0)
38  {
39  for (size_t i=0; i<4; i++)
40  sign_[i] = (s[i]) ? -1.0 : 1.0;
41  }
42 
44  unsigned int size () const
45  {
46  return 12;
47  }
48 
55  inline void evaluateFunction (const typename Traits::DomainType& in,
56  std::vector<typename Traits::RangeType>& out) const
57  {
58  out.resize(12);
59 
60  out[0][0] = sign_[0]*(-1.0 + 4.0*in[0]-3*in[0]*in[0]);
61  out[0][1] = 0.0;
62  out[1][0] = 3.0 - 12.0*in[0] - 6.0*in[1] + 24.0*in[0]*in[1]+9*in[0]*in[0] - 18.0*in[0]*in[0]*in[1];
63  out[1][1] = 0.0;
64  out[2][0] = sign_[1]*(-2.0*in[0] + 3.0*in[0]*in[0]);
65  out[2][1] = 0.0;
66  out[3][0] = -6.0*in[0] + 12.0*in[0]*in[1] + 9.0*in[0]*in[0] - 18.0*in[0]*in[0]*in[1];
67  out[3][1] = 0.0;
68  out[4][0] = 0.0;
69  out[4][1] = sign_[2]*(-1.0 + 4.0*in[1] - 3.0*in[1]*in[1]);
70  out[5][0] = 0.0;
71  out[5][1] = -3.0 + 6.0*in[0] + 12.0*in[1] - 24.0*in[0]*in[1] - 9.0*in[1]*in[1] + 18.0*in[0]*in[1]*in[1];
72  out[6][0] = 0.0;
73  out[6][1] = sign_[3]*(-2.0*in[1] + 3.0*in[1]*in[1]);
74  out[7][0] = 0.0;
75  out[7][1] = 6.0*in[1] - 12.0*in[0]*in[1] - 9.0*in[1]*in[1] + 18.0*in[0]*in[1]*in[1];
76  out[8][0] = 24.0*in[0] - 36.0*in[0]*in[1] - 24.0*in[0]*in[0] + 36.0*in[0]*in[0]*in[1];
77  out[8][1] = 0.0;
78  out[9][0] = 0.0;
79  out[9][1] = 24.0*in[1] - 36.0*in[0]*in[1] - 24.0*in[1]*in[1] + 36.0*in[0]*in[1]*in[1];
80  out[10][0] = -36.0*in[0] + 72.0*in[0]*in[1] + 36.0*in[0]*in[0] - 72.0*in[0]*in[0]*in[1];
81  out[10][1] = 0.0;
82  out[11][0] = 0.0;
83  out[11][1] = -36.0*in[1] + 72.0*in[0]*in[1] + 36*in[1]*in[1] - 72.0*in[0]*in[1]*in[1];
84  }
85 
92  inline void evaluateJacobian (const typename Traits::DomainType& in,
93  std::vector<typename Traits::JacobianType>& out) const
94  {
95  out.resize(12);
96 
97  out[0][0][0] = sign_[0]*(4.0 - 6.0*in[0]);
98  out[0][0][1] = 0.0;
99  out[0][1][0] = 0.0;
100  out[0][1][1] = 0.0;
101 
102  out[1][0][0] = -12.0 + 24.0*in[1] + 18.0*in[0] - 36.0*in[0]*in[1];
103  out[1][0][1] = -6 + 24.0*in[0] - 18.0*in[0]*in[0];
104  out[1][1][0] = 0.0;
105  out[1][1][1] = 0.0;
106 
107  out[2][0][0] = sign_[1]*(-2.0 + 6.0*in[0]);
108  out[2][0][1] = 0.0;
109  out[2][1][0] = 0.0;
110  out[2][1][1] = 0.0;
111 
112  out[3][0][0] = -6.0 + 12.0*in[1] + 18.0*in[0] - 36.0*in[0]*in[1];
113  out[3][0][1] = 12.0*in[0] - 18.0*in[0]*in[0];
114  out[3][1][0] = 0.0;
115  out[3][1][1] = 0.0;
116 
117  out[4][0][0] = 0.0;
118  out[4][0][1] = 0.0;
119  out[4][1][0] = 0.0;
120  out[4][1][1] = sign_[2]*(4.0 - 6.0*in[1]);
121 
122  out[5][0][0] = 0.0;
123  out[5][0][1] = 0.0;
124  out[5][1][0] = 6.0 - 24.0*in[1] + 18.0*in[1]*in[1];
125  out[5][1][1] = 12.0 - 24.0*in[0] - 18.0*in[1] + 36.0*in[0]*in[1];
126 
127  out[6][0][0] = 0.0;
128  out[6][0][1] = 0.0;
129  out[6][1][0] = 0.0;
130  out[6][1][1] = sign_[3]*(-2.0 + 6.0*in[1]);
131 
132  out[7][0][0] = 0.0;
133  out[7][0][1] = 0.0;
134  out[7][1][0] = -12.0*in[1] + 18.0*in[1]*in[1];
135  out[7][1][1] = 6.0 - 12.0*in[0] - 18.0*in[1] + 36.0*in[1]*in[0];
136 
137  out[8][0][0] = 24.0 - 36.0*in[1] - 48.0*in[0] + 72.0*in[0]*in[1];
138  out[8][0][1] = -36.0*in[0] + 36.0*in[0]*in[0];
139  out[8][1][0] = 0.0;
140  out[8][1][1] = 0.0;
141 
142  out[9][0][0] = 0.0;
143  out[9][0][1] = 0.0;
144  out[9][1][0] = -36.0*in[1] + 36.0*in[1]*in[1];
145  out[9][1][1] = 24.0 - 36.0*in[0] - 48.0*in[1] + 72.0*in[0]*in[1];
146 
147  out[10][0][0] = -36.0 + 72.0*in[1] + 72.0*in[0] - 144.0*in[0]*in[1];
148  out[10][0][1] = 72.0*in[0] - 72.0*in[0]*in[0];
149  out[10][1][0] = 0.0;
150  out[10][1][1] = 0.0;
151 
152  out[11][0][0] = 0.0;
153  out[11][0][1] = 0.0;
154  out[11][1][0] = 72.0*in[1] - 72.0*in[1]*in[1];
155  out[11][1][1] = -36.0 + 72.0*in[0] + 72.0*in[1] - 144.0*in[0]*in[1];
156  }
157 
159  void partial (const std::array<unsigned int, 2>& order,
160  const typename Traits::DomainType& in, // position
161  std::vector<typename Traits::RangeType>& out) const // return value
162  {
163  auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
164  if (totalOrder == 0) {
165  evaluateFunction(in, out);
166  } else {
167  DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
168  }
169  }
170 
172  unsigned int order () const
173  {
174  return 3;
175  }
176 
177  private:
178  std::array<R,4> sign_;
179  };
180 }
181 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS1_CUBE2D_LOCALBASIS_HH
Definition: bdfmcube.hh:16
Type traits for LocalBasisVirtualInterface.
Definition: common/localbasis.hh:32
D DomainType
domain type
Definition: common/localbasis.hh:43
First order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas1cube2dlocalbasis.hh:26
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas1cube2dlocalbasis.hh:92
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: raviartthomas1cube2dlocalbasis.hh:30
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas1cube2dlocalbasis.hh:172
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas1cube2dlocalbasis.hh:55
unsigned int size() const
number of shape functions
Definition: raviartthomas1cube2dlocalbasis.hh:44
RT1Cube2DLocalBasis(std::bitset< 4 > s=0)
Make set number s, where 0 <= s < 16.
Definition: raviartthomas1cube2dlocalbasis.hh:37
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: raviartthomas1cube2dlocalbasis.hh:159