Dune Core Modules (2.6.0)

pq22d.hh
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
7
8#include <dune/localfunctions/common/virtualinterface.hh>
9#include <dune/localfunctions/common/virtualwrappers.hh>
10#include "qk.hh"
11#include "pk2d.hh"
12
13namespace Dune
14{
15 template<class D, class R>
16 class PQ22DLocalFiniteElement
17 {
18 typedef Dune::FieldVector<D,2> Domain;
19 typedef Dune::FieldVector<R,1> Range;
20 typedef LocalBasisTraits<D,2,Domain, R,1,Range, Dune::FieldMatrix<R,1,2> > BasisTraits;
21
22 typedef typename Dune::LocalFiniteElementVirtualInterface<BasisTraits> LocalFiniteElementBase;
23 public:
24 typedef LocalFiniteElementTraits<
25 LocalBasisVirtualInterface<BasisTraits>,
26 LocalCoefficientsVirtualInterface,
27 LocalInterpolationVirtualInterface< Domain, Range >
28 > Traits;
29 typedef typename Traits::LocalBasisType LocalBasis;
30 typedef typename Traits::LocalCoefficientsType LocalCoefficients;
31 typedef typename Traits::LocalInterpolationType LocalInterpolation;
32
33 PQ22DLocalFiniteElement ( const GeometryType &gt )
34 : gt_(gt)
35 {
36 if ( gt.isTriangle() )
37 setup( Pk2DLocalFiniteElement<D,R,2>() );
38 else if ( gt.isQuadrilateral() )
39 setup( QkLocalFiniteElement<D,R,2,2>() );
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() )
48 setup( QkLocalFiniteElement<D,R,2,2>() );
49 }
50
51 PQ22DLocalFiniteElement ( const PQ22DLocalFiniteElement<D, R>& other )
52 : gt_(other.gt_)
53 {
54 fe_ = other.fe_->clone();
55 }
56
57 ~PQ22DLocalFiniteElement ( )
58 {
59 delete fe_;
60 }
61
62 const LocalBasis& localBasis () const
63 {
64 return fe_->localBasis();
65 }
66
67 const LocalCoefficients& localCoefficients () const
68 {
69 return fe_->localCoefficients();
70 }
71
72 const LocalInterpolation& localInterpolation () const
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 {
93 fe_ = new LocalFiniteElementVirtualImp<FE>(fe);
94 }
95
96 const GeometryType gt_;
97 const LocalFiniteElementBase *fe_;
98 };
99
100}
101
102#endif
vector space out of a tensor product of fields.
Definition: fvector.hh:93
virtual base class for local finite elements with functions
Definition: virtualinterface.hh:266
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:178
Implements a matrix constructed from a given type representing a field and compile-time given number ...
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:147
Dune namespace.
Definition: alignedallocator.hh:10
LB LocalBasisType
Definition: localfiniteelementtraits.hh:14
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:18
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:22
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 24, 23:30, 2024)