3 #ifndef DUNE_PQ22DLOCALFINITEELEMENT_HH
4 #define DUNE_PQ22DLOCALFINITEELEMENT_HH
6 #include <dune/common/fmatrix.hh>
15 template<
class D,
class R>
18 typedef Dune::FieldVector<D,2> Domain;
19 typedef Dune::FieldVector<R,1> Range;
36 if ( gt.isTriangle() )
38 else if ( gt.isQuadrilateral() )
45 if ( gt.isTriangle() )
47 else if ( gt.isQuadrilateral() )
54 fe_ = other.fe_->
clone();
69 return fe_->localCoefficients();
74 return fe_->localInterpolation();
83 const GeometryType &
type ()
const
91 void setup(
const FE& fe)
96 const GeometryType gt_;
97 const LocalFiniteElementBase *fe_;
virtual base class for local finite elements with functions
Definition: virtualinterface.hh:381
virtual base class for a local basis
Definition: virtualinterface.hh:24
virtual LocalFiniteElementVirtualInterface< T > * clone() const =0
virtual const Traits::LocalBasisType & localBasis() const =0
class for wrapping a finite element using the virtual interface
Definition: virtualwrappers.hh:19
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:14
~PQ22DLocalFiniteElement()
Definition: pq22d.hh:57
const LocalCoefficients & localCoefficients() const
Definition: pq22d.hh:67
Traits::LocalBasisType LocalBasis
Definition: pq22d.hh:29
PQ22DLocalFiniteElement(const GeometryType >)
Definition: pq22d.hh:33
General Lagrange finite element for cubes with arbitrary dimension and polynomial order...
Definition: qk.hh:22
virtual base class for local coefficients
Definition: virtualinterface.hh:354
const GeometryType & type() const
Definition: pq22d.hh:83
Traits::LocalInterpolationType LocalInterpolation
Definition: pq22d.hh:31
unsigned int size() const
Number of shape functions in this finite element.
Definition: pq22d.hh:78
traits helper struct
Definition: localfiniteelementtraits.hh:10
Traits::LocalCoefficientsType LocalCoefficients
Definition: pq22d.hh:30
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:22
PQ22DLocalFiniteElement(const GeometryType >, const std::vector< unsigned int > vertexmap)
Definition: pq22d.hh:42
PQ22DLocalFiniteElement(const PQ22DLocalFiniteElement< D, R > &other)
Definition: pq22d.hh:51
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:18
LocalFiniteElementTraits< LocalBasisVirtualInterface< BasisTraits >, LocalCoefficientsVirtualInterface, LocalInterpolationVirtualInterface< Domain, Range > > Traits
Definition: pq22d.hh:28
const LocalInterpolation & localInterpolation() const
Definition: pq22d.hh:72
const LocalBasis & localBasis() const
Definition: pq22d.hh:62
LB LocalBasisType
Definition: localfiniteelementtraits.hh:14
virtual base class for a local interpolation
Definition: virtualinterface.hh:21