DUNE PDELab (2.8)

localfiniteelementvariant.hh
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_COMMON_LOCALFINITEELEMENTVARIANT_HH
4#define DUNE_LOCALFUNCTIONS_COMMON_LOCALFINITEELEMENTVARIANT_HH
5
6#include <cstddef>
7#include <type_traits>
8#include <variant>
9
11#include <dune/common/std/type_traits.hh>
12#include <dune/common/overloadset.hh>
13
14#include <dune/geometry/type.hh>
15
16#include <dune/localfunctions/common/localfiniteelementtraits.hh>
17#include <dune/localfunctions/common/localbasis.hh>
18#include <dune/localfunctions/common/localkey.hh>
19
20
21namespace Dune {
22
23namespace Impl {
24
25 // Helper for visiting a variant containing monostate.
26 // Since a generic lambda will in most cases not compile
27 // for monostate, we add special empty overloads for monostate.
28 // Hence visitIf will simply do nothing in the case of a
29 // monostate value.
30 template<class Visitor, class Variant>
31 void visitIf(Visitor&& visitor, Variant&& variant)
32 {
33 auto visitorWithFallback = overload([&](std::monostate& impl) {}, [&](const std::monostate& impl) {}, visitor);
34 std::visit(visitorWithFallback, variant);
35 }
36
37 template<class... Implementations>
38 class LocalBasisVariant
39 {
40
41 template<class I0, class... II>
42 struct FirstType
43 { using type = I0; };
44
45 using FirstImpTraits = typename FirstType<Implementations...>::type::Traits;
46
47 public:
48
49 // We do not simply copy Implementation::LocalBasisTraits because this
50 // may be implementation specific. To stay clean, we simply put all its
51 // data into the default LocalBasisTraits.
52 using Traits = typename Dune::LocalBasisTraits<
53 typename FirstImpTraits::DomainFieldType,
54 FirstImpTraits::dimDomain,
55 typename FirstImpTraits::DomainType,
56 typename FirstImpTraits::RangeFieldType,
57 FirstImpTraits::dimRange,
58 typename FirstImpTraits::RangeType,
59 typename FirstImpTraits::JacobianType>;
60
61 template<class Implementation>
62 LocalBasisVariant(const Implementation& impl) :
63 impl_(&impl),
64 size_(impl.size()),
65 order_(impl.order())
66 {}
67
68 LocalBasisVariant() = default;
69 LocalBasisVariant(const LocalBasisVariant& other) = default;
70 LocalBasisVariant(LocalBasisVariant&& other) = default;
71 LocalBasisVariant& operator=(const LocalBasisVariant& other) = default;
72 LocalBasisVariant& operator=(LocalBasisVariant&& other) = default;
73
77 unsigned int size() const
78 {
79 return size_;
80 }
81
85 unsigned int order() const
86 {
87 return order_;
88 }
89
93 inline void evaluateFunction(
94 const typename Traits::DomainType& x,
95 std::vector<typename Traits::RangeType>& out) const
96 {
97 Impl::visitIf([&](const auto* impl) { impl->evaluateFunction(x, out); }, impl_);
98 }
99
103 inline void evaluateJacobian(
104 const typename Traits::DomainType& x,
105 std::vector<typename Traits::JacobianType>& out) const
106 {
107 Impl::visitIf([&](const auto* impl) { impl->evaluateJacobian(x, out); }, impl_);
108 }
109
117 void partial(
118 const std::array<unsigned int,Traits::dimDomain>& order,
119 const typename Traits::DomainType& x,
120 std::vector<typename Traits::RangeType>& out) const
121 {
122 Impl::visitIf([&](const auto* impl) { impl->partial(order, x, out); }, impl_);
123 }
124
125 private:
126 std::variant<std::monostate, const Implementations*...> impl_;
127 std::size_t size_;
128 std::size_t order_;
129 };
130
131
132 template<class... Implementations>
133 class LocalCoefficientsVariant
134 {
135 public:
136
137 template<class Implementation>
138 LocalCoefficientsVariant(const Implementation& impl) :
139 impl_(&impl),
140 size_(impl.size())
141 {}
142
143 LocalCoefficientsVariant() = default;
144 LocalCoefficientsVariant(const LocalCoefficientsVariant& other) = default;
145 LocalCoefficientsVariant(LocalCoefficientsVariant&& other) = default;
146 LocalCoefficientsVariant& operator=(const LocalCoefficientsVariant& other) = default;
147 LocalCoefficientsVariant& operator=(LocalCoefficientsVariant&& other) = default;
148
152 unsigned int size() const
153 {
154 return size_;
155 }
156
157 const Dune::LocalKey& localKey (std::size_t i) const
158 {
159 // We can't use visitIf since we have to return something
160 // even for a monostate value. Since the return type is
161 // an l-value reference, we use a default constructed
162 // dummy LocalKey value.
163 static const Dune::LocalKey dummyLocalKey;
164 return std::visit(overload(
165 [&](const std::monostate& impl) -> decltype(auto) { return (dummyLocalKey);},
166 [&](const auto* impl) -> decltype(auto) { return impl->localKey(i); }), impl_);
167 }
168
169 private:
170 std::variant<std::monostate, const Implementations*...> impl_;
171 std::size_t size_;
172 };
173
174
175 template<class... Implementations>
176 class LocalInterpolationVariant
177 {
178 public:
179
180 template<class Implementation>
181 LocalInterpolationVariant(const Implementation& impl) :
182 impl_(&impl)
183 {}
184
185 LocalInterpolationVariant() = default;
186 LocalInterpolationVariant(const LocalInterpolationVariant& other) = default;
187 LocalInterpolationVariant(LocalInterpolationVariant&& other) = default;
188 LocalInterpolationVariant& operator=(const LocalInterpolationVariant& other) = default;
189 LocalInterpolationVariant& operator=(LocalInterpolationVariant&& other) = default;
190
191 template<typename F, typename C>
192 void interpolate (const F& ff, std::vector<C>& out) const
193 {
194 Impl::visitIf([&](const auto* impl) { impl->interpolate(ff, out); }, impl_);
195 }
196
197 private:
198 std::variant<std::monostate, const Implementations*...> impl_;
199 };
200
201} // namespace Impl
202
203
232 template<class... Implementations>
234 {
235
236 // In each LocalFooVariant we store a std::variant<std::monostate, const FooImpl*...>, i.e. a std::variant
237 // with the pointer to the Foo implementation unless LocalFiniteElementVariant stores a monostate. In this
238 // case each LocalFooVariant also stores a monostate (and not a monostate*).
239 using LocalBasis = Impl::LocalBasisVariant<typename Implementations::Traits::LocalBasisType...>;
240 using LocalCoefficients = Impl::LocalCoefficientsVariant<typename Implementations::Traits::LocalCoefficientsType...>;
241 using LocalInterpolation = Impl::LocalInterpolationVariant<typename Implementations::Traits::LocalInterpolationType...>;
242
243 // Update members after changing impl_
244 void updateMembers()
245 {
246 std::visit(overload(
247 [&](std::monostate&) {
248 localBasis_ = LocalBasis();
249 localCoefficients_ = LocalCoefficients();
250 localInterpolation_ = LocalInterpolation();
251 size_ = 0;
252 geometryType_ = GeometryType{};
253 }, [&](auto&& impl) {
254 localBasis_ = LocalBasis(impl.localBasis());
255 localCoefficients_ = LocalCoefficients(impl.localCoefficients());
256 localInterpolation_ = LocalInterpolation(impl.localInterpolation());
257 size_ = impl.size();
258 geometryType_ = impl.type();
259 }), impl_);
260 }
261
262 public:
263
268
273
277 LocalFiniteElementVariant(const std::monostate& monostate)
278 {}
279
286 template<class Implementation,
287 std::enable_if_t<std::disjunction<std::is_same<std::decay_t<Implementation>, Implementations>...>::value, int> = 0>
288 LocalFiniteElementVariant(Implementation&& impl) :
289 impl_(std::forward<Implementation>(impl))
290 {
291 updateMembers();
292 }
293
298 impl_(other.impl_)
299 {
300 updateMembers();
301 }
302
307 impl_(std::move(other.impl_))
308 {
309 updateMembers();
310 }
311
316 {
317 impl_ = other.impl_;
318 updateMembers();
319 return *this;
320 }
321
326 {
327 impl_ = std::move(other.impl_);
328 updateMembers();
329 return *this;
330 }
331
335 template<class Implementation,
336 std::enable_if_t<std::disjunction<std::is_same<std::decay_t<Implementation>, Implementations>...>::value, int> = 0>
337 LocalFiniteElementVariant& operator=(Implementation&& impl)
338 {
339 impl_ = std::forward<Implementation>(impl);
340 updateMembers();
341 return *this;
342 }
343
344
348 const typename Traits::LocalBasisType& localBasis() const
349 {
350 return localBasis_;
351 }
352
356 const typename Traits::LocalCoefficientsType& localCoefficients() const
357 {
358 return localCoefficients_;
359 }
360
364 const typename Traits::LocalInterpolationType& localInterpolation() const
365 {
366 return localInterpolation_;
367 }
368
372 unsigned int size() const
373 {
374 return size_;
375 }
376
380 constexpr GeometryType type() const
381 {
382 return geometryType_;
383 }
384
396 const auto& variant() const
397 {
398 return impl_;
399 }
400
406 operator bool () const
407 {
408 return not(std::holds_alternative<std::monostate>(variant()));
409 }
410
411 private:
412 std::variant<std::monostate, Implementations...> impl_;
413 std::size_t size_;
414 GeometryType geometryType_;
415 LocalBasis localBasis_;
416 LocalCoefficients localCoefficients_;
417 LocalInterpolation localInterpolation_;
418 };
419
420} // end namespace Dune
421
422#endif // DUNE_LOCALFUNCTIONS_COMMON_LOCALFINITEELEMENTVARIANT_HH
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:123
Type erasure class for wrapping LocalFiniteElement classes.
Definition: localfiniteelementvariant.hh:234
typename Dune::LocalFiniteElementTraits< LocalBasis, LocalCoefficients, LocalInterpolation > Traits
Export LocalFiniteElementTraits.
Definition: localfiniteelementvariant.hh:267
const auto & variant() const
Provide access to underlying std::variant.
Definition: localfiniteelementvariant.hh:396
unsigned int size() const
Number of shape functions.
Definition: localfiniteelementvariant.hh:372
constexpr GeometryType type() const
Number of shape functions.
Definition: localfiniteelementvariant.hh:380
LocalFiniteElementVariant(LocalFiniteElementVariant &&other)
Move constructor.
Definition: localfiniteelementvariant.hh:306
LocalFiniteElementVariant & operator=(Implementation &&impl)
Assignment from implementation.
Definition: localfiniteelementvariant.hh:337
const Traits::LocalBasisType & localBasis() const
Provide access to LocalBasis implementation of this LocalFiniteElement.
Definition: localfiniteelementvariant.hh:348
LocalFiniteElementVariant & operator=(const LocalFiniteElementVariant &other)
Copy assignment.
Definition: localfiniteelementvariant.hh:315
LocalFiniteElementVariant(const std::monostate &monostate)
Construct empty LocalFiniteElementVariant.
Definition: localfiniteelementvariant.hh:277
LocalFiniteElementVariant(const LocalFiniteElementVariant &other)
Copy constructor.
Definition: localfiniteelementvariant.hh:297
const Traits::LocalCoefficientsType & localCoefficients() const
Provide access to LocalCoefficients implementation of this LocalFiniteElement.
Definition: localfiniteelementvariant.hh:356
LocalFiniteElementVariant(Implementation &&impl)
Construct LocalFiniteElementVariant.
Definition: localfiniteelementvariant.hh:288
LocalFiniteElementVariant & operator=(LocalFiniteElementVariant &&other)
Move assignment.
Definition: localfiniteelementvariant.hh:325
const Traits::LocalInterpolationType & localInterpolation() const
Provide access to LocalInterpolation implementation of this LocalFiniteElement.
Definition: localfiniteelementvariant.hh:364
LocalFiniteElementVariant()=default
Construct empty LocalFiniteElementVariant.
Describe position of one degree of freedom.
Definition: localkey.hh:21
auto overload(F &&... f)
Create an overload set.
Definition: overloadset.hh:59
void interpolate(const F &f, const GFS &gfs, XG &xg)
interpolation from a given grid function
Definition: interpolate.hh:177
Dune namespace.
Definition: alignedallocator.hh:11
STL namespace.
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:32
traits helper struct
Definition: localfiniteelementtraits.hh:11
A unique label for each type of element that can occur in a grid.
Utilities for type computations, constraining overloads, ...
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 27, 22:29, 2024)