3#ifndef DUNE_DUAL_Q1_LOCALFINITEELEMENT_HH
4#define DUNE_DUAL_Q1_LOCALFINITEELEMENT_HH
12#include <dune/geometry/referenceelements.hh>
15#include <dune/localfunctions/common/localfiniteelementtraits.hh>
16#include <dune/localfunctions/lagrange/lagrangecube.hh>
17#include "dualq1/dualq1localbasis.hh"
18#include "dualq1/dualq1localcoefficients.hh"
19#include "dualq1/dualq1localinterpolation.hh"
38 template<
class D,
class R,
int dim,
bool faceDual=false>
52 setupFaceDualCoefficients();
54 setupDualCoefficients();
93 void setupFaceDualCoefficients();
96 void setupDualCoefficients();
103 template<
class D,
class R,
int dim,
bool faceDual>
104 void DualQ1LocalFiniteElement<D,R,dim,faceDual>::setupDualCoefficients()
107 const int size = 1 <<dim;
108 std::array<Dune::FieldVector<R, size>, size> coeffs;
119 std::vector<Dune::FieldVector<R,1> > integral(size);
120 for (
int i=0; i<size; i++)
123 Dune::Impl::LagrangeCubeLocalBasis<D,R,dim,1> q1Basis;
124 for(
size_t pt=0; pt<quad.size(); pt++) {
127 std::vector<Dune::FieldVector<R,1> > q1Values(size);
128 q1Basis.evaluateFunction(pos,q1Values);
130 D weight = quad[pt].weight();
132 for (
int k=0; k<size; k++) {
133 integral[k] += q1Values[k]*weight;
135 for (
int l=0; l<=k; l++)
136 massMat[k][l] += weight*(q1Values[k]*q1Values[l]);
141 for (
int i=0; i<size-1; i++)
142 for (
int j=i+1; j<size; j++)
143 massMat[i][j] = massMat[j][i];
147 for (
int i=0; i<size; i++) {
150 rhs[i] = integral[i];
153 massMat.
solve(coeffs[i] ,rhs);
157 basis.setCoefficients(coeffs);
158 interpolation.setCoefficients(coeffs);
161 template<
class D,
class R,
int dim,
bool faceDual>
162 void DualQ1LocalFiniteElement<D,R,dim,faceDual>::setupFaceDualCoefficients()
165 const int size = 1 <<dim;
166 std::array<Dune::FieldVector<R, size>, size> coeffs;
169 Dune::Impl::LagrangeCubeLocalBasis<D,R,dim,1> q1Basis;
174 for (
int i=0; i<refElement.size(1);i++) {
185 const auto& geometry = refElement.template geometry<1>(i);
188 std::vector<Dune::FieldVector<R,1> > integral(size/2);
189 for (
int k=0; k<size/2; k++)
192 for(
size_t pt=0; pt<quad.size(); pt++) {
194 const auto& pos = quad[pt].position();
195 const auto& elementPos = geometry.global(pos);
197 std::vector<Dune::FieldVector<R,1> > q1Values(size);
198 q1Basis.evaluateFunction(elementPos,q1Values);
200 D weight = quad[pt].weight();
202 for (
int k=0; k<refElement.size(i,1,dim); k++) {
203 int row = refElement.subEntity(i,1,k,dim);
204 integral[k] += q1Values[row]*weight;
206 for (
int l=0; l<refElement.size(i,1,dim); l++) {
207 int col = refElement.subEntity(i,1,l,dim);
208 massMat[k][l] += weight*(q1Values[row]*q1Values[col]);
216 for (
int l=0; l<refElement.size(i,1,dim); l++) {
218 int row = refElement.subEntity(i,1,l,dim);
220 rhs[l] = integral[l];
223 massMat.
solve(x ,rhs);
225 for (
int k=0; k<refElement.size(i,1,dim); k++) {
226 int col = refElement.subEntity(i,1,k,dim);
227 coeffs[row][col]=x[k];
232 basis.setCoefficients(coeffs);
233 interpolation.setCoefficients(coeffs);
void solve(V1 &x, const V2 &b, bool doPivoting=true) const
Solve system A x = b.
Dual Lagrange shape functions of order 1 on the reference cube.
Definition: dualq1localbasis.hh:27
Layout map for dual Q1 elements.
Definition: dualq1localcoefficients.hh:23
The local dual Q1 finite element on cubes.
Definition: dualq1.hh:40
LocalFiniteElementTraits< DualQ1LocalBasis< D, R, dim >, DualQ1LocalCoefficients< dim >, DualQ1LocalInterpolation< dim, DualQ1LocalBasis< D, R, dim > > > Traits
Definition: dualq1.hh:45
unsigned int size() const
Number of shape functions in this finite element.
Definition: dualq1.hh:79
const Traits::LocalInterpolationType & localInterpolation() const
Definition: dualq1.hh:73
const Traits::LocalBasisType & localBasis() const
Definition: dualq1.hh:59
const Traits::LocalCoefficientsType & localCoefficients() const
Definition: dualq1.hh:66
DualQ1LocalFiniteElement()
Definition: dualq1.hh:49
static constexpr GeometryType type()
Definition: dualq1.hh:86
Definition: dualq1localinterpolation.hh:19
A dense n x m matrix.
Definition: fmatrix.hh:69
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:123
static const QuadratureRule & rule(const GeometryType &t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
select the appropriate QuadratureRule for GeometryType t and order p
Definition: quadraturerules.hh:280
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:470
Dune namespace.
Definition: alignedallocator.hh:11
static const ReferenceElement & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:196
traits helper struct
Definition: localfiniteelementtraits.hh:11
LB LocalBasisType
Definition: localfiniteelementtraits.hh:14
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:18
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:22
A unique label for each type of element that can occur in a grid.