4#ifndef DUNE_LOCALFUNCTIONS_MONOMIAL_HH
5#define DUNE_LOCALFUNCTIONS_MONOMIAL_HH
13#include <dune/geometry/type.hh>
37 template<
class D,
class R,
int d,
int p>
53 : basis(), interpolation(gt_, basis), gt(gt_)
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;
117 std::vector<std::shared_ptr<const LocalFE> > localFEs;
119 void init(
const GeometryType >) {
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));
135 template<
class ForwardIterator>
137 const ForwardIterator &end)
139 for(ForwardIterator it = begin; it != end; ++it)
155 static_assert(dim <= 3,
"MonomFiniteElementFactory knows the "
156 "available geometry types only up to dimension 3");
161 gt = Dune::GeometryTypes::vertex; init(gt);
164 gt = Dune::GeometryTypes::line; init(gt);
167 gt = Dune::GeometryTypes::triangle; init(gt);
168 gt = Dune::GeometryTypes::quadrilateral; init(gt);
171 gt = Dune::GeometryTypes::tetrahedron; init(gt);
172 gt = Dune::GeometryTypes::pyramid; init(gt);
173 gt = Dune::GeometryTypes::prism; init(gt);
174 gt = Dune::GeometryTypes::hexahedron; init(gt);
195 std::size_t index = geometry.
type().id() >> 1;
196 assert(localFEs.size() > index && localFEs[index]);
Definition: bdfmcube.hh:16
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
GeometryType type() const
Definition: localtoglobaladaptors.hh:227
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 >_)
Construct a MonomLocalFiniteElement.
Definition: monomial.hh:52
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 >)
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
ScalarLocalToGlobalFiniteElementAdaptor< LocalFE, Geometry > FiniteElement
Definition: monomial.hh:128
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
Definition: monomiallocalinterpolation.hh:20