dune-localfunctions  2.8.0
raviartthomaslfecache.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMASLFECACHE_HH
4 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMASLFECACHE_HH
5 
6 #include <tuple>
7 #include <utility>
8 
9 #include <dune/geometry/type.hh>
10 #include <dune/geometry/typeindex.hh>
11 
14 
15 namespace Dune {
16 
17 namespace Impl {
18 
19  // Provide implemented Raviart-Thomas local finite elements
20 
21  template<class D, class R, std::size_t dim, std::size_t order>
22  struct ImplementedRaviartThomasLocalFiniteElements
23  {};
24 
25  template<class D, class R>
26  struct ImplementedRaviartThomasLocalFiniteElements<D,R,2,0> : public FixedDimLocalGeometryTypeIndex<2>
27  {
28  using FixedDimLocalGeometryTypeIndex<2>::index;
29  static auto getImplementations()
30  {
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>(); })
34  );
35  }
36  };
37 
38  template<class D, class R>
39  struct ImplementedRaviartThomasLocalFiniteElements<D,R,2,1> : public FixedDimLocalGeometryTypeIndex<2>
40  {
41  using FixedDimLocalGeometryTypeIndex<2>::index;
42  static auto getImplementations()
43  {
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>(); })
47  );
48  }
49  };
50 
51  template<class D, class R>
52  struct ImplementedRaviartThomasLocalFiniteElements<D,R,2,2> : public FixedDimLocalGeometryTypeIndex<2>
53  {
54  using FixedDimLocalGeometryTypeIndex<2>::index;
55  static auto getImplementations()
56  {
57  return std::make_tuple(
58  std::make_pair(index(GeometryTypes::quadrilateral), []() { return RT2Cube2DLocalFiniteElement<D,R>(); })
59  );
60  }
61  };
62 
63  template<class D, class R>
64  struct ImplementedRaviartThomasLocalFiniteElements<D,R,3,0> : public FixedDimLocalGeometryTypeIndex<3>
65  {
66  using FixedDimLocalGeometryTypeIndex<3>::index;
67  static auto getImplementations()
68  {
69  return std::make_tuple(
70  std::make_pair(index(GeometryTypes::tetrahedron), []() { return RT03DLocalFiniteElement<D,R>(); }),
71  std::make_pair(index(GeometryTypes::hexahedron), []() { return RT0Cube3DLocalFiniteElement<D,R>(); })
72  );
73  }
74  };
75 
76  template<class D, class R>
77  struct ImplementedRaviartThomasLocalFiniteElements<D,R,3,1> : public FixedDimLocalGeometryTypeIndex<3>
78  {
79  using FixedDimLocalGeometryTypeIndex<3>::index;
80  static auto getImplementations()
81  {
82  return std::make_tuple(
83  std::make_pair(index(GeometryTypes::hexahedron), []() { RT1Cube3DLocalFiniteElement<D,R>(); })
84  );
85  }
86  };
87 
88 } // namespace Impl
89 
90 
91 
101 template<class D, class R, std::size_t dim, std::size_t order>
103 
104 } // namespace Dune
105 
106 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMASLFECACHE_HH
Definition: bdfmcube.hh:16
A cache storing a compile time selection of local finite element implementations.
Definition: localfiniteelementvariantcache.hh:66