3 #ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMASLFECACHE_HH
4 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMASLFECACHE_HH
9 #include <dune/geometry/type.hh>
10 #include <dune/geometry/typeindex.hh>
21 template<
class D,
class R, std::
size_t dim, std::
size_t order>
22 struct ImplementedRaviartThomasLocalFiniteElements
25 template<
class D,
class R>
26 struct ImplementedRaviartThomasLocalFiniteElements<D,R,2,0> :
public FixedDimLocalGeometryTypeIndex<2>
28 using FixedDimLocalGeometryTypeIndex<2>::index;
29 static auto getImplementations()
31 return std::make_tuple(
32 std::make_pair(index(GeometryTypes::triangle), []() {
return RT02DLocalFiniteElement<D,R>(); }),
33 std::make_pair(index(GeometryTypes::quadrilateral), []() {
return RT0Cube2DLocalFiniteElement<D,R>(); })
38 template<
class D,
class R>
39 struct ImplementedRaviartThomasLocalFiniteElements<D,R,2,1> :
public FixedDimLocalGeometryTypeIndex<2>
41 using FixedDimLocalGeometryTypeIndex<2>::index;
42 static auto getImplementations()
44 return std::make_tuple(
45 std::make_pair(index(GeometryTypes::triangle), []() {
return RT12DLocalFiniteElement<D,R>(); }),
46 std::make_pair(index(GeometryTypes::quadrilateral), []() {
return RT1Cube2DLocalFiniteElement<D,R>(); })
51 template<
class D,
class R>
52 struct ImplementedRaviartThomasLocalFiniteElements<D,R,2,2> :
public FixedDimLocalGeometryTypeIndex<2>
54 using FixedDimLocalGeometryTypeIndex<2>::index;
55 static auto getImplementations()
57 return std::make_tuple(
58 std::make_pair(index(GeometryTypes::quadrilateral), []() {
return RT2Cube2DLocalFiniteElement<D,R>(); })
63 template<
class D,
class R>
64 struct ImplementedRaviartThomasLocalFiniteElements<D,R,3,0> :
public FixedDimLocalGeometryTypeIndex<3>
66 using FixedDimLocalGeometryTypeIndex<3>::index;
67 static auto getImplementations()
69 return std::make_tuple(
70 std::make_pair(index(GeometryTypes::hexahedron), []() {
return RT0Cube3DLocalFiniteElement<D,R>(); })
75 template<
class D,
class R>
76 struct ImplementedRaviartThomasLocalFiniteElements<D,R,3,1> :
public FixedDimLocalGeometryTypeIndex<3>
78 using FixedDimLocalGeometryTypeIndex<3>::index;
79 static auto getImplementations()
81 return std::make_tuple(
82 std::make_pair(index(GeometryTypes::hexahedron), []() { RT1Cube3DLocalFiniteElement<D,R>(); })
100 template<
class D,
class R, std::
size_t dim, std::
size_t order>
105 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMASLFECACHE_HH