Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
template_constraints.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2003 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_template_constraints_h
16#define dealii_template_constraints_h
17
18
19#include <deal.II/base/config.h>
20
23
24#include <complex>
25#include <type_traits>
26#include <utility>
27
28
30
31// Detection idiom adapted from Version 2 of the C++ Extensions for Library
32// Fundamentals, ISO/IEC TS 19568:2017
33namespace internal
34{
40 namespace SupportsOperation
41 {
42 template <class...>
43 using void_t = void;
44
51 template <class Default,
52 class AlwaysVoid,
53 template <class...>
54 class Op,
55 class... Args>
56 struct detector
57 {
58 using value_t = std::false_type;
59 using type = Default;
60 };
61
74 template <class Default, template <class...> class Op, class... Args>
75 struct detector<Default, void_t<Op<Args...>>, Op, Args...>
76 {
77 using value_t = std::true_type;
78 using type = Op<Args...>;
79 };
80
81
87 {};
88
93 struct nonesuch : private nonesuch_base
94 {
95 ~nonesuch() = delete;
96 nonesuch(const nonesuch &) = delete;
97 void
98 operator=(const nonesuch &) = delete;
99 };
100
101 template <class Default, template <class...> class Op, class... Args>
102 using detected_or = detector<Default, void, Op, Args...>;
103
104 template <template <class...> class Op, class... Args>
105 using is_detected = typename detected_or<nonesuch, Op, Args...>::value_t;
106
107 template <template <class...> class Op, class... Args>
108 using detected_t = typename detected_or<nonesuch, Op, Args...>::type;
109
110 template <class Default, template <class...> class Op, class... Args>
111 using detected_or_t = typename detected_or<Default, Op, Args...>::type;
112
113 template <class Expected, template <class...> class Op, class... Args>
114 using is_detected_exact = std::is_same<Expected, detected_t<Op, Args...>>;
115
116 template <class To, template <class...> class Op, class... Args>
118 std::is_convertible<detected_t<Op, Args...>, To>;
119 } // namespace SupportsOperation
120
121
156 template <template <class...> class Op, class... Args>
157 constexpr bool is_supported_operation =
159} // namespace internal
160
161
162
163namespace internal
164{
165 namespace TemplateConstraints
166 {
175 template <bool... Values>
176 struct all_true
177 {
178 static constexpr bool value = (Values && ...);
179 };
180
181
186 template <bool... Values>
187 struct any_true
188 {
189 static constexpr bool value = (Values || ...);
190 };
191 } // namespace TemplateConstraints
192} // namespace internal
193
200template <class Base, class... Derived>
202{
204 std::is_base_of_v<Base, Derived>...>::value;
205};
206
207
208
215template <typename Type, class... Types>
217{
219 std::is_same_v<Type, Types>...>::value;
220};
221
222
223
230template <typename Type, class... Types>
232{
234 std::is_same_v<Type, Types>...>::value;
235};
236
237
238
239/*
240 * A generalization of `std::enable_if` that only works if
241 * <i>all</i> of the given boolean template parameters are
242 * true. See [here](https://en.cppreference.com/w/cpp/types/enable_if)
243 * for what `std::enable_if` does.
244 *
245 * @note
246 * In contrast to `std::enable_if`, this template has no additional
247 * template type (which for `std::enable_if` is defaulted to `void`).
248 * As a consequence, this structure cannot be used for anything other
249 * than enabling or disabling a template declaration; in particular
250 * it cannot be used to set the return type of a function as one might
251 * do with something like
252 * @code
253 * template <typename T>
254 * typename std::enable_if<std::is_floating_point_v<T>, T>::type
255 * abs (const T t);
256 * @endcode
257 * which declares a function template `abs()` that can only be
258 * instantiated if `T` is a floating point type; this function then
259 * returns an object of type `T` as indicated by the last argument to
260 * `std::enable_if`. The reason `enable_if_all` does not allow providing
261 * this additional type is that variadic templates (here, the list of
262 * `bool` arguments) must be the last template argument.
263 */
264template <bool... Values>
266 : std::enable_if<internal::TemplateConstraints::all_true<Values...>::value>
267{};
268
269
270
271/*
272 * A generalization of `std::enable_if_t` that only works if
273 * <i>all</i> of the given boolean template parameters are
274 * true. See [here](https://en.cppreference.com/w/cpp/types/enable_if)
275 * for what `std::enable_if_t` does.
276 *
277 * @note
278 * In contrast to `std::enable_if_t`, this template has no additional
279 * template type (which for `std::enable_if` is defaulted to `void`).
280 * As a consequence, this structure cannot be used for anything other
281 * than enabling or disabling a template declaration; in particular
282 * it cannot be used to set the return type of a function as one might
283 * do with something like
284 * @code
285 * template <typename T>
286 * std::enable_if_t<std::is_floating_point_v<T>, T>
287 * abs (const T t);
288 * @endcode
289 * which declares a function template `abs()` that can only be
290 * instantiated if `T` is a floating point type; this function then
291 * returns an object of type `T` as indicated by the last argument to
292 * `std::enable_if`. The reason `enable_if_all` does not allow providing
293 * this additional type is that variadic templates (here, the list of
294 * `bool` arguments) must be the last template argument.
295 */
296template <bool... Values>
297using enable_if_all_t = typename enable_if_all<Values...>::type;
298
299
305template <typename T>
307 decltype(std::begin(std::declval<T>()), std::end(std::declval<T>()));
308
309template <typename T>
310constexpr bool has_begin_and_end =
311 internal::is_supported_operation<begin_and_end_t, T>;
312
313
321template <typename T>
323
324
330template <typename ArgType, typename ValueType>
332{
333 ValueType value;
334 ValueType
335 operator()(const ArgType &)
336 {
337 return value;
338 }
339};
340
341
342
360{
365 template <typename T>
366 static bool
367 equal(const T *p1, const T *p2)
368 {
369 return (p1 == p2);
370 }
371
372
379 template <typename T, typename U>
380 static bool
381 equal(const T *, const U *)
382 {
383 return false;
384 }
385};
386
387
388
389namespace internal
390{
401 template <typename T, typename U>
403 {
404 using type = decltype(std::declval<T>() * std::declval<U>());
405 };
406
407} // namespace internal
408
409
410
457template <typename T, typename U>
459{
460 using type =
461 typename internal::ProductTypeImpl<std::decay_t<T>, std::decay_t<U>>::type;
462};
463
464namespace internal
465{
466 // Annoyingly, there is no std::complex<T>::operator*(U) for scalars U
467 // other than T (not even in C++11, or C++14). We provide our own overloads
468 // in base/complex_overloads.h, but in order for them to work, we have to
469 // manually specify all products we want to allow:
470
471 template <typename T>
472 struct ProductTypeImpl<std::complex<T>, std::complex<T>>
473 {
474 using type = std::complex<T>;
475 };
476
477 template <typename T, typename U>
478 struct ProductTypeImpl<std::complex<T>, std::complex<U>>
479 {
480 using type = std::complex<typename ProductType<T, U>::type>;
481 };
482
483 template <typename U>
484 struct ProductTypeImpl<double, std::complex<U>>
485 {
486 using type = std::complex<typename ProductType<double, U>::type>;
487 };
488
489 template <typename T>
490 struct ProductTypeImpl<std::complex<T>, double>
491 {
492 using type = std::complex<typename ProductType<T, double>::type>;
493 };
494
495 template <typename U>
496 struct ProductTypeImpl<float, std::complex<U>>
497 {
498 using type = std::complex<typename ProductType<float, U>::type>;
499 };
500
501 template <typename T>
502 struct ProductTypeImpl<std::complex<T>, float>
503 {
504 using type = std::complex<typename ProductType<T, float>::type>;
505 };
506
507} // namespace internal
508
509
510
566template <typename T>
568
569
570template <>
571struct EnableIfScalar<double>
572{
573 using type = double;
574};
575
576template <>
577struct EnableIfScalar<float>
578{
579 using type = float;
580};
581
582template <>
583struct EnableIfScalar<long double>
584{
585 using type = long double;
586};
587
588template <>
590{
591 using type = int;
592};
593
594template <>
595struct EnableIfScalar<unsigned int>
596{
597 using type = unsigned int;
598};
599
600template <typename T>
601struct EnableIfScalar<std::complex<T>>
602{
603 using type = std::complex<T>;
604};
605
606
607// Forward declarations of vector types
608template <typename Number>
609class Vector;
610
611template <typename Number>
612class BlockVector;
613
615{
616 template <typename Number>
617 class Vector;
618
619 template <typename Number>
621
622 namespace distributed
623 {
624 template <typename Number, typename MemorySpace>
625 class Vector;
626
627 template <typename Number>
628 class BlockVector;
629 } // namespace distributed
630} // namespace LinearAlgebra
631
632#ifdef DEAL_II_WITH_PETSC
634{
635 class VectorBase;
636 class Vector;
637 class BlockVector;
638
639 namespace MPI
640 {
641 class Vector;
642 class BlockVector;
643
644 class SparseMatrix;
645 class BlockSparseMatrix;
646 } // namespace MPI
647} // namespace PETScWrappers
648#endif
649
650#ifdef DEAL_II_WITH_TRILINOS
652{
653 namespace MPI
654 {
655 class Vector;
656 class BlockVector;
657 } // namespace MPI
658} // namespace TrilinosWrappers
659
660namespace LinearAlgebra
661{
662 namespace EpetraWrappers
663 {
664 class Vector;
665 }
666
667# ifdef DEAL_II_TRILINOS_WITH_TPETRA
668 namespace TpetraWrappers
669 {
670 template <typename Number, typename MemorySpace>
671 class Vector;
672
673 template <typename Number, typename MemorySpace>
674 class SparseMatrix;
675 } // namespace TpetraWrappers
676# endif
677} // namespace LinearAlgebra
678#endif
679
680
681
686namespace concepts
687{
688#if defined(DEAL_II_HAVE_CXX20) || defined(DOXYGEN)
699 template <typename C>
700 concept is_contiguous_container = requires(C &c) {
701 {
702 std::data(c)
703 };
704 {
705 std::size(c)
706 };
707 };
708
709
719 template <int dim, int spacedim>
721 (dim >= 1 && spacedim <= 3 && dim <= spacedim);
722
723 namespace internal
724 {
731 template <typename T>
732 inline constexpr bool is_dealii_vector_type = false;
733
734 template <typename Number>
735 inline constexpr bool is_dealii_vector_type<::Vector<Number>> = true;
736
737 template <typename Number>
738 inline constexpr bool is_dealii_vector_type<::BlockVector<Number>> =
739 true;
740
741 template <typename Number>
742 inline constexpr bool
744
745 template <typename Number, typename MemorySpace>
746 inline constexpr bool is_dealii_vector_type<
748
749 template <typename Number>
750 inline constexpr bool is_dealii_vector_type<
752
753# ifdef DEAL_II_WITH_PETSC
754 template <>
756 true;
757
758 template <>
759 inline constexpr bool
761
762 template <>
763 inline constexpr bool
765
766 template <>
767 inline constexpr bool
769# endif
770
771# ifdef DEAL_II_WITH_TRILINOS
772 template <>
773 inline constexpr bool
775
776 template <>
777 inline constexpr bool
779
780 template <>
781 inline constexpr bool
783 true;
784
785# ifdef DEAL_II_TRILINOS_WITH_TPETRA
786 template <typename Number, typename MemorySpace>
787 inline constexpr bool is_dealii_vector_type<
789 true;
790# endif
791# endif
792
793
801 template <typename T>
802 inline constexpr bool is_dealii_petsc_vector_type = false;
803
804# ifdef DEAL_II_WITH_PETSC
805 template <>
806 inline constexpr bool
808
809 template <>
810 inline constexpr bool
812
813 template <>
814 inline constexpr bool
816
817 template <>
818 inline constexpr bool
820
821 template <>
822 inline constexpr bool
824 true;
825# endif
826
827
835 template <typename T>
836 inline constexpr bool is_dealii_petsc_matrix_type = false;
837
838# ifdef DEAL_II_WITH_PETSC
839 template <>
840 inline constexpr bool
842 true;
843
844 template <>
845 inline constexpr bool is_dealii_petsc_matrix_type<
847# endif
848 } // namespace internal
849
850
862 template <typename VectorType>
864 internal::is_dealii_vector_type<std::remove_cv_t<VectorType>>;
865
875 template <typename VectorType>
877 is_dealii_vector_type<VectorType> && (std::is_const_v<VectorType> == false);
878
887 template <typename VectorType>
889 internal::is_dealii_petsc_vector_type<VectorType>;
890
899 template <typename VectorType>
901 internal::is_dealii_petsc_matrix_type<VectorType>;
902#endif
903} // namespace concepts
904
905
906template <int dim, int spacedim>
908class Triangulation;
909
910template <int dim, int spacedim>
912class DoFHandler;
913
914namespace parallel
915{
916 namespace distributed
917 {
918 template <int dim, int spacedim>
920 class Triangulation;
921 }
922 namespace shared
923 {
924 template <int dim, int spacedim>
926 class Triangulation;
927 }
928 namespace fullydistributed
929 {
930 template <int dim, int spacedim>
932 class Triangulation;
933 }
934} // namespace parallel
935
936namespace concepts
937{
938#if defined(DEAL_II_HAVE_CXX20) || defined(DOXYGEN)
939 namespace internal
940 {
941 template <typename T>
942 inline constexpr bool is_triangulation_or_dof_handler = false;
943
944 template <int dim, int spacedim>
945 inline constexpr bool
947
948 template <int dim, int spacedim>
949 inline constexpr bool is_triangulation_or_dof_handler<
951
952 template <int dim, int spacedim>
953 inline constexpr bool is_triangulation_or_dof_handler<
955
956 template <int dim, int spacedim>
957 inline constexpr bool is_triangulation_or_dof_handler<
959
960 template <int dim, int spacedim>
961 inline constexpr bool
963 } // namespace internal
964
965
973 template <typename MeshType>
975 internal::is_triangulation_or_dof_handler<MeshType>;
976
983 template <typename VectorType>
984 concept is_vector_space_vector = requires(VectorType U,
985 VectorType V,
986 VectorType W,
987 typename VectorType::value_type a,
988 typename VectorType::value_type b,
989 typename VectorType::value_type s) {
990 // Check local type requirements:
991 typename VectorType::value_type;
992 typename VectorType::size_type;
993 typename VectorType::real_type;
994
995 // Check some assignment and reinit operations;
996 U.reinit(V);
997 U.reinit(V, /* omit_zeroing_entries= */ true);
998
999 U = V;
1000 U = a; // assignment of scalar
1001
1002 U.equ(a, V);
1003
1004 // Check scaling operations
1005 U *= a;
1006 U /= a;
1007
1008 U.scale(V);
1009
1010 // Vector additions:
1011 U += V;
1012 U -= V;
1013
1014 U.add(a);
1015 U.add(a, V);
1016 U.add(a, V, b, W);
1017
1018 U.sadd(s, a, V);
1019
1020 // Norms and similar stuff:
1021 {
1022 U.mean_value()
1023 } -> std::convertible_to<typename VectorType::value_type>;
1024
1025 {
1026 U.l1_norm()
1027 } -> std::convertible_to<typename VectorType::real_type>;
1028
1029 {
1030 U.l2_norm()
1031 } -> std::convertible_to<typename VectorType::real_type>;
1032
1033 {
1034 U.linfty_norm()
1035 } -> std::convertible_to<typename VectorType::real_type>;
1036
1037 // Dot products:
1038 {
1039 U *V
1040 } -> std::convertible_to<typename VectorType::value_type>;
1041
1042 {
1043 U.add_and_dot(a, V, W)
1044 } -> std::convertible_to<typename VectorType::value_type>;
1045
1046 // Some queries:
1047 {
1048 U.size()
1049 } -> std::convertible_to<typename VectorType::size_type>;
1050
1051 {
1052 U.all_zero()
1053 } -> std::same_as<bool>;
1054 };
1055
1056#endif
1057
1058} // namespace concepts
1059
1060
1061
1063
1064#endif
#define DEAL_II_DEPRECATED
Definition config.h:207
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:177
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
typename detected_or< Default, Op, Args... >::type detected_or_t
std::is_convertible< detected_t< Op, Args... >, To > is_detected_convertible
typename detected_or< nonesuch, Op, Args... >::value_t is_detected
typename detected_or< nonesuch, Op, Args... >::type detected_t
std::is_same< Expected, detected_t< Op, Args... > > is_detected_exact
constexpr bool is_supported_operation
STL namespace.
static bool equal(const T *p1, const T *p2)
static bool equal(const T *, const U *)
typename internal::ProductTypeImpl< std::decay_t< T >, std::decay_t< U > >::type type
static constexpr bool value
ValueType operator()(const ArgType &)
std::complex< typename ProductType< double, U >::type > type
std::complex< typename ProductType< float, U >::type > type
std::complex< typename ProductType< T, double >::type > type
std::complex< typename ProductType< T, float >::type > type
std::complex< typename ProductType< T, U >::type > type
decltype(std::declval< T >() *std::declval< U >()) type
void operator=(const nonesuch &)=delete
nonesuch(const nonesuch &)=delete
static constexpr bool value
static constexpr bool value
typename enable_if_all< Values... >::type enable_if_all_t
decltype(std::begin(std::declval< T >()), std::end(std::declval< T >())) begin_and_end_t
constexpr bool has_begin_and_end