dune-localfunctions  2.3beta2
pk3dlocalbasis.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_PK3DLOCALBASIS_HH
4 #define DUNE_PK3DLOCALBASIS_HH
5 
6 #include <dune/common/fmatrix.hh>
7 
9 
10 namespace Dune
11 {
24  template<class D, class R, unsigned int k>
26  {
27  public:
28  enum {N = (k+1)*(k+2)*(k+3)/6};
29  enum {O = k};
30 
31  typedef LocalBasisTraits<D,3,Dune::FieldVector<D,3>,R,1,Dune::FieldVector<R,1>,
32  Dune::FieldMatrix<R,1,3> > Traits;
33 
36 
38  unsigned int size () const
39  {
40  return N;
41  }
42 
44  inline void evaluateFunction (const typename Traits::DomainType& x,
45  std::vector<typename Traits::RangeType>& out) const
46  {
47  out.resize(N);
48  typename Traits::DomainType kx = x;
49  kx *= k;
50  unsigned int n = 0;
51  unsigned int i[4];
52  R factor[4];
53  for (i[2] = 0; i[2] <= k; ++i[2])
54  {
55  factor[2] = 1.0;
56  for (unsigned int j = 0; j < i[2]; ++j)
57  factor[2] *= (kx[2]-j) / (i[2]-j);
58  for (i[1] = 0; i[1] <= k - i[2]; ++i[1])
59  {
60  factor[1] = 1.0;
61  for (unsigned int j = 0; j < i[1]; ++j)
62  factor[1] *= (kx[1]-j) / (i[1]-j);
63  for (i[0] = 0; i[0] <= k - i[1] - i[2]; ++i[0])
64  {
65  factor[0] = 1.0;
66  for (unsigned int j = 0; j < i[0]; ++j)
67  factor[0] *= (kx[0]-j) / (i[0]-j);
68  i[3] = k - i[0] - i[1] - i[2];
69  D kx3 = k - kx[0] - kx[1] - kx[2];
70  factor[3] = 1.0;
71  for (unsigned int j = 0; j < i[3]; ++j)
72  factor[3] *= (kx3-j) / (i[3]-j);
73  out[n++] = factor[0] * factor[1] * factor[2] * factor[3];
74  }
75  }
76  }
77  }
78 
80  inline void
81  evaluateJacobian (const typename Traits::DomainType& x, // position
82  std::vector<typename Traits::JacobianType>& out) const // return value
83  {
84  out.resize(N);
85  typename Traits::DomainType kx = x;
86  kx *= k;
87  unsigned int n = 0;
88  unsigned int i[4];
89  R factor[4];
90  for (i[2] = 0; i[2] <= k; ++i[2])
91  {
92  factor[2] = 1.0;
93  for (unsigned int j = 0; j < i[2]; ++j)
94  factor[2] *= (kx[2]-j) / (i[2]-j);
95  for (i[1] = 0; i[1] <= k - i[2]; ++i[1])
96  {
97  factor[1] = 1.0;
98  for (unsigned int j = 0; j < i[1]; ++j)
99  factor[1] *= (kx[1]-j) / (i[1]-j);
100  for (i[0] = 0; i[0] <= k - i[1] - i[2]; ++i[0])
101  {
102  factor[0] = 1.0;
103  for (unsigned int j = 0; j < i[0]; ++j)
104  factor[0] *= (kx[0]-j) / (i[0]-j);
105  i[3] = k - i[0] - i[1] - i[2];
106  D kx3 = k - kx[0] - kx[1] - kx[2];
107  R sum3 = 0.0;
108  factor[3] = 1.0;
109  for (unsigned int j = 0; j < i[3]; ++j)
110  factor[3] /= i[3] - j;
111  R prod_all = factor[0] * factor[1] * factor[2] * factor[3];
112  for (unsigned int j = 0; j < i[3]; ++j)
113  {
114  R prod = prod_all;
115  for (unsigned int l = 0; l < i[3]; ++l)
116  if (j == l)
117  prod *= -R(k);
118  else
119  prod *= kx3 - l;
120  sum3 += prod;
121  }
122  for (unsigned int j = 0; j < i[3]; ++j)
123  factor[3] *= kx3 - j;
124  for (unsigned int m = 0; m < 3; ++m)
125  {
126  out[n][0][m] = sum3;
127  for (unsigned int j = 0; j < i[m]; ++j)
128  {
129  R prod = factor[3];
130  for (unsigned int p = 0; p < 3; ++p)
131  {
132  if (m == p)
133  for (unsigned int l = 0; l < i[p]; ++l)
134  if (j == l)
135  prod *= R(k) / (i[p]-l);
136  else
137  prod *= (kx[p]-l) / (i[p]-l);
138  else
139  prod *= factor[p];
140  }
141  out[n][0][m] += prod;
142  }
143  }
144  n++;
145  }
146  }
147  }
148  }
149 
151  unsigned int order () const
152  {
153  return k;
154  }
155  };
156 
157 
158  //Specialization for k=0
159  template<class D, class R>
160  class Pk3DLocalBasis<D,R,0>
161  {
162  public:
163  typedef LocalBasisTraits<D,3,Dune::FieldVector<D,3>,R,1,Dune::FieldVector<R,1>,
164  Dune::FieldMatrix<R,1,3> > Traits;
165 
167  enum {N = 1};
168 
170  enum {O = 0};
171 
172  unsigned int size () const
173  {
174  return 1;
175  }
176 
177  inline void evaluateFunction (const typename Traits::DomainType& in,
178  std::vector<typename Traits::RangeType>& out) const
179  {
180  out.resize(1);
181  out[0] = 1;
182  }
183 
184  // evaluate derivative of a single component
185  inline void
186  evaluateJacobian (const typename Traits::DomainType& in, // position
187  std::vector<typename Traits::JacobianType>& out) const // return value
188  {
189  out.resize(1);
190  out[0][0][0] = 0;
191  out[0][0][1] = 0;
192  out[0][0][2] = 0;
193  }
194 
195  // local interpolation of a function
196  template<typename E, typename F, typename C>
197  void interpolate (const E& e, const F& f, std::vector<C>& out) const
198  {
199  typename Traits::DomainType x;
200  typename Traits::RangeType y;
201  x[0] = 1.0/4.0;
202  x[1] = 1.0/4.0;
203  x[2] = 1.0/4.0;
204  f.eval_local(e,x,y);
205  out[0] = y;
206  }
207 
208  unsigned int order () const
209  {
210  return 0;
211  }
212  };
213 }
214 #endif