dune-localfunctions 2.8.0
pq22d.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_PQ22DLOCALFINITEELEMENT_HH
4#define DUNE_PQ22DLOCALFINITEELEMENT_HH
5
6#include <dune/common/fmatrix.hh>
7
9
12
13namespace Dune
14{
15 template<class D, class R>
17 {
20 public:
21 using Traits = typename LFEVariant::Traits;
22
23 PQ22DLocalFiniteElement ( const GeometryType &gt )
24 {
25 if ( gt.isTriangle() )
27 else if ( gt.isQuadrilateral() )
29 }
30
31 PQ22DLocalFiniteElement ( const GeometryType &gt, const std::vector<unsigned int> vertexmap )
32 {
33 if ( gt.isTriangle() )
34 lfeVariant_ = LagrangeSimplexLocalFiniteElement<D,R,2,2>(vertexmap);
35 else if ( gt.isQuadrilateral() )
37 }
38
39 const typename Traits::LocalBasisType& localBasis () const
40 {
41 return lfeVariant_.localBasis();
42 }
43
44 const typename Traits::LocalCoefficientsType& localCoefficients () const
45 {
46 return lfeVariant_.localCoefficients();
47 }
48
49 const typename Traits::LocalInterpolationType& localInterpolation () const
50 {
51 return lfeVariant_.localInterpolation();
52 }
53
55 unsigned int size () const
56 {
57 return lfeVariant_.size();
58 }
59
60 GeometryType type () const
61 {
62 return lfeVariant_.type();
63 }
64
65 private:
66
67 LFEVariant lfeVariant_;
68 };
69
70}
71
72#endif
Definition: bdfmcube.hh:16
typename Dune::LocalFiniteElementTraits< LocalBasis, LocalCoefficients, LocalInterpolation > Traits
Export LocalFiniteElementTraits.
Definition: localfiniteelementvariant.hh:267
unsigned int size() const
Number of shape functions.
Definition: localfiniteelementvariant.hh:372
constexpr GeometryType type() const
Number of shape functions.
Definition: localfiniteelementvariant.hh:380
const Traits::LocalBasisType & localBasis() const
Provide access to LocalBasis implementation of this LocalFiniteElement.
Definition: localfiniteelementvariant.hh:348
const Traits::LocalCoefficientsType & localCoefficients() const
Provide access to LocalCoefficients implementation of this LocalFiniteElement.
Definition: localfiniteelementvariant.hh:356
const Traits::LocalInterpolationType & localInterpolation() const
Provide access to LocalInterpolation implementation of this LocalFiniteElement.
Definition: localfiniteelementvariant.hh:364
Lagrange finite element for cubes with arbitrary compile-time dimension and polynomial order.
Definition: lagrangecube.hh:709
Lagrange finite element for simplices with arbitrary compile-time dimension and polynomial order.
Definition: lagrangesimplex.hh:836
Definition: pq22d.hh:17
typename LFEVariant::Traits Traits
Definition: pq22d.hh:21
const Traits::LocalCoefficientsType & localCoefficients() const
Definition: pq22d.hh:44
PQ22DLocalFiniteElement(const GeometryType &gt, const std::vector< unsigned int > vertexmap)
Definition: pq22d.hh:31
PQ22DLocalFiniteElement(const GeometryType &gt)
Definition: pq22d.hh:23
unsigned int size() const
Number of shape functions in this finite element.
Definition: pq22d.hh:55
const Traits::LocalInterpolationType & localInterpolation() const
Definition: pq22d.hh:49
GeometryType type() const
Definition: pq22d.hh:60
const Traits::LocalBasisType & localBasis() const
Definition: pq22d.hh:39