5#ifndef DUNE_DUAL_Q1_LOCALFINITEELEMENT_HH
6#define DUNE_DUAL_Q1_LOCALFINITEELEMENT_HH
14#include <dune/geometry/referenceelements.hh>
17#include <dune/localfunctions/common/localfiniteelementtraits.hh>
18#include <dune/localfunctions/lagrange/lagrangecube.hh>
19#include "dualq1/dualq1localbasis.hh"
20#include "dualq1/dualq1localcoefficients.hh"
21#include "dualq1/dualq1localinterpolation.hh"
40 template<
class D,
class R,
int dim,
bool faceDual=false>
54 setupFaceDualCoefficients();
56 setupDualCoefficients();
95 void setupFaceDualCoefficients();
98 void setupDualCoefficients();
105 template<
class D,
class R,
int dim,
bool faceDual>
106 void DualQ1LocalFiniteElement<D,R,dim,faceDual>::setupDualCoefficients()
109 const int size = 1 <<dim;
110 std::array<Dune::FieldVector<R, size>, size> coeffs;
121 std::vector<Dune::FieldVector<R,1> > integral(size);
122 for (
int i=0; i<size; i++)
125 Dune::Impl::LagrangeCubeLocalBasis<D,R,dim,1> q1Basis;
126 for(
size_t pt=0; pt<quad.size(); pt++) {
129 std::vector<Dune::FieldVector<R,1> > q1Values(size);
130 q1Basis.evaluateFunction(pos,q1Values);
132 D weight = quad[pt].weight();
134 for (
int k=0; k<size; k++) {
135 integral[k] += q1Values[k]*weight;
137 for (
int l=0; l<=k; l++)
138 massMat[k][l] += weight*(q1Values[k]*q1Values[l]);
143 for (
int i=0; i<size-1; i++)
144 for (
int j=i+1; j<size; j++)
145 massMat[i][j] = massMat[j][i];
149 for (
int i=0; i<size; i++) {
152 rhs[i] = integral[i];
155 massMat.
solve(coeffs[i] ,rhs);
159 basis.setCoefficients(coeffs);
160 interpolation.setCoefficients(coeffs);
163 template<
class D,
class R,
int dim,
bool faceDual>
164 void DualQ1LocalFiniteElement<D,R,dim,faceDual>::setupFaceDualCoefficients()
167 const int size = 1 <<dim;
168 std::array<Dune::FieldVector<R, size>, size> coeffs;
171 Dune::Impl::LagrangeCubeLocalBasis<D,R,dim,1> q1Basis;
176 for (
int i=0; i<refElement.size(1);i++) {
187 const auto& geometry = refElement.template geometry<1>(i);
190 std::vector<Dune::FieldVector<R,1> > integral(size/2);
191 for (
int k=0; k<size/2; k++)
194 for(
size_t pt=0; pt<quad.size(); pt++) {
196 const auto& pos = quad[pt].position();
197 const auto& elementPos = geometry.global(pos);
199 std::vector<Dune::FieldVector<R,1> > q1Values(size);
200 q1Basis.evaluateFunction(elementPos,q1Values);
202 D weight = quad[pt].weight();
204 for (
int k=0; k<refElement.size(i,1,dim); k++) {
205 int row = refElement.subEntity(i,1,k,dim);
206 integral[k] += q1Values[row]*weight;
208 for (
int l=0; l<refElement.size(i,1,dim); l++) {
209 int col = refElement.subEntity(i,1,l,dim);
210 massMat[k][l] += weight*(q1Values[row]*q1Values[col]);
218 for (
int l=0; l<refElement.size(i,1,dim); l++) {
220 int row = refElement.subEntity(i,1,l,dim);
222 rhs[l] = integral[l];
225 massMat.
solve(x ,rhs);
227 for (
int k=0; k<refElement.size(i,1,dim); k++) {
228 int col = refElement.subEntity(i,1,k,dim);
229 coeffs[row][col]=x[k];
234 basis.setCoefficients(coeffs);
235 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:29
Layout map for dual Q1 elements.
Definition: dualq1localcoefficients.hh:25
The local dual Q1 finite element on cubes.
Definition: dualq1.hh:42
LocalFiniteElementTraits< DualQ1LocalBasis< D, R, dim >, DualQ1LocalCoefficients< dim >, DualQ1LocalInterpolation< dim, DualQ1LocalBasis< D, R, dim > > > Traits
Definition: dualq1.hh:47
unsigned int size() const
Number of shape functions in this finite element.
Definition: dualq1.hh:81
const Traits::LocalInterpolationType & localInterpolation() const
Definition: dualq1.hh:75
const Traits::LocalBasisType & localBasis() const
Definition: dualq1.hh:61
const Traits::LocalCoefficientsType & localCoefficients() const
Definition: dualq1.hh:68
DualQ1LocalFiniteElement()
Definition: dualq1.hh:51
static constexpr GeometryType type()
Definition: dualq1.hh:88
Definition: dualq1localinterpolation.hh:21
A dense n x m matrix.
Definition: fmatrix.hh:117
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:126
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:266
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:473
Dune namespace.
Definition: alignedallocator.hh:13
static const ReferenceElement & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:198
traits helper struct
Definition: localfiniteelementtraits.hh:13
LB LocalBasisType
Definition: localfiniteelementtraits.hh:16
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:20
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:24
A unique label for each type of element that can occur in a grid.