DUNE PDELab (git)

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
37template <class Domain, class Range, int dim>
39{
40public:
42
44 explicit DynamicLagrangeLocalFiniteElementCache (unsigned int order)
45 : order_(order)
46 , data_()
47 {}
48
56 {
57 auto [it,_] = data_.try_emplace(type,type,order_);
58 return it->second;
59 }
60
61private:
62 unsigned int order_;
63 mutable std::map<GeometryType, FiniteElementType> data_;
64};
65
66
78template <GeometryType::Id id, class Domain, class Range, std::size_t dim, std::size_t order>
80{
81 struct UnknownToplogy {};
82
83 static constexpr bool isSimplex = GeometryType(id).isSimplex();
84 static constexpr bool isCube = GeometryType(id).isCube();
85 static constexpr bool isPrism = GeometryType(id).isPrism();
86 static constexpr bool isPyramid = GeometryType(id).isPyramid();
87
88public:
89 using FiniteElementType
90 = std::conditional_t<isSimplex, LagrangeSimplexLocalFiniteElement<Domain,Range,dim,order>,
91 std::conditional_t<isCube, LagrangeCubeLocalFiniteElement<Domain,Range,dim,order>,
92 std::conditional_t<isPrism, LagrangePrismLocalFiniteElement<Domain,Range,order>,
93 std::conditional_t<isPyramid, LagrangePyramidLocalFiniteElement<Domain,Range,order>, UnknownToplogy> > > >;
94
96 explicit StaticLagrangeLocalFiniteElementCache (std::integral_constant<std::size_t,order> = {})
97 {
98 lfe_.emplace();
99 }
100
102 const FiniteElementType& get ([[maybe_unused]] GeometryType type) const
103 {
104 assert(GeometryType::Id(type) == id);
105 assert(!!lfe_);
106 return *lfe_;
107 }
108
109private:
110 std::optional<FiniteElementType> lfe_{};
111};
112
113
127template <class Domain, class Range, std::size_t dim, std::size_t order>
128class StaticLagrangeLocalFiniteElementCache<GeometryType::Id(~0u), Domain, Range, dim, order>
129 : public LagrangeLocalFiniteElementCache<Domain,Range,dim,order>
130{
132public:
133 using Base::Base;
134};
135
136} // namespace Dune
137
138#endif // DUNE_LOCALFUNCTIONS_LAGRANGE_CACHE_HH
A cache that stores Lagrange finite elements for the given dimension and order.
Definition: cache.hh:39
DynamicLagrangeLocalFiniteElementCache(unsigned int order)
Construct an empty cache.
Definition: cache.hh:44
const FiniteElementType & get(GeometryType type) const
Obtain the cached local finite-element.
Definition: cache.hh:55
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:66
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:80
const FiniteElementType & get(GeometryType type) const
Obtain the cached local finite-element.
Definition: cache.hh:102
StaticLagrangeLocalFiniteElementCache(std::integral_constant< std::size_t, order >={})
Construct the local-finite element for the order specified as template parameter.
Definition: cache.hh:96
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  |  generated with Hugo v0.111.3 (Nov 24, 23:30, 2024)