dune-localfunctions  2.4.1-rc2
dualp1localbasis.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_DUAL_P1_LOCALBASIS_HH
4 #define DUNE_DUAL_P1_LOCALBASIS_HH
5 
6 #include <dune/common/fvector.hh>
7 #include <dune/common/fmatrix.hh>
9 
10 namespace Dune
11 {
28  template<class D, class R, int dim, bool faceDualT=false>
30  {
31  public:
33  static const bool faceDual = faceDualT;
35  typedef LocalBasisTraits<D,dim,Dune::FieldVector<D,dim>,R,1,Dune::FieldVector<R,1>,
36  Dune::FieldMatrix<R,1,dim> > Traits;
37 
39  unsigned int size () const
40  {
41  return dim+1;
42  }
43 
45  inline void evaluateFunction (const typename Traits::DomainType& in,
46  std::vector<typename Traits::RangeType>& out) const
47  {
48  // evaluate P1 basis functions
49  std::vector<typename Traits::RangeType> p1Values(size());
50 
51  p1Values[0] = 1.0;
52 
53  for (int i=0; i<dim; i++) {
54  p1Values[0] -= in[i];
55  p1Values[i+1] = in[i];
56  }
57 
58  // compute dual basis function values as a linear combination of the Lagrange values
59  out.resize(size());
60 
61  for (int i=0; i<=dim; i++) {
62  out[i] = (dim+!faceDual)*p1Values[i];
63  for (int j=0; j<i; j++)
64  out[i] -= p1Values[j];
65 
66  for (int j=i+1; j<=dim; j++)
67  out[i] -= p1Values[j];
68  }
69  }
70 
72  inline void
73  evaluateJacobian (const typename Traits::DomainType& in,
74  std::vector<typename Traits::JacobianType>& out) const
75  {
76  // evaluate P1 jacobians
77  std::vector<typename Traits::JacobianType> p1Jacs(size());
78 
79  for (int i=0; i<dim; i++)
80  p1Jacs[0][0][i] = -1;
81 
82  for (int i=0; i<dim; i++)
83  for (int j=0; j<dim; j++)
84  p1Jacs[i+1][0][j] = (i==j);
85 
86  // compute dual basis jacobians as linear combination of the Lagrange jacobians
87  out.resize(size());
88 
89  for (size_t i=0; i<=dim; i++) {
90  out[i][0] = 0;
91  out[i][0].axpy(dim+!faceDual,p1Jacs[i][0]);
92 
93  for (size_t j=0; j<i; j++)
94  out[i][0] -= p1Jacs[j][0];
95 
96  for (int j=i+1; j<=dim; j++)
97  out[i][0] -= p1Jacs[j][0];
98  }
99  }
100 
102  unsigned int order () const
103  {
104  return 1;
105  }
106  };
107 }
108 #endif
Dual Lagrange shape functions on the simplex.
Definition: dualp1localbasis.hh:29
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: dualp1localbasis.hh:73
unsigned int size() const
number of shape functions
Definition: dualp1localbasis.hh:39
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: dualp1localbasis.hh:45
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim > > Traits
export type traits for function signature
Definition: dualp1localbasis.hh:36
unsigned int order() const
Polynomial order of the shape functions.
Definition: dualp1localbasis.hh:102
D DomainType
domain type
Definition: localbasis.hh:49
static const bool faceDual
Determines if the basis is only biorthogonal on adjacent faces.
Definition: dualp1localbasis.hh:33