3 #ifndef DUNE_DUAL_Q1_LOCALFINITEELEMENT_HH 4 #define DUNE_DUAL_Q1_LOCALFINITEELEMENT_HH 8 #include <dune/common/fvector.hh> 9 #include <dune/common/fmatrix.hh> 11 #include <dune/geometry/type.hh> 12 #include <dune/geometry/quadraturerules.hh> 31 template<
class D,
class R,
int dim>
49 const Dune::QuadratureRule<D,dim>& quad = Dune::QuadratureRules<D,dim>::rule(gt, 2*dim);
51 const int size = 1 <<dim;
54 Dune::FieldMatrix<R, size, size> massMat;
58 std::vector<Dune::FieldVector<R,1> > integral(size);
59 for (
int i=0; i<
size; i++)
62 for(
size_t pt=0; pt<quad.size(); pt++) {
64 const Dune::FieldVector<D ,dim>& pos = quad[pt].position();
65 std::vector<typename Traits::LocalBasisType::Traits::RangeType> q1Values(size);
68 for (
int i=0; i<
size; i++) {
72 for (
int j=0; j<dim; j++)
74 q1Values[i] *= (i & (1<<j)) ? pos[j] : 1-pos[j];
77 double weight = quad[pt].weight();
79 for (
int k=0; k<
size; k++) {
80 integral[k] += q1Values[k]*weight;
82 for (
int l=0; l<=k; l++)
83 massMat[k][l] += weight*(q1Values[k]*q1Values[l]);
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];
93 std::array<Dune::FieldVector<R, size>, size> coefficients;
95 for (
int i=0; i<
size; i++) {
97 Dune::FieldVector<R, size> rhs(0);
101 massMat.solve(coefficients[i] ,rhs);
105 basis.setCoefficients(coefficients);
106 interpolation.setCoefficients(coefficients);
127 return interpolation;
150 DualQ1LocalCoefficients<dim> coefficients;
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