6#ifndef DUNE_LOCALFUNCTIONS_MONOMIAL_HH
7#define DUNE_LOCALFUNCTIONS_MONOMIAL_HH
15#include <dune/geometry/type.hh>
39 template<
class D,
class R,
int d,
int p>
55 : basis(), interpolation(gt_, basis), gt(gt_)
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;
119 std::vector<std::shared_ptr<const LocalFE> > localFEs;
121 void init(
const GeometryType >) {
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));
137 template<
class ForwardIterator>
139 const ForwardIterator &end)
141 for(ForwardIterator it = begin; it != end; ++it)
157 static_assert(dim <= 3,
"MonomFiniteElementFactory knows the "
158 "available geometry types only up to dimension 3");
163 gt = Dune::GeometryTypes::vertex; init(gt);
166 gt = Dune::GeometryTypes::line; init(gt);
169 gt = Dune::GeometryTypes::triangle; init(gt);
170 gt = Dune::GeometryTypes::quadrilateral; init(gt);
173 gt = Dune::GeometryTypes::tetrahedron; init(gt);
174 gt = Dune::GeometryTypes::pyramid; init(gt);
175 gt = Dune::GeometryTypes::prism; init(gt);
176 gt = Dune::GeometryTypes::hexahedron; init(gt);
197 std::size_t index = geometry.
type().id() >> 1;
198 assert(localFEs.size() > index && localFEs[index]);
Definition bdfmcube.hh:18
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
GeometryType type() const
Definition localtoglobaladaptors.hh:229
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 >_)
Construct a MonomLocalFiniteElement.
Definition monomial.hh:54
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 >)
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
ScalarLocalToGlobalFiniteElementAdaptor< LocalFE, Geometry > FiniteElement
Definition monomial.hh:130
MonomialFiniteElementFactory()
construct a MonomFiniteElementFactory for all applicable GeometryType's
Definition monomial.hh:156
Definition monomiallocalbasis.hh:202
static constexpr unsigned int size()
Number of shape functions.
Definition monomiallocalbasis.hh:217
Definition monomiallocalcoefficients.hh:25
Definition monomiallocalinterpolation.hh:24