DUNE PDELab (2.8)

monomial.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4#ifndef DUNE_LOCALFUNCTIONS_MONOMIAL_HH
5#define DUNE_LOCALFUNCTIONS_MONOMIAL_HH
6
7#include <cassert>
8#include <cstddef>
9#include <cstdlib>
10#include <memory>
11#include <vector>
12
13#include <dune/geometry/type.hh>
14
15#include "common/localfiniteelementtraits.hh"
16#include "common/localtoglobaladaptors.hh"
17#include "monomial/monomiallocalbasis.hh"
18#include "monomial/monomiallocalcoefficients.hh"
19#include "monomial/monomiallocalinterpolation.hh"
20
21namespace Dune
22{
23
24
37 template<class D, class R, int d, int p>
39 {
40 enum { static_size = MonomialLocalBasis<D,R,d,p>::size() };
41
42 public:
48 MonomialLocalInterpolation<MonomialLocalBasis<D,R,d,p>,static_size>
50
53 : basis(), interpolation(gt_, basis), gt(gt_)
54 {}
55
58 const typename Traits::LocalBasisType& localBasis () const
59 {
60 return basis;
61 }
62
66 {
67 return coefficients;
68 }
69
73 {
74 return interpolation;
75 }
76
78 unsigned int size () const
79 {
80 return basis.size();
81 }
82
86 {
87 return gt;
88 }
89
90 private:
93 MonomialLocalInterpolation<MonomialLocalBasis<D,R,d,p>,static_size> interpolation;
94 GeometryType gt;
95 };
96
98
110 template<class Geometry, class RF, std::size_t p>
112 typedef typename Geometry::ctype DF;
113 static const std::size_t dim = Geometry::mydimension;
114
116
117 std::vector<std::shared_ptr<const LocalFE> > localFEs;
118
119 void init(const GeometryType &gt) {
120 std::size_t index = gt.id() >> 1;
121 if(localFEs.size() <= index)
122 localFEs.resize(index+1);
123 localFEs[index].reset(new LocalFE(gt));
124 }
125
126 public:
129
131
135 template<class ForwardIterator>
136 MonomialFiniteElementFactory(const ForwardIterator &begin,
137 const ForwardIterator &end)
138 {
139 for(ForwardIterator it = begin; it != end; ++it)
140 init(*it);
141 }
142
144
148 { init(gt); }
149
151
155 static_assert(dim <= 3, "MonomFiniteElementFactory knows the "
156 "available geometry types only up to dimension 3");
157
159 switch(dim) {
160 case 0 :
162 break;
163 case 1 :
165 break;
166 case 2 :
169 break;
170 case 3 :
175 break;
176 default :
177 // this should never happen -- it should be caught by the static
178 // assert above.
179 std::abort();
180 };
181 }
182
184
194 const FiniteElement make(const Geometry& geometry) {
195 std::size_t index = geometry.type().id() >> 1;
196 assert(localFEs.size() > index && localFEs[index]);
197 return FiniteElement(*localFEs[index], geometry);
198 }
199 };
200}
201
202#endif // DUNE_LOCALFUNCTIONS_MONOMIAL_HH
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:123
constexpr unsigned int id() const
Return the topology id of the type.
Definition: type.hh:374
Wrapper class for geometries.
Definition: geometry.hh:67
GeometryType type() const
Return the type of the reference element. The type can be used to access the Dune::ReferenceElement.
Definition: geometry.hh:131
GridImp::ctype ctype
define type used for coordinates in grid module
Definition: geometry.hh:95
@ mydimension
Definition: geometry.hh:90
Factory for global-valued MonomFiniteElement objects.
Definition: monomial.hh:111
MonomialFiniteElementFactory(const ForwardIterator &begin, const ForwardIterator &end)
construct a MonomialFiniteElementFactory from a list of GeometryType's
Definition: monomial.hh:136
MonomialFiniteElementFactory(const GeometryType &gt)
construct a MonomialFiniteElementFactory from a single GeometryType
Definition: monomial.hh:147
const FiniteElement make(const Geometry &geometry)
construct a global-valued MonomFiniteElement
Definition: monomial.hh:194
MonomialFiniteElementFactory()
construct a MonomFiniteElementFactory for all applicable GeometryType's
Definition: monomial.hh:154
Constant shape function.
Definition: monomiallocalbasis.hh:200
static constexpr unsigned int size()
Number of shape functions.
Definition: monomiallocalbasis.hh:215
Layout map for monomial finite elements.
Definition: monomiallocalcoefficients.hh:22
Monomial basis for discontinuous Galerkin methods.
Definition: monomial.hh:39
unsigned int size() const
Number of shape functions in this finite element.
Definition: monomial.hh:78
GeometryType type() const
Definition: monomial.hh:85
const Traits::LocalInterpolationType & localInterpolation() const
Definition: monomial.hh:72
LocalFiniteElementTraits< MonomialLocalBasis< D, R, d, p >, MonomialLocalCoefficients< static_size >, MonomialLocalInterpolation< MonomialLocalBasis< D, R, d, p >, static_size > > Traits
Definition: monomial.hh:49
const Traits::LocalCoefficientsType & localCoefficients() const
Definition: monomial.hh:65
const Traits::LocalBasisType & localBasis() const
Definition: monomial.hh:58
MonomialLocalFiniteElement(const GeometryType &gt_)
Construct a MonomLocalFiniteElement.
Definition: monomial.hh:52
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:156
constexpr GeometryType line
GeometryType representing a line.
Definition: type.hh:510
constexpr GeometryType prism
GeometryType representing a 3D prism.
Definition: type.hh:540
constexpr GeometryType triangle
GeometryType representing a triangle.
Definition: type.hh:516
constexpr GeometryType quadrilateral
GeometryType representing a quadrilateral (a square).
Definition: type.hh:522
constexpr GeometryType hexahedron
GeometryType representing a hexahedron.
Definition: type.hh:546
constexpr GeometryType pyramid
GeometryType representing a 3D pyramid.
Definition: type.hh:534
constexpr GeometryType tetrahedron
GeometryType representing a tetrahedron.
Definition: type.hh:528
constexpr GeometryType vertex
GeometryType representing a vertex.
Definition: type.hh:504
Dune namespace.
Definition: alignedallocator.hh:11
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
Convert a simple scalar local finite element into a global finite element.
Definition: localtoglobaladaptors.hh:185
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 (Dec 21, 23:30, 2024)