dune-localfunctions  2.7.1
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 (unsigned int s = 0)
38  {
39  sign0 = sign1 = sign2 = sign3 = 1.0;
40  if (s & 1)
41  {
42  sign0 = -1.0;
43  }
44  if (s & 2)
45  {
46  sign1 = -1.0;
47  }
48  if (s & 4)
49  {
50  sign2 = -1.0;
51  }
52  if (s & 8)
53  {
54  sign3 = -1.0;
55  }
56  }
57 
59  unsigned int size () const
60  {
61  return 12;
62  }
63 
70  inline void evaluateFunction (const typename Traits::DomainType& in,
71  std::vector<typename Traits::RangeType>& out) const
72  {
73  out.resize(12);
74 
75  out[0][0] = sign0*(-1.0 + 4.0*in[0]-3*in[0]*in[0]);
76  out[0][1] = 0.0;
77  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];
78  out[1][1] = 0.0;
79  out[2][0] = sign1*(-2.0*in[0] + 3.0*in[0]*in[0]);
80  out[2][1] = 0.0;
81  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];
82  out[3][1] = 0.0;
83  out[4][0] = 0.0;
84  out[4][1] = sign2*(-1.0 + 4.0*in[1] - 3.0*in[1]*in[1]);
85  out[5][0] = 0.0;
86  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];
87  out[6][0] = 0.0;
88  out[6][1] = sign3*(-2.0*in[1] + 3.0*in[1]*in[1]);
89  out[7][0] = 0.0;
90  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];
91  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];
92  out[8][1] = 0.0;
93  out[9][0] = 0.0;
94  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];
95  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];
96  out[10][1] = 0.0;
97  out[11][0] = 0.0;
98  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];
99  }
100 
107  inline void evaluateJacobian (const typename Traits::DomainType& in,
108  std::vector<typename Traits::JacobianType>& out) const
109  {
110  out.resize(12);
111 
112  out[0][0][0] = sign0*(4.0 - 6.0*in[0]);
113  out[0][0][1] = 0.0;
114  out[0][1][0] = 0.0;
115  out[0][1][1] = 0.0;
116 
117  out[1][0][0] = -12.0 + 24.0*in[1] + 18.0*in[0] - 36.0*in[0]*in[1];
118  out[1][0][1] = -6 + 24.0*in[0] - 18.0*in[0]*in[0];
119  out[1][1][0] = 0.0;
120  out[1][1][1] = 0.0;
121 
122  out[2][0][0] = sign1*(-2.0 + 6.0*in[0]);
123  out[2][0][1] = 0.0;
124  out[2][1][0] = 0.0;
125  out[2][1][1] = 0.0;
126 
127  out[3][0][0] = -6.0 + 12.0*in[1] + 18.0*in[0] - 36.0*in[0]*in[1];
128  out[3][0][1] = 12.0*in[0] - 18.0*in[0]*in[0];
129  out[3][1][0] = 0.0;
130  out[3][1][1] = 0.0;
131 
132  out[4][0][0] = 0.0;
133  out[4][0][1] = 0.0;
134  out[4][1][0] = 0.0;
135  out[4][1][1] = sign2*(4.0 - 6.0*in[1]);
136 
137  out[5][0][0] = 0.0;
138  out[5][0][1] = 0.0;
139  out[5][1][0] = 6.0 - 24.0*in[1] + 18.0*in[1]*in[1];
140  out[5][1][1] = 12.0 - 24.0*in[0] - 18.0*in[1] + 36.0*in[0]*in[1];
141 
142  out[6][0][0] = 0.0;
143  out[6][0][1] = 0.0;
144  out[6][1][0] = 0.0;
145  out[6][1][1] = sign3*(-2.0 + 6.0*in[1]);
146 
147  out[7][0][0] = 0.0;
148  out[7][0][1] = 0.0;
149  out[7][1][0] = -12.0*in[1] + 18.0*in[1]*in[1];
150  out[7][1][1] = 6.0 - 12.0*in[0] - 18.0*in[1] + 36.0*in[1]*in[0];
151 
152  out[8][0][0] = 24.0 - 36.0*in[1] - 48.0*in[0] + 72.0*in[0]*in[1];
153  out[8][0][1] = -36.0*in[0] + 36.0*in[0]*in[0];
154  out[8][1][0] = 0.0;
155  out[8][1][1] = 0.0;
156 
157  out[9][0][0] = 0.0;
158  out[9][0][1] = 0.0;
159  out[9][1][0] = -36.0*in[1] + 36.0*in[1]*in[1];
160  out[9][1][1] = 24.0 - 36.0*in[0] - 48.0*in[1] + 72.0*in[0]*in[1];
161 
162  out[10][0][0] = -36.0 + 72.0*in[1] + 72.0*in[0] - 144.0*in[0]*in[1];
163  out[10][0][1] = 72.0*in[0] - 72.0*in[0]*in[0];
164  out[10][1][0] = 0.0;
165  out[10][1][1] = 0.0;
166 
167  out[11][0][0] = 0.0;
168  out[11][0][1] = 0.0;
169  out[11][1][0] = 72.0*in[1] - 72.0*in[1]*in[1];
170  out[11][1][1] = -36.0 + 72.0*in[0] + 72.0*in[1] - 144.0*in[0]*in[1];
171  }
172 
174  void partial (const std::array<unsigned int, 2>& order,
175  const typename Traits::DomainType& in, // position
176  std::vector<typename Traits::RangeType>& out) const // return value
177  {
178  auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
179  if (totalOrder == 0) {
180  evaluateFunction(in, out);
181  } else {
182  DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
183  }
184  }
185 
187  unsigned int order () const
188  {
189  return 3;
190  }
191 
192  private:
193  R sign0, sign1, sign2, sign3;
194  };
195 }
196 #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:107
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:187
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas1cube2dlocalbasis.hh:70
unsigned int size() const
number of shape functions
Definition: raviartthomas1cube2dlocalbasis.hh:59
RT1Cube2DLocalBasis(unsigned int 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:174