Dune Core Modules (2.9.0)

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// SPDX-FileCopyrightInfo: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5
6#ifndef DUNE_LOCALFUNCTIONS_MONOMIAL_HH
7#define DUNE_LOCALFUNCTIONS_MONOMIAL_HH
8
9#include <cassert>
10#include <cstddef>
11#include <cstdlib>
12#include <memory>
13#include <vector>
14
15#include <dune/geometry/type.hh>
16
17#include "common/localfiniteelementtraits.hh"
18#include "common/localtoglobaladaptors.hh"
19#include "monomial/monomiallocalbasis.hh"
20#include "monomial/monomiallocalcoefficients.hh"
21#include "monomial/monomiallocalinterpolation.hh"
22
23namespace Dune
24{
25
26
39 template<class D, class R, int d, int p>
41 {
42 constexpr static int static_size = MonomialLocalBasis<D,R,d,p>::size();
43
44 public:
50 MonomialLocalInterpolation<MonomialLocalBasis<D,R,d,p>,static_size>
52
55 : basis(), interpolation(gt_, basis), gt(gt_)
56 {}
57
60 const typename Traits::LocalBasisType& localBasis () const
61 {
62 return basis;
63 }
64
68 {
69 return coefficients;
70 }
71
75 {
76 return interpolation;
77 }
78
80 unsigned int size () const
81 {
82 return basis.size();
83 }
84
88 {
89 return gt;
90 }
91
92 private:
95 MonomialLocalInterpolation<MonomialLocalBasis<D,R,d,p>,static_size> interpolation;
96 GeometryType gt;
97 };
98
100
112 template<class Geometry, class RF, std::size_t p>
114 typedef typename Geometry::ctype DF;
115 static const std::size_t dim = Geometry::mydimension;
116
118
119 std::vector<std::shared_ptr<const LocalFE> > localFEs;
120
121 void init(const GeometryType &gt) {
122 std::size_t index = gt.id() >> 1;
123 if(localFEs.size() <= index)
124 localFEs.resize(index+1);
125 localFEs[index].reset(new LocalFE(gt));
126 }
127
128 public:
131
133
137 template<class ForwardIterator>
138 MonomialFiniteElementFactory(const ForwardIterator &begin,
139 const ForwardIterator &end)
140 {
141 for(ForwardIterator it = begin; it != end; ++it)
142 init(*it);
143 }
144
146
150 { init(gt); }
151
153
157 static_assert(dim <= 3, "MonomFiniteElementFactory knows the "
158 "available geometry types only up to dimension 3");
159
161 switch(dim) {
162 case 0 :
164 break;
165 case 1 :
167 break;
168 case 2 :
171 break;
172 case 3 :
177 break;
178 default :
179 // this should never happen -- it should be caught by the static
180 // assert above.
181 std::abort();
182 };
183 }
184
186
196 const FiniteElement make(const Geometry& geometry) {
197 std::size_t index = geometry.type().id() >> 1;
198 assert(localFEs.size() > index && localFEs[index]);
199 return FiniteElement(*localFEs[index], geometry);
200 }
201 };
202}
203
204#endif // DUNE_LOCALFUNCTIONS_MONOMIAL_HH
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:125
constexpr unsigned int id() const
Return the topology id of the type.
Definition: type.hh:376
Wrapper class for geometries.
Definition: geometry.hh:71
GeometryType type() const
Return the type of the reference element. The type can be used to access the Dune::ReferenceElement.
Definition: geometry.hh:194
GridImp::ctype ctype
define type used for coordinates in grid module
Definition: geometry.hh:100
static constexpr int mydimension
geometry dimension
Definition: geometry.hh:94
Factory for global-valued MonomFiniteElement objects.
Definition: monomial.hh:113
MonomialFiniteElementFactory(const ForwardIterator &begin, const ForwardIterator &end)
construct a MonomialFiniteElementFactory from a list of GeometryType's
Definition: monomial.hh:138
MonomialFiniteElementFactory(const GeometryType &gt)
construct a MonomialFiniteElementFactory from a single GeometryType
Definition: monomial.hh:149
const FiniteElement make(const Geometry &geometry)
construct a global-valued MonomFiniteElement
Definition: monomial.hh:196
MonomialFiniteElementFactory()
construct a MonomFiniteElementFactory for all applicable GeometryType's
Definition: monomial.hh:156
Constant shape function.
Definition: monomiallocalbasis.hh:201
static constexpr unsigned int size()
Number of shape functions.
Definition: monomiallocalbasis.hh:216
Layout map for monomial finite elements.
Definition: monomiallocalcoefficients.hh:24
Monomial basis for discontinuous Galerkin methods.
Definition: monomial.hh:41
unsigned int size() const
Number of shape functions in this finite element.
Definition: monomial.hh:80
GeometryType type() const
Definition: monomial.hh:87
const Traits::LocalInterpolationType & localInterpolation() const
Definition: monomial.hh:74
LocalFiniteElementTraits< MonomialLocalBasis< D, R, d, p >, MonomialLocalCoefficients< static_size >, MonomialLocalInterpolation< MonomialLocalBasis< D, R, d, p >, static_size > > Traits
Definition: monomial.hh:51
const Traits::LocalCoefficientsType & localCoefficients() const
Definition: monomial.hh:67
const Traits::LocalBasisType & localBasis() const
Definition: monomial.hh:60
MonomialLocalFiniteElement(const GeometryType &gt_)
Construct a MonomLocalFiniteElement.
Definition: monomial.hh:54
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:158
constexpr GeometryType line
GeometryType representing a line.
Definition: type.hh:512
constexpr GeometryType prism
GeometryType representing a 3D prism.
Definition: type.hh:542
constexpr GeometryType triangle
GeometryType representing a triangle.
Definition: type.hh:518
constexpr GeometryType quadrilateral
GeometryType representing a quadrilateral (a square).
Definition: type.hh:524
constexpr GeometryType hexahedron
GeometryType representing a hexahedron.
Definition: type.hh:548
constexpr GeometryType pyramid
GeometryType representing a 3D pyramid.
Definition: type.hh:536
constexpr GeometryType tetrahedron
GeometryType representing a tetrahedron.
Definition: type.hh:530
constexpr GeometryType vertex
GeometryType representing a vertex.
Definition: type.hh:506
Dune namespace.
Definition: alignedallocator.hh:13
Dummy struct used for documentation purposes.
Definition: documentation.hh:42
traits helper struct
Definition: localfiniteelementtraits.hh:13
LB LocalBasisType
Definition: localfiniteelementtraits.hh:16
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:20
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:24
Convert a simple scalar local finite element into a global finite element.
Definition: localtoglobaladaptors.hh:187
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)