4 #ifndef DUNE_LOCALFUNCTIONS_QKLOCALBASIS_HH
5 #define DUNE_LOCALFUNCTIONS_QKLOCALBASIS_HH
9 #include <dune/common/fvector.hh>
10 #include <dune/common/fmatrix.hh>
11 #include <dune/common/power.hh>
13 #include <dune/geometry/type.hh>
33 template<
class D,
class R,
int k,
int d>
37 static R p (
int i, D x)
40 for (
int j=0; j<=k; j++)
41 if (j!=i) result *= (k*x-j)/(i-j);
46 static R dp (
int i, D x)
50 for (
int j=0; j<=k; j++)
53 R prod( (k*1.0)/(i-j) );
54 for (
int l=0; l<=k; l++)
56 prod *= (k*x-l)/(i-l);
64 static R ddp (
int j, D x)
68 for (
int i=0; i<=k; i++)
75 for (
int m=0; m<=k; m++)
80 R prod( (k*1.0)/(j-m) );
81 for (
int l=0; l<=k; l++)
82 if (l!=i && l!=j && l!=m)
83 prod *= (k*x-l)/(j-l);
87 result += sum * ( (k*1.0)/(j-i) );
94 static Dune::FieldVector<int,d> multiindex (
int i)
96 Dune::FieldVector<int,d> alpha;
97 for (
int j=0; j<d; j++)
106 typedef LocalBasisTraits<D,d,Dune::FieldVector<D,d>,R,1,Dune::FieldVector<R,1>,Dune::FieldMatrix<R,1,d> >
Traits;
111 return StaticPower<k+1,d>::power;
116 std::vector<typename Traits::RangeType>& out)
const
119 for (
size_t i=0; i<
size(); i++)
122 Dune::FieldVector<int,d> alpha(multiindex(i));
128 for (
int j=0; j<d; j++)
129 out[i] *= p(alpha[j],in[j]);
139 std::vector<typename Traits::JacobianType>& out)
const
144 for (
size_t i=0; i<
size(); i++)
147 Dune::FieldVector<int,d> alpha(multiindex(i));
150 for (
int j=0; j<d; j++)
154 out[i][0][j] = dp(alpha[j],in[j]);
157 for (
int l=0; l<d; l++)
159 out[i][0][j] *= p(alpha[l],in[l]);
171 std::vector<typename Traits::RangeType>& out)
const
176 for (
size_t i=0; i<
size(); i++)
179 Dune::FieldVector<int,d> alpha(multiindex(i));
185 for (std::size_t l=0; l<d; l++)
190 out[i][0] *= p(alpha[l],in[l]);
193 out[i][0] *= dp(alpha[l],in[l]);
196 out[i][0] *= ddp(alpha[l],in[l]);
199 DUNE_THROW(NotImplemented,
"Desired derivative order is not implemented");