dune-localfunctions  2.4
qklocalbasis.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 
4 #ifndef DUNE_LOCALFUNCTIONS_QKLOCALBASIS_HH
5 #define DUNE_LOCALFUNCTIONS_QKLOCALBASIS_HH
6 
7 #include <dune/common/fvector.hh>
8 #include <dune/common/fmatrix.hh>
9 #include <dune/common/power.hh>
10 
11 #include <dune/geometry/type.hh>
12 
15 
16 
17 namespace Dune
18 {
31  template<class D, class R, int k, int d>
33  {
34  enum { n = StaticPower<k+1,d>::power };
35 
36  // ith Lagrange polynomial of degree k in one dimension
37  static R p (int i, D x)
38  {
39  R result(1.0);
40  for (int j=0; j<=k; j++)
41  if (j!=i) result *= (k*x-j)/(i-j);
42  return result;
43  }
44 
45  // derivative of ith Lagrange polynomial of degree k in one dimension
46  static R dp (int i, D x)
47  {
48  R result(0.0);
49 
50  for (int j=0; j<=k; j++)
51  if (j!=i)
52  {
53  R prod( (k*1.0)/(i-j) );
54  for (int l=0; l<=k; l++)
55  if (l!=i && l!=j)
56  prod *= (k*x-l)/(i-l);
57  result += prod;
58  }
59  return result;
60  }
61 
62  // Return i as a d-digit number in the (k+1)-nary system
63  static Dune::FieldVector<int,d> multiindex (int i)
64  {
65  Dune::FieldVector<int,d> alpha;
66  for (int j=0; j<d; j++)
67  {
68  alpha[j] = i % (k+1);
69  i = i/(k+1);
70  }
71  return alpha;
72  }
73 
74  public:
75  typedef LocalBasisTraits<D,d,Dune::FieldVector<D,d>,R,1,Dune::FieldVector<R,1>,Dune::FieldMatrix<R,1,d>, 1> Traits;
76 
78  unsigned int size () const
79  {
80  return StaticPower<k+1,d>::power;
81  }
82 
84  inline void evaluateFunction (const typename Traits::DomainType& in,
85  std::vector<typename Traits::RangeType>& out) const
86  {
87  out.resize(size());
88  for (size_t i=0; i<size(); i++)
89  {
90  // convert index i to multiindex
91  Dune::FieldVector<int,d> alpha(multiindex(i));
92 
93  // initialize product
94  out[i] = 1.0;
95 
96  // dimension by dimension
97  for (int j=0; j<d; j++)
98  out[i] *= p(alpha[j],in[j]);
99  }
100  }
101 
106  inline void
108  std::vector<typename Traits::JacobianType>& out) const
109  {
110  out.resize(size());
111 
112  // Loop over all shape functions
113  for (size_t i=0; i<size(); i++)
114  {
115  // convert index i to multiindex
116  Dune::FieldVector<int,d> alpha(multiindex(i));
117 
118  // Loop over all coordinate directions
119  for (int j=0; j<d; j++)
120  {
121  // Initialize: the overall expression is a product
122  // if j-th bit of i is set to -1, else 1
123  out[i][0][j] = dp(alpha[j],in[j]);
124 
125  // rest of the product
126  for (int l=0; l<d; l++)
127  if (l!=j)
128  out[i][0][j] *= p(alpha[l],in[l]);
129  }
130  }
131  }
132 
138  template<int diffOrder>
139  inline void evaluate(
140  const std::array<int,1>& direction,
141  const typename Traits::DomainType& in,
142  std::vector<typename Traits::RangeType>& out) const
143  {
144  static_assert(diffOrder == 1, "We only can compute first derivatives");
145  out.resize(size());
146 
147  // Loop over all shape functions
148  for (size_t i=0; i<size(); i++)
149  {
150  // convert index i to multiindex
151  Dune::FieldVector<int,d> alpha(multiindex(i));
152 
153  // Loop over all coordinate directions
154  std::size_t j = direction[0];
155 
156  // Initialize: the overall expression is a product
157  // if j-th bit of i is set to -1, else 1
158  out[i][0] = dp(alpha[j],in[j]);
159 
160  // rest of the product
161  for (std::size_t l=0; l<d; l++)
162  if (l!=j)
163  out[i][0] *= p(alpha[l],in[l]);
164  }
165  }
166 
168  unsigned int order () const
169  {
170  return k;
171  }
172  };
173 }
174 
175 #endif
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: qklocalbasis.hh:107
LocalBasisTraits< D, d, Dune::FieldVector< D, d >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, d >, 1 > Traits
Definition: qklocalbasis.hh:75
void evaluate(const std::array< int, 1 > &direction, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate derivative in a given direction.
Definition: qklocalbasis.hh:139
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
unsigned int size() const
number of shape functions
Definition: qklocalbasis.hh:78
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:14
Lagrange shape functions of order k on the reference cube.
Definition: qklocalbasis.hh:32
unsigned int order() const
Polynomial order of the shape functions.
Definition: qklocalbasis.hh:168
D DomainType
domain type
Definition: localbasis.hh:49
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: qklocalbasis.hh:84