dune-localfunctions 3.0-git
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
10#include "qk.hh"
11#include "pk2d.hh"
12
13namespace Dune
14{
15 template<class D, class R>
17 {
18 typedef Dune::FieldVector<D,2> Domain;
19 typedef Dune::FieldVector<R,1> Range;
21
23 public:
32
33 PQ22DLocalFiniteElement ( const GeometryType &gt )
34 : gt_(gt)
35 {
36 if ( gt.isTriangle() )
38 else if ( gt.isQuadrilateral() )
40 }
41
42 PQ22DLocalFiniteElement ( const GeometryType &gt, const std::vector<unsigned int> vertexmap )
43 : gt_(gt)
44 {
45 if ( gt.isTriangle() )
46 setup( Pk2DLocalFiniteElement<D,R,2>(vertexmap) );
47 else if ( gt.isQuadrilateral() )
49 }
50
52 : gt_(other.gt_)
53 {
54 fe_ = other.fe_->clone();
55 }
56
58 {
59 delete fe_;
60 }
61
62 const LocalBasis& localBasis () const
63 {
64 return fe_->localBasis();
65 }
66
68 {
69 return fe_->localCoefficients();
70 }
71
73 {
74 return fe_->localInterpolation();
75 }
76
78 unsigned int size () const
79 {
80 return fe_->localBasis().size();
81 }
82
83 const GeometryType &type () const
84 {
85 return gt_;
86 }
87
88 private:
89
90 template <class FE>
91 void setup(const FE& fe)
92 {
94 }
95
96 const GeometryType gt_;
97 const LocalFiniteElementBase *fe_;
98 };
99
100}
101
102#endif
Definition brezzidouglasmarini1cube2d.hh:14
Type traits for LocalBasisVirtualInterface.
Definition localbasis.hh:38
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
virtual base class for a local interpolation
Definition virtualinterface.hh:288
virtual base class for a local basis
Definition virtualinterface.hh:207
virtual base class for local coefficients
Definition virtualinterface.hh:365
virtual base class for local finite elements with functions
Definition virtualinterface.hh:393
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:308
Definition pk2d.hh:23
Definition pq22d.hh:17
LocalFiniteElementTraits< LocalBasisVirtualInterface< BasisTraits >, LocalCoefficientsVirtualInterface, LocalInterpolationVirtualInterface< Domain, Range > > Traits
Definition pq22d.hh:28
const LocalInterpolation & localInterpolation() const
Definition pq22d.hh:72
const LocalCoefficients & localCoefficients() const
Definition pq22d.hh:67
PQ22DLocalFiniteElement(const GeometryType &gt, const std::vector< unsigned int > vertexmap)
Definition pq22d.hh:42
Traits::LocalBasisType LocalBasis
Definition pq22d.hh:29
PQ22DLocalFiniteElement(const GeometryType &gt)
Definition pq22d.hh:33
~PQ22DLocalFiniteElement()
Definition pq22d.hh:57
Traits::LocalInterpolationType LocalInterpolation
Definition pq22d.hh:31
const GeometryType & type() const
Definition pq22d.hh:83
unsigned int size() const
Number of shape functions in this finite element.
Definition pq22d.hh:78
Traits::LocalCoefficientsType LocalCoefficients
Definition pq22d.hh:30
PQ22DLocalFiniteElement(const PQ22DLocalFiniteElement< D, R > &other)
Definition pq22d.hh:51
const LocalBasis & localBasis() const
Definition pq22d.hh:62
General Lagrange finite element for cubes with arbitrary dimension and polynomial order.
Definition qk.hh:22