dune-localfunctions  2.8.0
raviartthomas12dlocalbasis.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_RAVIARTTHOMAS12DLOCALBASIS_HH
4 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS12DLOCALBASIS_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  RT12DLocalBasis (std::bitset<3> s = 0)
38  {
39  for (size_t i=0; i<3; i++)
40  sign_[i] = (s[i]) ? -1.0 : 1.0;
41  }
42 
44  unsigned int size () const
45  {
46  return 8;
47  }
48 
55  inline void evaluateFunction (const typename Traits::DomainType& in,
56  std::vector<typename Traits::RangeType>& out) const
57  {
58  out.resize(8);
59  out[0][0] = sign_[0]*(in[0] - 4.0*in[0]*in[1]);
60  out[0][1] = sign_[0]*(-1.0 + 5.0*in[1] - 4.0*in[1]*in[1]);
61  out[1][0] = sign_[1]*(-1.0 + 5.0*in[0] - 4.0*in[0]*in[0]);
62  out[1][1] = sign_[1]*(in[1] - 4.0*in[0]*in[1]);
63  out[2][0] = sign_[2]*(-3.0*in[0] + 4.0*in[0]*in[0] + 4.0*in[1]*in[0]);
64  out[2][1] = sign_[2]*(-3.0*in[1] + 4.0*in[0]*in[1] + 4.0*in[1]*in[1]);
65  out[3][0] = -5.0*in[0] + 8.0*in[0]*in[0] + 4.0*in[1]*in[0];
66  out[3][1] = 3.0 - 6.0*in[0] - 7.0*in[1] + 8.0*in[0]*in[1] + 4.0*in[1]*in[1];
67  out[4][0] = -3.0 + 7.0*in[0] + 6.0*in[1] - 4.0*in[0]*in[0] - 8.0*in[1]*in[0];
68  out[4][1] = 5.0*in[1] - 4.0*in[0]*in[1] - 8.0*in[1]*in[1];
69  out[5][0] = in[0] - 4.0*in[0]*in[0] + 4.0*in[1]*in[0];
70  out[5][1] = -1.0*in[1] - 4.0*in[0]*in[1] + 4.0*in[1]*in[1];
71  out[6][0] = 16.0*in[0] - 16.0*in[0]*in[0] - 8.0*in[1]*in[0];
72  out[6][1] = 8.0*in[1] - 16.0*in[0]*in[1] - 8.0*in[1]*in[1];
73  out[7][0] = 8.0*in[0] - 8.0*in[0]*in[0] - 16.0*in[1]*in[0];
74  out[7][1] = 16.0*in[1] - 8.0*in[0]*in[1] - 16.0*in[1]*in[1];
75  }
76 
83  inline void evaluateJacobian (const typename Traits::DomainType& in,
84  std::vector<typename Traits::JacobianType>& out) const
85  {
86  out.resize(8);
87 
88  out[0][0][0] = sign_[0]*(1.0 - 4.0*in[1]);
89  out[0][0][1] = sign_[0]*(-4.0*in[0]);
90  out[0][1][0] = 0.0;
91  out[0][1][1] = sign_[0]*(5.0 - 8.0*in[1]);
92 
93  out[1][0][0] = sign_[1]*(5.0 - 8.0*in[0]);
94  out[1][0][1] = 0.0;
95  out[1][1][0] = sign_[1]*(-4.0*in[1]);
96  out[1][1][1] = sign_[1]*(1.0 - 4.0*in[0]);
97 
98  out[2][0][0] = sign_[2]*(-3.0 + 8.0*in[0] + 4.0*in[1]);
99  out[2][0][1] = sign_[2]*(4.0*in[0]);
100  out[2][1][0] = sign_[2]*(4.0*in[1]);
101  out[2][1][1] = sign_[2]*(-3.0 + 4.0*in[0] + 8.0*in[1]);
102 
103  out[3][0][0] = -5.0 + 16.0*in[0] + 4.0*in[1];
104  out[3][0][1] = 4.0*in[0];
105  out[3][1][0] = -6.0 + 8.0*in[1];
106  out[3][1][1] = -7.0 + 8.0*in[0] + 8.0*in[1];
107 
108  out[4][0][0] = 7.0 - 8.0*in[0] - 8.0*in[1];
109  out[4][0][1] = 6.0 - 8.0*in[0];
110  out[4][1][0] = -4.0*in[1];
111  out[4][1][1] = 5.0 - 4.0*in[0] - 16.0*in[1];
112 
113  out[5][0][0] = 1.0 - 8.0*in[0] + 4*in[1];
114  out[5][0][1] = 4.0*in[0];
115  out[5][1][0] = -4.0*in[1];
116  out[5][1][1] = -1.0 - 4.0*in[0] + 8.0*in[1];
117 
118  out[6][0][0] = 16.0 - 32.0*in[0] - 8.0*in[1];
119  out[6][0][1] = -8.0*in[0];
120  out[6][1][0] = -16.0*in[1];
121  out[6][1][1] = 8.0 - 16.0*in[0] - 16.0*in[1];
122 
123  out[7][0][0] = 8.0 - 16.0*in[0] - 16.0*in[1];
124  out[7][0][1] = -16.0*in[0];
125  out[7][1][0] = -8.0*in[1];
126  out[7][1][1] = 16.0 - 8.0*in[0] - 32.0*in[1];
127  }
128 
130  void partial (const std::array<unsigned int, 2>& order,
131  const typename Traits::DomainType& in, // position
132  std::vector<typename Traits::RangeType>& out) const // return value
133  {
134  auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
135  if (totalOrder == 0) {
136  evaluateFunction(in, out);
137  } else if (totalOrder == 1) {
138  auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
139  out.resize(size());
140 
141  switch (direction) {
142  case 0:
143  out[0][0] = sign_[0]*(1.0 - 4.0*in[1]);
144  out[0][1] = 0.0;
145  out[1][0] = sign_[1]*(5.0 - 8.0*in[0]);
146  out[1][1] = sign_[1]*(-4.0*in[1]);
147  out[2][0] = sign_[2]*(-3.0 + 8.0*in[0] + 4.0*in[1]);
148  out[2][1] = sign_[2]*(4.0*in[1]);
149  out[3][0] = -5.0 + 16.0*in[0] + 4.0*in[1];
150  out[3][1] = -6.0 + 8.0*in[1];
151  out[4][0] = 7.0 - 8.0*in[0] - 8.0*in[1];
152  out[4][1] = -4.0*in[1];
153  out[5][0] = 1.0 - 8.0*in[0] + 4*in[1];
154  out[5][1] = -4.0*in[1];
155  out[6][0] = 16.0 - 32.0*in[0] - 8.0*in[1];
156  out[6][1] = -16.0*in[1];
157  out[7][0] = 8.0 - 16.0*in[0] - 16.0*in[1];
158  out[7][1] = -8.0*in[1];
159  break;
160  case 1:
161  out[2][1] = sign_[2]*(-3.0 + 4.0*in[0] + 8.0*in[1]);
162  out[2][0] = sign_[2]*(4.0*in[0]);
163  out[1][1] = sign_[1]*(1.0 - 4.0*in[0]);
164  out[1][0] = 0.0;
165  out[0][0] = sign_[0]*(-4.0*in[0]);
166  out[0][1] = sign_[0]*(5.0 - 8.0*in[1]);
167  out[3][0] = 4.0*in[0];
168  out[3][1] = -7.0 + 8.0*in[0] + 8.0*in[1];
169  out[4][0] = 6.0 - 8.0*in[0];
170  out[4][1] = 5.0 - 4.0*in[0] - 16.0*in[1];
171  out[5][0] = 4.0*in[0];
172  out[5][1] = -1.0 - 4.0*in[0] + 8.0*in[1];
173  out[6][0] = -8.0*in[0];
174  out[6][1] = 8.0 - 16.0*in[0] - 16.0*in[1];
175  out[7][0] = -16.0*in[0];
176  out[7][1] = 16.0 - 8.0*in[0] - 32.0*in[1];
177  break;
178  default:
179  DUNE_THROW(RangeError, "Component out of range.");
180  }
181  } else {
182  DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
183  }
184  }
185 
187  unsigned int order () const
188  {
189  return 2;
190  }
191 
192  private:
193  std::array<R,3> sign_;
194  };
195 }
196 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS12DLOCALBASIS_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 triangle.
Definition: raviartthomas12dlocalbasis.hh:26
RT12DLocalBasis(std::bitset< 3 > s=0)
Make set number s, where 0 <= s < 8.
Definition: raviartthomas12dlocalbasis.hh:37
unsigned int size() const
number of shape functions
Definition: raviartthomas12dlocalbasis.hh:44
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: raviartthomas12dlocalbasis.hh:130
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas12dlocalbasis.hh:83
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas12dlocalbasis.hh:187
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: raviartthomas12dlocalbasis.hh:30
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas12dlocalbasis.hh:55