dune-localfunctions  2.3.1-rc1
pk2dlocalbasis.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_PK2DLOCALBASIS_HH
4 #define DUNE_PK2DLOCALBASIS_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 
30  enum {N = (k+1)*(k+2)/2};
31 
35  enum {O = k};
36 
37  typedef LocalBasisTraits<D,2,Dune::FieldVector<D,2>,R,1,Dune::FieldVector<R,1>,
38  Dune::FieldMatrix<R,1,2> > Traits;
39 
42  {
43  for (unsigned int i=0; i<=k; i++)
44  pos[i] = (1.0*i)/std::max(k,(unsigned int)1);
45  }
46 
48  unsigned int size () const
49  {
50  return N;
51  }
52 
54  inline void evaluateFunction (const typename Traits::DomainType& x,
55  std::vector<typename Traits::RangeType>& out) const
56  {
57  out.resize(N);
58  // specialization for k==0, not clear whether that is needed
59  if (k==0) {
60  out[0] = 1;
61  return;
62  }
63 
64  int n=0;
65  for (unsigned int j=0; j<=k; j++)
66  for (unsigned int i=0; i<=k-j; i++)
67  {
68  out[n] = 1.0;
69  for (unsigned int alpha=0; alpha<i; alpha++)
70  out[n] *= (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
71  for (unsigned int beta=0; beta<j; beta++)
72  out[n] *= (x[1]-pos[beta])/(pos[j]-pos[beta]);
73  for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
74  out[n] *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
75  n++;
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 
86  // specialization for k==0, not clear whether that is needed
87  if (k==0) {
88  out[0][0][0] = 0; out[0][0][1] = 0;
89  return;
90  }
91 
92  int n=0;
93  for (unsigned int j=0; j<=k; j++)
94  for (unsigned int i=0; i<=k-j; i++)
95  {
96  // x_0 derivative
97  out[n][0][0] = 0.0;
98  R factor=1.0;
99  for (unsigned int beta=0; beta<j; beta++)
100  factor *= (x[1]-pos[beta])/(pos[j]-pos[beta]);
101  for (unsigned int a=0; a<i; a++)
102  {
103  R product=factor;
104  for (unsigned int alpha=0; alpha<i; alpha++)
105  if (alpha==a)
106  product *= 1.0/(pos[i]-pos[alpha]);
107  else
108  product *= (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
109  for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
110  product *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
111  out[n][0][0] += product;
112  }
113  for (unsigned int c=i+j+1; c<=k; c++)
114  {
115  R product=factor;
116  for (unsigned int alpha=0; alpha<i; alpha++)
117  product *= (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
118  for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
119  if (gamma==c)
120  product *= -1.0/(pos[gamma]-pos[i]-pos[j]);
121  else
122  product *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
123  out[n][0][0] += product;
124  }
125 
126  // x_1 derivative
127  out[n][0][1] = 0.0;
128  factor = 1.0;
129  for (unsigned int alpha=0; alpha<i; alpha++)
130  factor *= (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
131  for (unsigned int b=0; b<j; b++)
132  {
133  R product=factor;
134  for (unsigned int beta=0; beta<j; beta++)
135  if (beta==b)
136  product *= 1.0/(pos[j]-pos[beta]);
137  else
138  product *= (x[1]-pos[beta])/(pos[j]-pos[beta]);
139  for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
140  product *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
141  out[n][0][1] += product;
142  }
143  for (unsigned int c=i+j+1; c<=k; c++)
144  {
145  R product=factor;
146  for (unsigned int beta=0; beta<j; beta++)
147  product *= (x[1]-pos[beta])/(pos[j]-pos[beta]);
148  for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
149  if (gamma==c)
150  product *= -1.0/(pos[gamma]-pos[i]-pos[j]);
151  else
152  product *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
153  out[n][0][1] += product;
154  }
155 
156  n++;
157  }
158 
159  }
160 
162  unsigned int order () const
163  {
164  return k;
165  }
166 
167  private:
168  R pos[k+1]; // positions on the interval
169  };
170 
171 }
172 #endif
Definition: pk2dlocalbasis.hh:35
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:39
Definition: pk2dlocalbasis.hh:30
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, 2 > > Traits
Definition: pk2dlocalbasis.hh:38
void evaluateFunction(const typename Traits::DomainType &x, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: pk2dlocalbasis.hh:54
unsigned int size() const
number of shape functions
Definition: pk2dlocalbasis.hh:48
Lagrange shape functions of arbitrary order on the reference triangle.
Definition: pk2dlocalbasis.hh:25
Pk2DLocalBasis()
Standard constructor.
Definition: pk2dlocalbasis.hh:41
D DomainType
domain type
Definition: localbasis.hh:51
void evaluateJacobian(const typename Traits::DomainType &x, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: pk2dlocalbasis.hh:81
unsigned int order() const
Polynomial order of the shape functions.
Definition: pk2dlocalbasis.hh:162