dune-localfunctions  2.6-git
raviartthomas0cube2dall.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_RAVIARTTHOMAS0_CUBE2D_ALL_HH
4 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_CUBE2D_ALL_HH
5 
6 #include <cstddef>
7 #include <numeric>
8 #include <vector>
9 
10 #include <dune/common/fmatrix.hh>
11 
14 
15 namespace Dune
16 {
25  template<class D, class R>
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  std::fill(sign_.begin(), sign_.end(), 1.0);
36  }
37 
39  RT0Cube2DLocalBasis (std::bitset<4> s)
40  {
41  for (int i=0; i<4; i++)
42  sign_[i] = s[i] ? -1.0 : 1.0;
43  }
44 
46  unsigned int size () const
47  {
48  return 4;
49  }
50 
52  inline void evaluateFunction (const typename Traits::DomainType& in,
53  std::vector<typename Traits::RangeType>& out) const
54  {
55  out.resize(4);
56  out[0] = {sign_[0]*(in[0]-1.0), 0.0};
57  out[1] = {sign_[1]*(in[0]), 0.0};
58  out[2] = {0.0, sign_[2]*(in[1]-1.0)};
59  out[3] = {0.0, sign_[3]*(in[1])};
60  }
61 
63  inline void
64  evaluateJacobian (const typename Traits::DomainType& in, // position
65  std::vector<typename Traits::JacobianType>& out) const // return value
66  {
67  out.resize(4);
68  out[0][0] = {sign_[0], 0};
69  out[0][1] = {0, 0};
70 
71  out[1][0] = {sign_[1], 0};
72  out[1][1] = {0, 0};
73 
74  out[2][0] = {0, 0};
75  out[2][1] = {0, sign_[2]};
76 
77  out[3][0] = {0, 0};
78  out[3][1] = {0, sign_[3]};
79  }
80 
82  void partial (const std::array<unsigned int, 2>& order,
83  const typename Traits::DomainType& in, // position
84  std::vector<typename Traits::RangeType>& out) const // return value
85  {
86  auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
87  if (totalOrder == 0) {
88  evaluateFunction(in, out);
89  } else if (totalOrder == 1) {
90  auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
91  out.resize(size());
92 
93  for (std::size_t i = 0; i < size(); ++i)
94  out[i] = {0, 0};
95 
96  switch (direction) {
97  case 0:
98  out[0][0] = sign_[0];
99  out[1][0] = sign_[1];
100  break;
101  case 1:
102  out[2][1] = sign_[2];
103  out[3][1] = sign_[3];
104  break;
105  default:
106  DUNE_THROW(RangeError, "Component out of range.");
107  }
108  } else {
109  out.resize(size());
110  for (std::size_t i = 0; i < size(); ++i)
111  for (std::size_t j = 0; j < 2; ++j)
112  out[i][j] = 0;
113  }
114 
115  }
116 
118  unsigned int order () const
119  {
120  return 1;
121  }
122 
123  private:
124  std::array<R,4> sign_;
125  };
126 
127 
135  template<class LB>
137  {
138  public:
139 
141  RT0Cube2DLocalInterpolation (std::bitset<4> s = 0)
142  {
143  for (int i=0; i<4; i++)
144  sign_[i] = s[i] ? -1.0 : 1.0;
145 
146  m0 = {0.0, 0.5};
147  m1 = {1.0, 0.5};
148  m2 = {0.5, 0.0};
149  m3 = {0.5, 1.0};
150 
151  n0 = {-1.0, 0.0};
152  n1 = { 1.0, 0.0};
153  n2 = { 0.0, -1.0};
154  n3 = { 0.0, 1.0};
155  }
156 
157  template<typename F, typename C>
158  void interpolate (const F& f, std::vector<C>& out) const
159  {
160  // f gives v*outer normal at a point on the edge!
161  typename F::Traits::RangeType y;
162 
163  out.resize(4);
164 
165  // Evaluate the normal components at the edge midpoints
166  f.evaluate(m0,y); out[0] = (y[0]*n0[0]+y[1]*n0[1])*sign_[0];
167  f.evaluate(m1,y); out[1] = (y[0]*n1[0]+y[1]*n1[1])*sign_[1];
168  f.evaluate(m2,y); out[2] = (y[0]*n2[0]+y[1]*n2[1])*sign_[2];
169  f.evaluate(m3,y); out[3] = (y[0]*n3[0]+y[1]*n3[1])*sign_[3];
170  }
171 
172  private:
173  std::array<typename LB::Traits::RangeFieldType,4> sign_;
174 
175  // The four edge midpoints of the reference quadrilateral
176  typename LB::Traits::DomainType m0,m1,m2,m3;
177 
178  // The four edge normals of the reference quadrilateral
179  typename LB::Traits::DomainType n0,n1,n2,n3;
180  };
181 
189  {
190  public:
193  {
194  for (std::size_t i=0; i<4; i++)
195  li[i] = LocalKey(i,1,0);
196  }
197 
199  std::size_t size () const
200  {
201  return 4;
202  }
203 
205  const LocalKey& localKey (std::size_t i) const
206  {
207  return li[i];
208  }
209 
210  private:
211  std::vector<LocalKey> li;
212  };
213 
214 }
215 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_CUBE2D_ALL_HH
Dune::RT0Cube2DLocalBasis::size
unsigned int size() const
number of shape functions
Definition: raviartthomas0cube2dall.hh:46
localbasis.hh
Dune::RT0Cube2DLocalBasis::evaluateJacobian
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas0cube2dall.hh:64
Dune::RT0Cube2DLocalCoefficients
Layout map for RT0 elements on quadrilaterals.
Definition: raviartthomas0cube2dall.hh:188
localkey.hh
Dune::RT0Cube2DLocalInterpolation::RT0Cube2DLocalInterpolation
RT0Cube2DLocalInterpolation(std::bitset< 4 > s=0)
Constructor with explicitly given edge orientations.
Definition: raviartthomas0cube2dall.hh:141
Dune::RT0Cube2DLocalBasis::Traits
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: raviartthomas0cube2dall.hh:30
Dune::RT0Cube2DLocalBasis::RT0Cube2DLocalBasis
RT0Cube2DLocalBasis()
Standard constructor.
Definition: raviartthomas0cube2dall.hh:33
Dune::RT0Cube2DLocalBasis
Lowest order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas0cube2dall.hh:26
Dune::RT0Cube2DLocalBasis::partial
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: raviartthomas0cube2dall.hh:82
Dune::RT0Cube2DLocalInterpolation::interpolate
void interpolate(const F &f, std::vector< C > &out) const
Definition: raviartthomas0cube2dall.hh:158
Dune::RT0Cube2DLocalBasis::RT0Cube2DLocalBasis
RT0Cube2DLocalBasis(std::bitset< 4 > s)
Constructor with a set of edge orientations.
Definition: raviartthomas0cube2dall.hh:39
Dune::RT0Cube2DLocalBasis::order
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas0cube2dall.hh:118
Dune::RT0Cube2DLocalCoefficients::RT0Cube2DLocalCoefficients
RT0Cube2DLocalCoefficients()
Standard constructor.
Definition: raviartthomas0cube2dall.hh:192
Dune::LocalBasisTraits::DomainType
D DomainType
domain type
Definition: localbasis.hh:43
Dune::RT0Cube2DLocalCoefficients::size
std::size_t size() const
number of coefficients
Definition: raviartthomas0cube2dall.hh:199
Dune::LocalBasisTraits
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:31
Dune::RT0Cube2DLocalInterpolation
Lowest order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas0cube2dall.hh:136
Dune::RT0Cube2DLocalBasis::evaluateFunction
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas0cube2dall.hh:52
Dune::LocalKey
Describe position of one degree of freedom.
Definition: localkey.hh:20
Dune::RT0Cube2DLocalCoefficients::localKey
const LocalKey & localKey(std::size_t i) const
get i'th index
Definition: raviartthomas0cube2dall.hh:205
Dune
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:15