dune-localfunctions  2.4.1
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 <vector>
8 
9 #include <dune/common/fmatrix.hh>
10 
13 
14 namespace Dune
15 {
24  template<class D, class R>
26  {
27  public:
28  typedef LocalBasisTraits<D,2,Dune::FieldVector<D,2>,R,2,Dune::FieldVector<R,2>,
29  Dune::FieldMatrix<R,2,2> > Traits;
30 
33  {
34  sign0 = sign1 = sign2 = sign3 = 1.0;
35  }
36 
38  RT0Cube2DLocalBasis (unsigned int s)
39  {
40  sign0 = sign1 = sign2 = sign3 = 1.0;
41  if (s&1) sign0 = -1.0;
42  if (s&2) sign1 = -1.0;
43  if (s&4) sign2 = -1.0;
44  if (s&8) sign3 = -1.0;
45  }
46 
48  unsigned int size () const
49  {
50  return 4;
51  }
52 
54  inline void evaluateFunction (const typename Traits::DomainType& in,
55  std::vector<typename Traits::RangeType>& out) const
56  {
57  out.resize(4);
58  out[0][0] = sign0*(in[0]-1.0); out[0][1]=0.0;
59  out[1][0] = sign1*(in[0]); out[1][1]=0.0;
60  out[2][0] = 0.0; out[2][1]=sign2*(in[1]-1.0);
61  out[3][0] = 0.0; out[3][1]=sign3*(in[1]);
62  }
63 
65  inline void
66  evaluateJacobian (const typename Traits::DomainType& in, // position
67  std::vector<typename Traits::JacobianType>& out) const // return value
68  {
69  out.resize(4);
70  out[0][0][0] = sign0; out[0][0][1] = 0;
71  out[0][1][0] = 0; out[0][1][1] = 0;
72 
73  out[1][0][0] = sign1; out[1][0][1] = 0;
74  out[1][1][0] = 0; out[1][1][1] = 0;
75 
76  out[2][0][0] = 0; out[2][0][1] = 0;
77  out[2][1][0] = 0; out[2][1][1] = sign2;
78 
79  out[3][0][0] = 0; out[3][0][1] = 0;
80  out[3][1][0] = 0; out[3][1][1] = sign3;
81  }
82 
84  unsigned int order () const
85  {
86  return 1;
87  }
88 
89  private:
90  R sign0, sign1, sign2, sign3;
91  };
92 
93 
101  template<class LB>
103  {
104  public:
105 
108  {
109  sign0 = sign1 = sign2 = sign3 = 1.0;
110  }
111 
114  {
115  sign0 = sign1 = sign2 = sign3 = 1.0;
116  if (s&1) sign0 *= -1.0;
117  if (s&2) sign1 *= -1.0;
118  if (s&4) sign2 *= -1.0;
119  if (s&8) sign3 *= -1.0;
120 
121  m0[0] = 0.0; m0[1] = 0.5;
122  m1[0] = 1.0; m1[1] = 0.5;
123  m2[0] = 0.5; m2[1] = 0.0;
124  m3[0] = 0.5; m3[1] = 1.0;
125 
126  n0[0] = -1.0; n0[1] = 0.0;
127  n1[0] = 1.0; n1[1] = 0.0;
128  n2[0] = 0.0; n2[1] = -1.0;
129  n3[0] = 0.0; n3[1] = 1.0;
130  }
131 
132  template<typename F, typename C>
133  void interpolate (const F& f, std::vector<C>& out) const
134  {
135  // f gives v*outer normal at a point on the edge!
136  typename F::Traits::RangeType y;
137 
138  out.resize(4);
139 
140  f.evaluate(m0,y); out[0] = (y[0]*n0[0]+y[1]*n0[1])*sign0;
141  f.evaluate(m1,y); out[1] = (y[0]*n1[0]+y[1]*n1[1])*sign1;
142  f.evaluate(m2,y); out[2] = (y[0]*n2[0]+y[1]*n2[1])*sign2;
143  f.evaluate(m3,y); out[3] = (y[0]*n3[0]+y[1]*n3[1])*sign3;
144  }
145 
146  private:
147  typename LB::Traits::RangeFieldType sign0,sign1,sign2,sign3;
148  typename LB::Traits::DomainType m0,m1,m2,m3;
149  typename LB::Traits::DomainType n0,n1,n2,n3;
150  };
151 
159  {
160  public:
163  {
164  for (std::size_t i=0; i<4; i++)
165  li[i] = LocalKey(i,1,0);
166  }
167 
169  std::size_t size () const
170  {
171  return 4;
172  }
173 
175  const LocalKey& localKey (std::size_t i) const
176  {
177  return li[i];
178  }
179 
180  private:
181  std::vector<LocalKey> li;
182  };
183 
184 }
185 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_CUBE2D_ALL_HH
void interpolate(const F &f, std::vector< C > &out) const
Definition: raviartthomas0cube2dall.hh:133
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas0cube2dall.hh:84
Lowest order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas0cube2dall.hh:102
Layout map for RT0 elements on quadrilaterals.
Definition: raviartthomas0cube2dall.hh:158
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition: raviartthomas0cube2dall.hh:29
RT0Cube2DLocalInterpolation()
Standard constructor.
Definition: raviartthomas0cube2dall.hh:107
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:14
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas0cube2dall.hh:54
Lowest order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas0cube2dall.hh:25
unsigned int size() const
number of shape functions
Definition: raviartthomas0cube2dall.hh:48
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas0cube2dall.hh:66
RT0Cube2DLocalInterpolation(unsigned int s)
Make set numer s, where 0<=s<8.
Definition: raviartthomas0cube2dall.hh:113
RT0Cube2DLocalCoefficients()
Standard constructor.
Definition: raviartthomas0cube2dall.hh:162
RT0Cube2DLocalBasis(unsigned int s)
Make set numer s, where 0<=s<16.
Definition: raviartthomas0cube2dall.hh:38
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
RT0Cube2DLocalBasis()
Standard constructor.
Definition: raviartthomas0cube2dall.hh:32
D DomainType
domain type
Definition: localbasis.hh:49
std::size_t size() const
number of coefficients
Definition: raviartthomas0cube2dall.hh:169
Describe position of one degree of freedom.
Definition: localkey.hh:21
const LocalKey & localKey(std::size_t i) const
get i'th index
Definition: raviartthomas0cube2dall.hh:175