dune-localfunctions  2.4.1-rc2
dualq1.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_Q1_LOCALFINITEELEMENT_HH
4 #define DUNE_DUAL_Q1_LOCALFINITEELEMENT_HH
5 
6 #include <array>
7 
8 #include <dune/common/fvector.hh>
9 #include <dune/common/fmatrix.hh>
10 
11 #include <dune/geometry/type.hh>
12 #include <dune/geometry/quadraturerules.hh>
13 
18 
19 namespace Dune
20 {
21 
31  template<class D, class R, int dim>
33  {
34  public:
39 
43  {
44  gt.makeCube(dim);
45 
46  // dual basis functions are linear combinations of Lagrange elements
47  // compute these coefficients here because the basis and the local interpolation needs them
48 
49  const Dune::QuadratureRule<D,dim>& quad = Dune::QuadratureRules<D,dim>::rule(gt, 2*dim);
50 
51  const int size = 1 <<dim;
52 
53  // assemble mass matrix on the reference element
54  Dune::FieldMatrix<R, size, size> massMat;
55  massMat = 0;
56 
57  // and the integrals of the lagrange shape functions
58  std::vector<Dune::FieldVector<R,1> > integral(size);
59  for (int i=0; i<size; i++)
60  integral[i] = 0;
61 
62  for(size_t pt=0; pt<quad.size(); pt++) {
63 
64  const Dune::FieldVector<D ,dim>& pos = quad[pt].position();
65  std::vector<typename Traits::LocalBasisType::Traits::RangeType> q1Values(size);
66 
67  // evaluate q1 basis functions
68  for (int i=0; i<size; i++) {
69 
70  q1Values[i] = 1;
71 
72  for (int j=0; j<dim; j++)
73  // if j-th bit of i is set multiply with in[j], else with 1-in[j]
74  q1Values[i] *= (i & (1<<j)) ? pos[j] : 1-pos[j];
75  }
76 
77  double weight = quad[pt].weight();
78 
79  for (int k=0; k<size; k++) {
80  integral[k] += q1Values[k]*weight;
81 
82  for (int l=0; l<=k; l++)
83  massMat[k][l] += weight*(q1Values[k]*q1Values[l]);
84  }
85  }
86 
87  // make matrix symmetric
88  for (int i=0; i<size-1; i++)
89  for (int j=i+1; j<size; j++)
90  massMat[i][j] = massMat[j][i];
91 
92  //solve for the coefficients
93  std::array<Dune::FieldVector<R, size>, size> coefficients;
94 
95  for (int i=0; i<size; i++) {
96 
97  Dune::FieldVector<R, size> rhs(0);
98  rhs[i] = integral[i];
99 
100  coefficients[i] = 0;
101  massMat.solve(coefficients[i] ,rhs);
102 
103  }
104 
105  basis.setCoefficients(coefficients);
106  interpolation.setCoefficients(coefficients);
107  }
108 
111  const typename Traits::LocalBasisType& localBasis () const
112  {
113  return basis;
114  }
115 
119  {
120  return coefficients;
121  }
122 
126  {
127  return interpolation;
128  }
129 
131  unsigned int size () const
132  {
133  return basis.size();
134  }
135 
138  GeometryType type () const
139  {
140  return gt;
141  }
142 
144  {
145  return new DualQ1LocalFiniteElement(*this);
146  }
147 
148  private:
150  DualQ1LocalCoefficients<dim> coefficients;
152  GeometryType gt;
153  };
154 }
155 
156 #endif
traits helper struct
Definition: localfiniteelementtraits.hh:10
DualQ1LocalFiniteElement()
Definition: dualq1.hh:42
Layout map for dual Q1 elements.
Definition: dualq1localcoefficients.hh:22
const Traits::LocalInterpolationType & localInterpolation() const
Definition: dualq1.hh:125
DualQ1LocalFiniteElement * clone() const
Definition: dualq1.hh:143
Dual Lagrange shape functions of order 1 on the reference cube.
Definition: dualq1localbasis.hh:25
const Traits::LocalBasisType & localBasis() const
Definition: dualq1.hh:111
const Traits::LocalCoefficientsType & localCoefficients() const
Definition: dualq1.hh:118
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:14
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:22
LocalFiniteElementTraits< DualQ1LocalBasis< D, R, dim >, DualQ1LocalCoefficients< dim >, DualQ1LocalInterpolation< dim, DualQ1LocalBasis< D, R, dim > > > Traits
Definition: dualq1.hh:38
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:18
Definition: dualq1localinterpolation.hh:17
unsigned int size() const
Number of shape functions in this finite element.
Definition: dualq1.hh:131
GeometryType type() const
Definition: dualq1.hh:138
LB LocalBasisType
Definition: localfiniteelementtraits.hh:14
The local dual Q1 finite element on cubes.
Definition: dualq1.hh:32