Dune Core Modules (2.6.0)

pk1d.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_PK1DLOCALFINITEELEMENT_HH
4#define DUNE_PK1DLOCALFINITEELEMENT_HH
5
6#include <cstddef>
7
9
10#include <dune/localfunctions/common/localfiniteelementtraits.hh>
11#include <dune/localfunctions/common/localtoglobaladaptors.hh>
12#include "pk1d/pk1dlocalbasis.hh"
13#include "pk1d/pk1dlocalcoefficients.hh"
14#include "pk1d/pk1dlocalinterpolation.hh"
15
16namespace Dune
17{
18
21 template<class D, class R, unsigned int k>
23 {
24 public:
29 Pk1DLocalInterpolation<Pk1DLocalBasis<D,R,k> > > Traits;
30
34 {}
35
38 Pk1DLocalFiniteElement (int variant) :
39 coefficients(variant)
40 {}
41
46 Pk1DLocalFiniteElement (const unsigned int vertexmap[3]) :
47 coefficients(vertexmap)
48 {}
49
52 const typename Traits::LocalBasisType& localBasis () const
53 {
54 return basis;
55 }
56
60 {
61 return coefficients;
62 }
63
67 {
68 return interpolation;
69 }
70
72 unsigned int size () const
73 {
74 return basis.size();
75 }
76
79 static constexpr GeometryType type ()
80 {
82 }
83
84 private:
86 Pk1DLocalCoefficients<k> coefficients;
87 Pk1DLocalInterpolation<Pk1DLocalBasis<D,R,k> > interpolation;
88 };
89
91
98 template<class Geometry, class RF, std::size_t k>
100 typedef typename Geometry::ctype DF;
102 typedef Pk1DLocalInterpolation<LocalBasis> LocalInterpolation;
103
104 public:
108 struct Traits {
111 LocalInterpolation,
112 typename Basis::Traits
113 > Interpolation;
114 typedef Pk1DLocalCoefficients<k> Coefficients;
115 };
116
117 private:
118 static const GeometryType gt;
119 static const LocalBasis localBasis;
120 static const LocalInterpolation localInterpolation;
121
122 typename Traits::Basis basis_;
123 typename Traits::Interpolation interpolation_;
124 typename Traits::Coefficients coefficients_;
125
126 public:
128
141 template<class VertexOrder>
143 const VertexOrder& vertexOrder) :
144 basis_(localBasis, geometry), interpolation_(localInterpolation),
145 coefficients_(vertexOrder.begin(0, 0))
146 { }
147
148 const typename Traits::Basis& basis() const { return basis_; }
149 const typename Traits::Interpolation& interpolation() const
150 { return interpolation_; }
151 const typename Traits::Coefficients& coefficients() const
152 { return coefficients_; }
153 const GeometryType &type() const { return gt; }
154 };
155
156 template<class Geometry, class RF, std::size_t k>
157 const GeometryType
158 Pk1DFiniteElement<Geometry, RF, k>::gt(GeometryTypes::simplex(2));
159
160 template<class Geometry, class RF, std::size_t k>
161 const typename Pk1DFiniteElement<Geometry, RF, k>::LocalBasis
162 Pk1DFiniteElement<Geometry, RF, k>::localBasis = LocalBasis();
163
164 template<class Geometry, class RF, std::size_t k>
165 const typename Pk1DFiniteElement<Geometry, RF, k>::LocalInterpolation
166 Pk1DFiniteElement<Geometry, RF, k>::localInterpolation =
167 LocalInterpolation();
168
170
180 template<class Geometry, class RF, std::size_t k>
183
185
199 template<class VertexOrder>
200 const FiniteElement make(const Geometry& geometry,
201 const VertexOrder& vertexOrder)
202 { return FiniteElement(geometry, vertexOrder); }
203 };
204}
205
206#endif
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:277
Wrapper class for geometries.
Definition: geometry.hh:67
GridImp::ctype ctype
define type used for coordinates in grid module
Definition: geometry.hh:95
Convert a local interpolation into a global interpolation.
Definition: localtoglobaladaptors.hh:147
Langrange finite element of arbitrary order on triangles.
Definition: pk1d.hh:99
Pk1DFiniteElement(const Geometry &geometry, const VertexOrder &vertexOrder)
construct a Pk1DFiniteElement
Definition: pk1d.hh:142
Lagrange shape functions of arbitrary order on the 1D reference triangle.
Definition: pk1dlocalbasis.hh:26
Layout map for Pk elements.
Definition: pk1dlocalcoefficients.hh:23
Definition: pk1d.hh:23
Pk1DLocalFiniteElement(const unsigned int vertexmap[3])
Definition: pk1d.hh:46
Pk1DLocalFiniteElement()
Definition: pk1d.hh:33
LocalFiniteElementTraits< Pk1DLocalBasis< D, R, k >, Pk1DLocalCoefficients< k >, Pk1DLocalInterpolation< Pk1DLocalBasis< D, R, k > > > Traits
Definition: pk1d.hh:29
unsigned int size() const
Number of shape functions in this finite element.
Definition: pk1d.hh:72
static constexpr GeometryType type()
Definition: pk1d.hh:79
const Traits::LocalCoefficientsType & localCoefficients() const
Definition: pk1d.hh:59
Pk1DLocalFiniteElement(int variant)
Definition: pk1d.hh:38
const Traits::LocalBasisType & localBasis() const
Definition: pk1d.hh:52
const Traits::LocalInterpolationType & localInterpolation() const
Definition: pk1d.hh:66
Convert a simple scalar local basis into a global basis.
Definition: localtoglobaladaptors.hh:63
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:178
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:147
constexpr GeometryType line
GeometryType representing a line.
Definition: type.hh:733
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:696
Dune namespace.
Definition: alignedallocator.hh:10
Dummy struct used for documentation purposes.
Definition: documentation.hh:40
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
Factory for Pk1DFiniteElement objects.
Definition: pk1d.hh:181
const FiniteElement make(const Geometry &geometry, const VertexOrder &vertexOrder)
construct Pk1DFiniteElementFactory
Definition: pk1d.hh:200
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 24, 23:30, 2024)