dune-localfunctions  2.9.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 // 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 #ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMASLFECACHE_HH
6 #define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMASLFECACHE_HH
7 
8 #include <tuple>
9 #include <utility>
10 
11 #include <dune/geometry/type.hh>
12 #include <dune/geometry/typeindex.hh>
13 
16 
17 namespace Dune {
18 
19 namespace Impl {
20 
21  // Provide implemented Raviart-Thomas local finite elements
22 
23  template<class D, class R, std::size_t dim, std::size_t order>
24  struct ImplementedRaviartThomasLocalFiniteElements
25  {};
26 
27  template<class D, class R>
28  struct ImplementedRaviartThomasLocalFiniteElements<D,R,2,0> : public FixedDimLocalGeometryTypeIndex<2>
29  {
30  using FixedDimLocalGeometryTypeIndex<2>::index;
31  static auto getImplementations()
32  {
33  return std::make_tuple(
34  std::make_pair(index(GeometryTypes::triangle), []() { return RT02DLocalFiniteElement<D,R>(); }),
35  std::make_pair(index(GeometryTypes::quadrilateral), []() { return RT0Cube2DLocalFiniteElement<D,R>(); })
36  );
37  }
38  };
39 
40  template<class D, class R>
41  struct ImplementedRaviartThomasLocalFiniteElements<D,R,2,1> : public FixedDimLocalGeometryTypeIndex<2>
42  {
43  using FixedDimLocalGeometryTypeIndex<2>::index;
44  static auto getImplementations()
45  {
46  return std::make_tuple(
47  std::make_pair(index(GeometryTypes::triangle), []() { return RT12DLocalFiniteElement<D,R>(); }),
48  std::make_pair(index(GeometryTypes::quadrilateral), []() { return RT1Cube2DLocalFiniteElement<D,R>(); })
49  );
50  }
51  };
52 
53  template<class D, class R>
54  struct ImplementedRaviartThomasLocalFiniteElements<D,R,2,2> : public FixedDimLocalGeometryTypeIndex<2>
55  {
56  using FixedDimLocalGeometryTypeIndex<2>::index;
57  static auto getImplementations()
58  {
59  return std::make_tuple(
60  std::make_pair(index(GeometryTypes::quadrilateral), []() { return RT2Cube2DLocalFiniteElement<D,R>(); })
61  );
62  }
63  };
64 
65  template<class D, class R>
66  struct ImplementedRaviartThomasLocalFiniteElements<D,R,3,0> : public FixedDimLocalGeometryTypeIndex<3>
67  {
68  using FixedDimLocalGeometryTypeIndex<3>::index;
69  static auto getImplementations()
70  {
71  return std::make_tuple(
72  std::make_pair(index(GeometryTypes::tetrahedron), []() { return RT03DLocalFiniteElement<D,R>(); }),
73  std::make_pair(index(GeometryTypes::hexahedron), []() { return RT0Cube3DLocalFiniteElement<D,R>(); })
74  );
75  }
76  };
77 
78  template<class D, class R>
79  struct ImplementedRaviartThomasLocalFiniteElements<D,R,3,1> : public FixedDimLocalGeometryTypeIndex<3>
80  {
81  using FixedDimLocalGeometryTypeIndex<3>::index;
82  static auto getImplementations()
83  {
84  return std::make_tuple(
85  std::make_pair(index(GeometryTypes::hexahedron), []() { RT1Cube3DLocalFiniteElement<D,R>(); })
86  );
87  }
88  };
89 
90 } // namespace Impl
91 
92 
93 
103 template<class D, class R, std::size_t dim, std::size_t order>
105 
106 } // namespace Dune
107 
108 #endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS_RAVIARTTHOMASLFECACHE_HH
Definition: bdfmcube.hh:18
A cache storing a compile time selection of local finite element implementations.
Definition: localfiniteelementvariantcache.hh:68