4#ifndef DUNE_LOCALFUNCTIONS_SERENDIPITYLOCALBASIS_HH
5#define DUNE_LOCALFUNCTIONS_SERENDIPITYLOCALBASIS_HH
7#include <dune/common/fvector.hh>
8#include <dune/common/fmatrix.hh>
10#include <dune/geometry/type.hh>
12#include <dune/localfunctions/common/localbasis.hh>
13#include <dune/localfunctions/common/localfiniteelementtraits.hh>
28 template<
class D,
class R,
int k,
int d>
31 typedef typename Dune::FieldVector<int,d> DuneVector;
35 static R p (
int i,
int deg, D x)
38 for (
int j=0; j<=deg; j++)
39 if (j!=i) result *= (deg*x-j)/(i-j);
44 static R dp (
int i,
int deg, D x)
48 for (
int j=0; j<=deg; j++)
51 R prod( (deg*1.0)/(i-j) );
52 for (
int l=0; l<=deg; l++)
54 prod *= (deg*x-l)/(i-l);
60 static void multiindex (DuneVector& alpha)
63 for (
int j=0; j<d; j++){
64 i+=alpha[j]*pow(k+1,j);
67 for (
int j=0; j<d; j++)
73 for (
unsigned int j=0; j<d; j++){
74 if(alpha[j] > 0 && alpha[j] < k){ sum++; }
81 bool isEdge(DuneVector& alpha)
const{
82 for(
int i = 0; i < d; i++){
83 if(alpha[i]>0 && alpha[i] < k)
89 double edgeVal(
int alpha,
const double& in,
int der)
const{
91 if(alpha > 0 && alpha < k){
95 return dp(alpha,k,in);
99 return p(alpha/k,k-1,in);
101 return dp(alpha/k,k-1,in);
105 double nodeVal(
int alpha,
const double& in,
int der)
const{
108 return p(alpha/k,k-1,in);
110 return dp(alpha/k,k-1,in);
114 typedef LocalBasisTraits<D,d,Dune::FieldVector<D,d>,R,1,Dune::FieldVector<R,1>,Dune::FieldMatrix<R,1,d> > Traits;
132 std::vector<typename Traits::RangeType>& out)
const
136 for (
size_t i=0; i<
size(); i++)
139 if(i>0){ multiindex(alpha); }
142 for(
int j = 0; j < d; j++){
144 out[i] *= edgeVal(alpha[j], in[j], 0);
146 out[i] *= nodeVal(alpha[j], in[j], 0);
151 out[0] -= 0.5*(out[1]+out[3]);
152 out[2] -= 0.5*(out[1]+out[4]);
153 out[5] -= 0.5*(out[3]+out[6]);
154 out[7] -= 0.5*(out[4]+out[6]);
157 out[0] -= 0.5*(out[1] +out[3] +out[8]);
158 out[2] -= 0.5*(out[1] +out[4] +out[9]);
159 out[5] -= 0.5*(out[3] +out[6] +out[10]);
160 out[7] -= 0.5*(out[4] +out[6] +out[11]);
161 out[12] -= 0.5*(out[8] +out[13]+out[15]);
162 out[14] -= 0.5*(out[9] +out[13]+out[16]);
163 out[17] -= 0.5*(out[10]+out[15]+out[18]);
164 out[19] -= 0.5*(out[11]+out[16]+out[18]);
208 std::vector<typename Traits::JacobianType>& out)
const
214 for (
size_t i=0; i<
size(); i++)
217 if(i>0){ multiindex(alpha); }
221 for (
int j=0; j<d; j++)
226 out[i][0][j] *= edgeVal(alpha[j],in[j], 1);
228 out[i][0][j] *= nodeVal(alpha[j],in[j], 1);
230 for (
int l=0; l<d; l++){
233 out[i][0][j] *= edgeVal(alpha[l],in[l],0);
235 out[i][0][j] *= nodeVal(alpha[l],in[l],0);
241 for(
int j = 0; j < d; j++){
242 out[0][0][j] -= 0.5*(out[1][0][j]+out[3][0][j]);
243 out[2][0][j] -= 0.5*(out[1][0][j]+out[4][0][j]);
244 out[5][0][j] -= 0.5*(out[3][0][j]+out[6][0][j]);
245 out[7][0][j] -= 0.5*(out[4][0][j]+out[6][0][j]);
249 for(
int j = 0; j < d; j++){
250 out[0][0][j] -= 0.5*(out[1][0][j] +out[3][0][j] +out[8][0][j]);
251 out[2][0][j] -= 0.5*(out[1][0][j] +out[4][0][j] +out[9][0][j]);
252 out[5][0][j] -= 0.5*(out[3][0][j] +out[6][0][j] +out[10][0][j]);
253 out[7][0][j] -= 0.5*(out[4][0][j] +out[6][0][j] +out[11][0][j]);
254 out[12][0][j] -= 0.5*(out[8][0][j] +out[13][0][j]+out[15][0][j]);
255 out[14][0][j] -= 0.5*(out[9][0][j] +out[13][0][j]+out[16][0][j]);
256 out[17][0][j] -= 0.5*(out[10][0][j]+out[15][0][j]+out[18][0][j]);
257 out[19][0][j] -= 0.5*(out[11][0][j]+out[16][0][j]+out[18][0][j]);
267 template<
int diffOrder>
269 const std::array<int,1>& direction,
270 const typename Traits::DomainType& in,
271 std::vector<typename Traits::RangeType>& out)
const
273 static_assert(diffOrder == 1,
"We only can compute first derivatives");
278 for (
size_t i=0; i<
size(); i++)
281 if(i>0){ multiindex(alpha);}
284 std::size_t j = direction[0];
290 out[i][0] *= edgeVal(alpha[j],in[j],1);
292 out[i][0] *= nodeVal(alpha[j],in[j],1);
295 for (std::size_t l=0; l<d; l++){
298 out[i][0] *= edgeVal(alpha[l],in[l],0);
300 out[i][0] *= nodeVal(alpha[l],in[l],0);
305 out[0][0] -= 0.5*(out[1][0]+out[3][0]);
306 out[2][0] -= 0.5*(out[1][0]+out[4][0]);
307 out[5][0] -= 0.5*(out[3][0]+out[6][0]);
308 out[7][0] -= 0.5*(out[4][0]+out[6][0]);
311 out[0][0] -= 0.5*(out[1][0] +out[3][0] +out[8][0]);
312 out[2][0] -= 0.5*(out[1][0] +out[4][0] +out[9][0]);
313 out[5][0] -= 0.5*(out[3][0] +out[6][0] +out[10][0]);
314 out[7][0] -= 0.5*(out[4][0] +out[6][0] +out[11][0]);
315 out[12][0] -= 0.5*(out[8][0] +out[13][0]+out[15][0]);
316 out[14][0] -= 0.5*(out[9][0] +out[13][0]+out[16][0]);
317 out[17][0] -= 0.5*(out[10][0]+out[15][0]+out[18][0]);
318 out[19][0] -= 0.5*(out[11][0]+out[16][0]+out[18][0]);
Serendipity basis functions of order k on the reference cube.
Definition: serendipitylocalbasis.hh:30
unsigned int order() const
Polynomial order of the shape functions.
Definition: serendipitylocalbasis.hh:323
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: serendipitylocalbasis.hh:207
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: serendipitylocalbasis.hh:268
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: serendipitylocalbasis.hh:131
unsigned int size() const
number of shape functions
Definition: serendipitylocalbasis.hh:117