Dune Core Modules (unstable)

cache.hh
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 © 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_LAGRANGE_CACHE_HH
6#define DUNE_LOCALFUNCTIONS_LAGRANGE_CACHE_HH
7
8#include <map>
9#include <optional>
10#include <type_traits>
11
12#include <dune/geometry/type.hh>
14
16#include <dune/localfunctions/lagrange/equidistantpoints.hh>
17#include <dune/localfunctions/lagrange/lagrangecube.hh>
18#include <dune/localfunctions/lagrange/lagrangelfecache.hh>
19#include <dune/localfunctions/lagrange/lagrangeprism.hh>
20#include <dune/localfunctions/lagrange/lagrangepyramid.hh>
21#include <dune/localfunctions/lagrange/lagrangesimplex.hh>
22
23
24namespace Dune {
25
39template <class Domain, class Range, int dim,
40 class SF=Range, class CF=Range>
42{
43public:
44 using PointSet = EquidistantPointSet;
46
48 explicit DynamicLagrangeLocalFiniteElementCache (unsigned int order)
49 : order_(order)
50 , data_()
51 {}
52
60 {
61 auto [it,_] = data_.try_emplace(type,type,order_);
62 return it->second;
63 }
64
65private:
66 unsigned int order_;
67 mutable std::map<GeometryType, FiniteElementType> data_;
68};
69
70
82template <GeometryType::Id id, class Domain, class Range, std::size_t dim, std::size_t order>
84{
85 struct UnknownToplogy {};
86
87 static constexpr bool isSimplex = GeometryType(id).isSimplex();
88 static constexpr bool isCube = GeometryType(id).isCube();
89 static constexpr bool isPrism = GeometryType(id).isPrism();
90 static constexpr bool isPyramid = GeometryType(id).isPyramid();
91
92public:
93 using FiniteElementType
94 = std::conditional_t<isSimplex, LagrangeSimplexLocalFiniteElement<Domain,Range,dim,order>,
95 std::conditional_t<isCube, LagrangeCubeLocalFiniteElement<Domain,Range,dim,order>,
96 std::conditional_t<isPrism, LagrangePrismLocalFiniteElement<Domain,Range,order>,
97 std::conditional_t<isPyramid, LagrangePyramidLocalFiniteElement<Domain,Range,order>, UnknownToplogy> > > >;
98
100 explicit StaticLagrangeLocalFiniteElementCache (std::integral_constant<std::size_t,order> = {})
101 {
102 lfe_.emplace();
103 }
104
106 const FiniteElementType& get ([[maybe_unused]] GeometryType type) const
107 {
108 assert(GeometryType::Id(type) == id);
109 assert(!!lfe_);
110 return *lfe_;
111 }
112
113private:
114 std::optional<FiniteElementType> lfe_{};
115};
116
117
131template <class Domain, class Range, std::size_t dim, std::size_t order>
132class StaticLagrangeLocalFiniteElementCache<GeometryType::Id(~0u), Domain, Range, dim, order>
133 : public LagrangeLocalFiniteElementCache<Domain,Range,dim,order>
134{
136public:
137 using Base::Base;
138};
139
140} // namespace Dune
141
142#endif // DUNE_LOCALFUNCTIONS_LAGRANGE_CACHE_HH
A cache that stores Lagrange finite elements for the given dimension and order.
Definition: cache.hh:42
const FiniteElementType & get(GeometryType type) const
Obtain the cached local finite-element.
Definition: cache.hh:59
DynamicLagrangeLocalFiniteElementCache(unsigned int order)
Construct an empty cache.
Definition: cache.hh:48
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
constexpr bool isPyramid() const
Return true if entity is a pyramid.
Definition: type.hh:304
constexpr bool isPrism() const
Return true if entity is a prism.
Definition: type.hh:309
constexpr bool isCube() const
Return true if entity is a cube of any dimension.
Definition: type.hh:324
constexpr bool isSimplex() const
Return true if entity is a simplex of any dimension.
Definition: type.hh:319
Lagrange local finite elements for a given set of interpolation points.
Definition: lagrange.hh:69
A cache storing a compile time selection of local finite element implementations.
Definition: localfiniteelementvariantcache.hh:68
A cache that stores all available Pk/Qk like local finite elements for the given dimension and order ...
Definition: cache.hh:84
const FiniteElementType & get(GeometryType type) const
Obtain the cached local finite-element.
Definition: cache.hh:106
StaticLagrangeLocalFiniteElementCache(std::integral_constant< std::size_t, order >={})
Construct the local-finite element for the order specified as template parameter.
Definition: cache.hh:100
Convenience header that includes all implementations of Lagrange finite elements.
Dune namespace.
Definition: alignedallocator.hh:13
A unique label for each type of element that can occur in a grid.
Helper classes to provide indices for geometrytypes for use in a vector.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Sep 3, 22:42, 2025)