dune-localfunctions  2.2.1
raviartthomas0q3dall.hh
Go to the documentation of this file.
1 #ifndef DUNE_RT0Q3DALL_HH
2 #define DUNE_RT0Q3DALL_HH
3 
4 #include <cstddef>
5 #include <vector>
6 
7 #include <dune/common/fmatrix.hh>
8 
11 
12 namespace Dune
13 {
22  template<class D, class R>
24  {
25  public:
26  typedef LocalBasisTraits<D,3,Dune::FieldVector<D,3>,R,3,Dune::FieldVector<R,3>,
27  Dune::FieldMatrix<R,3,3> > Traits;
28 
31  {
32  sign0 = sign1 = sign2 = sign3 = sign4 = sign5 = 1.0;
33  }
34 
36  RT0Q3DLocalBasis (unsigned int s)
37  {
38  sign0 = sign1 = sign2 = sign3 = sign4 = sign5 = 1.0;
39  if (s&1) sign0 = -1.0;
40  if (s&2) sign1 = -1.0;
41  if (s&4) sign2 = -1.0;
42  if (s&8) sign3 = -1.0;
43  if (s&16) sign4 = -1.0;
44  if (s&32) sign5 = -1.0;
45  }
46 
48  unsigned int size () const
49  {
50  return 6;
51  }
52 
54  inline void evaluateFunction (const typename Traits::DomainType& in,
55  std::vector<typename Traits::RangeType>& out) const
56  {
57  out.resize(6);
58  out[0][0] = sign0*(in[0]-1.0); out[0][1]=0.0; out[0][2]=0.0;
59  out[1][0] = sign1*(in[0]); out[1][1]=0.0; out[1][2]=0.0;
60  out[2][0] = 0.0; out[2][1]=sign2*(in[1]-1.0); out[2][2]=0.0;
61  out[3][0] = 0.0; out[3][1]=sign3*(in[1]); out[3][2]=0.0;
62  out[4][0] = 0.0; out[4][1]=0.0; out[4][2]=sign4*(in[2]-1.0);
63  out[5][0] = 0.0; out[5][1]=0.0; out[5][2]=sign5*(in[2]);
64  }
65 
67  inline void
68  evaluateJacobian (const typename Traits::DomainType& in, // position
69  std::vector<typename Traits::JacobianType>& out) const // return value
70  {
71  out.resize(6);
72  out[0][0][0] = sign0; out[0][0][1] = 0; out[0][0][2] = 0;
73  out[0][1][0] = 0; out[0][1][1] = 0; out[0][1][2] = 0;
74  out[0][2][0] = 0; out[0][2][1] = 0; out[0][2][2] = 0;
75 
76  out[1][0][0] = sign1; out[1][0][1] = 0; out[1][0][2] = 0;
77  out[1][1][0] = 0; out[1][1][1] = 0; out[1][1][2] = 0;
78  out[1][2][0] = 0; out[1][2][1] = 0; out[1][2][2] = 0;
79 
80  out[2][0][0] = 0; out[2][0][1] = 0; out[2][0][2] = 0;
81  out[2][1][0] = 0; out[2][1][1] = sign2; out[2][1][2] = 0;
82  out[2][2][0] = 0; out[2][2][1] = 0; out[2][2][2] = 0;
83 
84  out[3][0][0] = 0; out[3][0][1] = 0; out[3][0][2] = 0;
85  out[3][1][0] = 0; out[3][1][1] = sign3; out[3][1][2] = 0;
86  out[3][2][0] = 0; out[3][2][1] = 0; out[3][2][2] = 0;
87 
88  out[4][0][0] = 0; out[4][0][1] = 0; out[4][0][2] = 0;
89  out[4][1][0] = 0; out[4][1][1] = 0; out[4][1][2] = 0;
90  out[4][2][0] = 0; out[4][2][1] = 0; out[4][2][2] = sign4;
91 
92  out[5][0][0] = 0; out[5][0][1] = 0; out[5][0][2] = 0;
93  out[5][1][0] = 0; out[5][1][1] = 0; out[5][1][2] = 0;
94  out[5][2][0] = 0; out[5][2][1] = 0; out[5][2][2] = sign5;
95  }
96 
98  unsigned int order () const
99  {
100  return 1;
101  }
102 
103  private:
104  R sign0, sign1, sign2, sign3, sign4, sign5;
105  };
106 
107 
115  template<class LB>
117  {
118  public:
119 
122  {
123  sign0 = sign1 = sign2 = sign3 = sign4 = sign5 = 1.0;
124  }
125 
127  RT0Q3DLocalInterpolation (unsigned int s)
128  {
129  sign0 = sign1 = sign2 = sign3 = sign4 = sign5 = 1.0;
130  if (s&1) sign0 *= -1.0;
131  if (s&2) sign1 *= -1.0;
132  if (s&4) sign2 *= -1.0;
133  if (s&8) sign3 *= -1.0;
134  if (s&16) sign4 *= -1.0;
135  if (s&32) sign5 *= -1.0;
136 
137  m0[0] = 0.0; m0[1] = 0.5; m0[2] = 0.5;
138  m1[0] = 1.0; m1[1] = 0.5; m1[2] = 0.5;
139  m2[0] = 0.5; m2[1] = 0.0; m2[2] = 0.5;
140  m3[0] = 0.5; m3[1] = 1.0; m3[2] = 0.5;
141  m4[0] = 0.5; m4[1] = 0.5; m4[2] = 0.0;
142  m5[0] = 0.5; m5[1] = 0.5; m5[2] = 1.0;
143 
144  n0[0] = -1.0; n0[1] = 0.0; n0[2] = 0.0;
145  n1[0] = 1.0; n1[1] = 0.0; n1[2] = 0.0;
146  n2[0] = 0.0; n2[1] = -1.0; n2[2] = 0.0;
147  n3[0] = 0.0; n3[1] = 1.0; n3[2] = 0.0;
148  n4[0] = 0.0; n4[1] = 0.0; n4[2] =-1.0;
149  n5[0] = 0.0; n5[1] = 0.0; n5[2] = 1.0;
150  }
151 
152  template<typename F, typename C>
153  void interpolate (const F& f, std::vector<C>& out) const
154  {
155  // f gives v*outer normal at a point on the edge!
156  typename F::Traits::RangeType y;
157 
158  out.resize(6);
159 
160  f.evaluate(m0,y); out[0] = (y[0]*n0[0]+y[1]*n0[1]+y[2]*n0[2])*sign0;
161  f.evaluate(m1,y); out[1] = (y[0]*n1[0]+y[1]*n1[1]+y[2]*n1[2])*sign1;
162  f.evaluate(m2,y); out[2] = (y[0]*n2[0]+y[1]*n2[1]+y[2]*n2[2])*sign2;
163  f.evaluate(m3,y); out[3] = (y[0]*n3[0]+y[1]*n3[1]+y[2]*n3[2])*sign3;
164  f.evaluate(m4,y); out[4] = (y[0]*n4[0]+y[1]*n4[1]+y[2]*n4[2])*sign4;
165  f.evaluate(m5,y); out[5] = (y[0]*n5[0]+y[1]*n5[1]+y[2]*n5[2])*sign5;
166  }
167 
168  private:
169  typename LB::Traits::RangeFieldType sign0,sign1,sign2,sign3,sign4,sign5;
170  typename LB::Traits::DomainType m0,m1,m2,m3,m4,m5;
171  typename LB::Traits::DomainType n0,n1,n2,n3,n4,n5;
172  };
173 
181  {
182  public:
185  {
186  for (std::size_t i=0; i<6; i++)
187  li[i] = LocalKey(i,1,0);
188  }
189 
191  std::size_t size () const
192  {
193  return 6;
194  }
195 
197  const LocalKey& localKey (std::size_t i) const
198  {
199  return li[i];
200  }
201 
202  private:
203  std::vector<LocalKey> li;
204  };
205 
206 }
207 #endif