deal.II version GIT relicensing-2013-g7f3fb24d6c 2024-10-21 19:30:00+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
24
25#include <complex>
26#include <type_traits>
27#include <utility>
28
29
31
32// Detection idiom adapted from Version 2 of the C++ Extensions for Library
33// Fundamentals, ISO/IEC TS 19568:2017
34namespace internal
35{
41 namespace SupportsOperation
42 {
43 template <class...>
44 using void_t = void;
45
52 template <class Default,
53 class AlwaysVoid,
54 template <class...>
55 class Op,
56 class... Args>
57 struct detector
58 {
59 using value_t = std::false_type;
60 using type = Default;
61 };
62
75 template <class Default, template <class...> class Op, class... Args>
76 struct detector<Default, void_t<Op<Args...>>, Op, Args...>
77 {
78 using value_t = std::true_type;
79 using type = Op<Args...>;
80 };
81
82
88 {};
89
94 struct nonesuch : private nonesuch_base
95 {
96 ~nonesuch() = delete;
97 nonesuch(const nonesuch &) = delete;
98 void
99 operator=(const nonesuch &) = delete;
100 };
101
102 template <class Default, template <class...> class Op, class... Args>
103 using detected_or = detector<Default, void, Op, Args...>;
104
105 template <template <class...> class Op, class... Args>
106 using is_detected = typename detected_or<nonesuch, Op, Args...>::value_t;
107
108 template <template <class...> class Op, class... Args>
109 using detected_t = typename detected_or<nonesuch, Op, Args...>::type;
110
111 template <class Default, template <class...> class Op, class... Args>
112 using detected_or_t = typename detected_or<Default, Op, Args...>::type;
113
114 template <class Expected, template <class...> class Op, class... Args>
115 using is_detected_exact = std::is_same<Expected, detected_t<Op, Args...>>;
116
117 template <class To, template <class...> class Op, class... Args>
119 std::is_convertible<detected_t<Op, Args...>, To>;
120 } // namespace SupportsOperation
121
122
157 template <template <class...> class Op, class... Args>
158 constexpr bool is_supported_operation =
160} // namespace internal
161
162
163
164namespace internal
165{
166 namespace TemplateConstraints
167 {
176 template <bool... Values>
177 struct all_true
178 {
179 static constexpr bool value = (Values && ...);
180 };
181
182
187 template <bool... Values>
188 struct any_true
189 {
190 static constexpr bool value = (Values || ...);
191 };
192 } // namespace TemplateConstraints
193} // namespace internal
194
201template <class Base, class... Derived>
203{
205 std::is_base_of_v<Base, Derived>...>::value;
206};
207
208
209
216template <typename Type, class... Types>
218{
220 std::is_same_v<Type, Types>...>::value;
221};
222
223
224
231template <typename Type, class... Types>
233{
235 std::is_same_v<Type, Types>...>::value;
236};
237
238
239
240/*
241 * A generalization of `std::enable_if` that only works if
242 * <i>all</i> of the given boolean template parameters are
243 * true. See [here](https://en.cppreference.com/w/cpp/types/enable_if)
244 * for what `std::enable_if` does.
245 *
246 * @note
247 * In contrast to `std::enable_if`, this template has no additional
248 * template type (which for `std::enable_if` is defaulted to `void`).
249 * As a consequence, this structure cannot be used for anything other
250 * than enabling or disabling a template declaration; in particular
251 * it cannot be used to set the return type of a function as one might
252 * do with something like
253 * @code
254 * template <typename T>
255 * typename std::enable_if<std::is_floating_point_v<T>, T>::type
256 * abs (const T t);
257 * @endcode
258 * which declares a function template `abs()` that can only be
259 * instantiated if `T` is a floating point type; this function then
260 * returns an object of type `T` as indicated by the last argument to
261 * `std::enable_if`. The reason `enable_if_all` does not allow providing
262 * this additional type is that variadic templates (here, the list of
263 * `bool` arguments) must be the last template argument.
264 */
265template <bool... Values>
267 : std::enable_if<internal::TemplateConstraints::all_true<Values...>::value>
268{};
269
270
271
272/*
273 * A generalization of `std::enable_if_t` that only works if
274 * <i>all</i> of the given boolean template parameters are
275 * true. See [here](https://en.cppreference.com/w/cpp/types/enable_if)
276 * for what `std::enable_if_t` does.
277 *
278 * @note
279 * In contrast to `std::enable_if_t`, this template has no additional
280 * template type (which for `std::enable_if` is defaulted to `void`).
281 * As a consequence, this structure cannot be used for anything other
282 * than enabling or disabling a template declaration; in particular
283 * it cannot be used to set the return type of a function as one might
284 * do with something like
285 * @code
286 * template <typename T>
287 * std::enable_if_t<std::is_floating_point_v<T>, T>
288 * abs (const T t);
289 * @endcode
290 * which declares a function template `abs()` that can only be
291 * instantiated if `T` is a floating point type; this function then
292 * returns an object of type `T` as indicated by the last argument to
293 * `std::enable_if`. The reason `enable_if_all` does not allow providing
294 * this additional type is that variadic templates (here, the list of
295 * `bool` arguments) must be the last template argument.
296 */
297template <bool... Values>
298using enable_if_all_t = typename enable_if_all<Values...>::type;
299
300
306template <typename T>
308 decltype(std::begin(std::declval<T>()), std::end(std::declval<T>()));
309
310template <typename T>
311constexpr bool has_begin_and_end =
312 internal::is_supported_operation<begin_and_end_t, T>;
313
314
322template <typename T>
324
325
331template <typename ArgType, typename ValueType>
333{
334 ValueType value;
335 ValueType
336 operator()(const ArgType &)
337 {
338 return value;
339 }
340};
341
342
343
361{
366 template <typename T>
367 static bool
368 equal(const T *p1, const T *p2)
369 {
370 return (p1 == p2);
371 }
372
373
380 template <typename T, typename U>
381 static bool
382 equal(const T *, const U *)
383 {
384 return false;
385 }
386};
387
388
389
390namespace internal
391{
402 template <typename T, typename U>
404 {
405 using type = decltype(std::declval<T>() * std::declval<U>());
406 };
407
408} // namespace internal
409
410
411
458template <typename T, typename U>
460{
461 using type =
462 typename internal::ProductTypeImpl<std::decay_t<T>, std::decay_t<U>>::type;
463};
464
465namespace internal
466{
467 // Annoyingly, there is no std::complex<T>::operator*(U) for scalars U
468 // other than T (not even in C++11, or C++14). We provide our own overloads
469 // in base/complex_overloads.h, but in order for them to work, we have to
470 // manually specify all products we want to allow:
471
472 template <typename T>
473 struct ProductTypeImpl<std::complex<T>, std::complex<T>>
474 {
475 using type = std::complex<T>;
476 };
477
478 template <typename T, typename U>
479 struct ProductTypeImpl<std::complex<T>, std::complex<U>>
480 {
481 using type = std::complex<typename ProductType<T, U>::type>;
482 };
483
484 template <typename U>
485 struct ProductTypeImpl<double, std::complex<U>>
486 {
487 using type = std::complex<typename ProductType<double, U>::type>;
488 };
489
490 template <typename T>
491 struct ProductTypeImpl<std::complex<T>, double>
492 {
493 using type = std::complex<typename ProductType<T, double>::type>;
494 };
495
496 template <typename U>
497 struct ProductTypeImpl<float, std::complex<U>>
498 {
499 using type = std::complex<typename ProductType<float, U>::type>;
500 };
501
502 template <typename T>
503 struct ProductTypeImpl<std::complex<T>, float>
504 {
505 using type = std::complex<typename ProductType<T, float>::type>;
506 };
507
508} // namespace internal
509
510
511
567template <typename T>
569
570
571template <>
572struct EnableIfScalar<double>
573{
574 using type = double;
575};
576
577template <>
578struct EnableIfScalar<float>
579{
580 using type = float;
581};
582
583template <>
584struct EnableIfScalar<long double>
585{
586 using type = long double;
587};
588
589template <>
591{
592 using type = int;
593};
594
595template <>
596struct EnableIfScalar<unsigned int>
597{
598 using type = unsigned int;
599};
600
601template <typename T>
602struct EnableIfScalar<std::complex<T>>
603{
604 using type = std::complex<T>;
605};
606
607
608// Forward declarations of vector types
609template <typename Number>
610class Vector;
611
612template <typename Number>
613class BlockVector;
614
616{
617 template <typename Number>
618 class Vector;
619
620 template <typename Number>
622
623 namespace distributed
624 {
625 template <typename Number, typename MemorySpace>
626 class Vector;
627
628 template <typename Number>
629 class BlockVector;
630 } // namespace distributed
631} // namespace LinearAlgebra
632
633#ifdef DEAL_II_WITH_PETSC
635{
636 class VectorBase;
637 class Vector;
638 class BlockVector;
639
640 namespace MPI
641 {
642 class Vector;
643 class BlockVector;
644
645 class SparseMatrix;
646 class BlockSparseMatrix;
647 } // namespace MPI
648} // namespace PETScWrappers
649#endif
650
651#ifdef DEAL_II_WITH_TRILINOS
653{
654 namespace MPI
655 {
656 class Vector;
657 class BlockVector;
658 } // namespace MPI
659} // namespace TrilinosWrappers
660
661namespace LinearAlgebra
662{
663 namespace EpetraWrappers
664 {
665 class Vector;
666 }
667
668# ifdef DEAL_II_TRILINOS_WITH_TPETRA
669 namespace TpetraWrappers
670 {
671 template <typename Number, typename MemorySpace>
672 class Vector;
673
674 template <typename Number, typename MemorySpace>
675 class SparseMatrix;
676 } // namespace TpetraWrappers
677# endif
678} // namespace LinearAlgebra
679#endif
680
681
682
687namespace concepts
688{
689#if defined(DEAL_II_HAVE_CXX20) || defined(DOXYGEN)
700 template <typename C>
701 concept is_contiguous_container = requires(C &c) {
702 {
703 std::data(c)
704 };
705 {
706 std::size(c)
707 };
708 };
709
710
720 template <int dim, int spacedim>
722 (dim >= 1 && spacedim <= 3 && dim <= spacedim);
723
724 namespace internal
725 {
732 template <typename T>
733 inline constexpr bool is_dealii_vector_type = false;
734
735 template <typename Number>
736 inline constexpr bool is_dealii_vector_type<::Vector<Number>> = true;
737
738 template <typename Number>
739 inline constexpr bool is_dealii_vector_type<::BlockVector<Number>> =
740 true;
741
742 template <typename Number>
743 inline constexpr bool
745
746 template <typename Number, typename MemorySpace>
747 inline constexpr bool is_dealii_vector_type<
749
750 template <typename Number>
751 inline constexpr bool is_dealii_vector_type<
753
754# ifdef DEAL_II_WITH_PETSC
755 template <>
757 true;
758
759 template <>
760 inline constexpr bool
762
763 template <>
764 inline constexpr bool
766
767 template <>
768 inline constexpr bool
770# endif
771
772# ifdef DEAL_II_WITH_TRILINOS
773 template <>
774 inline constexpr bool
776
777 template <>
778 inline constexpr bool
780
781 template <>
782 inline constexpr bool
784 true;
785
786# ifdef DEAL_II_TRILINOS_WITH_TPETRA
787 template <typename Number, typename MemorySpace>
788 inline constexpr bool is_dealii_vector_type<
790 true;
791# endif
792# endif
793
794
802 template <typename T>
803 inline constexpr bool is_dealii_petsc_vector_type = false;
804
805# ifdef DEAL_II_WITH_PETSC
806 template <>
807 inline constexpr bool
809
810 template <>
811 inline constexpr bool
813
814 template <>
815 inline constexpr bool
817
818 template <>
819 inline constexpr bool
821
822 template <>
823 inline constexpr bool
825 true;
826# endif
827
828
836 template <typename T>
837 inline constexpr bool is_dealii_petsc_matrix_type = false;
838
839# ifdef DEAL_II_WITH_PETSC
840 template <>
841 inline constexpr bool
843 true;
844
845 template <>
846 inline constexpr bool is_dealii_petsc_matrix_type<
848# endif
849 } // namespace internal
850
851
863 template <typename VectorType>
865 internal::is_dealii_vector_type<std::remove_cv_t<VectorType>>;
866
876 template <typename VectorType>
878 is_dealii_vector_type<VectorType> && (std::is_const_v<VectorType> == false);
879
888 template <typename VectorType>
890 internal::is_dealii_petsc_vector_type<VectorType>;
891
900 template <typename VectorType>
902 internal::is_dealii_petsc_matrix_type<VectorType>;
903#endif
904} // namespace concepts
905
906
907template <int dim, int spacedim>
909class Triangulation;
910
911template <int dim, int spacedim>
913class DoFHandler;
914
915namespace parallel
916{
917 namespace distributed
918 {
919 template <int dim, int spacedim>
921 class Triangulation;
922 }
923 namespace shared
924 {
925 template <int dim, int spacedim>
927 class Triangulation;
928 }
929 namespace fullydistributed
930 {
931 template <int dim, int spacedim>
933 class Triangulation;
934 }
935} // namespace parallel
936
937namespace concepts
938{
939#if defined(DEAL_II_HAVE_CXX20) || defined(DOXYGEN)
940 namespace internal
941 {
942 template <typename T>
943 inline constexpr bool is_triangulation_or_dof_handler = false;
944
945 template <int dim, int spacedim>
946 inline constexpr bool
948
949 template <int dim, int spacedim>
950 inline constexpr bool is_triangulation_or_dof_handler<
952
953 template <int dim, int spacedim>
954 inline constexpr bool is_triangulation_or_dof_handler<
956
957 template <int dim, int spacedim>
958 inline constexpr bool is_triangulation_or_dof_handler<
960
961 template <int dim, int spacedim>
962 inline constexpr bool
964 } // namespace internal
965
966
974 template <typename MeshType>
976 internal::is_triangulation_or_dof_handler<MeshType>;
977
984 template <typename VectorType>
985 concept is_vector_space_vector = requires(VectorType U,
986 VectorType V,
987 VectorType W,
988 typename VectorType::value_type a,
989 typename VectorType::value_type b,
990 typename VectorType::value_type s) {
991 // Check local type requirements:
992 typename VectorType::value_type;
993 typename VectorType::size_type;
994 typename VectorType::real_type;
995
996 // Check some assignment and reinit operations;
997 U.reinit(V);
998 U.reinit(V, /* omit_zeroing_entries= */ true);
999
1000 U = V;
1001 U = a; // assignment of scalar
1002
1003 U.equ(a, V);
1004
1005 // Check scaling operations
1006 U *= a;
1007 U /= a;
1008
1009 U.scale(V);
1010
1011 // Vector additions:
1012 U += V;
1013 U -= V;
1014
1015 U.add(a);
1016 U.add(a, V);
1017 U.add(a, V, b, W);
1018
1019 U.sadd(s, a, V);
1020
1021 // Norms and similar stuff:
1022 {
1023 U.mean_value()
1024 } -> std::convertible_to<typename VectorType::value_type>;
1025
1026 {
1027 U.l1_norm()
1028 } -> std::convertible_to<typename VectorType::real_type>;
1029
1030 {
1031 U.l2_norm()
1032 } -> std::convertible_to<typename VectorType::real_type>;
1033
1034 {
1035 U.linfty_norm()
1036 } -> std::convertible_to<typename VectorType::real_type>;
1037
1038 // Dot products:
1039 {
1040 U *V
1041 } -> std::convertible_to<typename VectorType::value_type>;
1042
1043 {
1044 U.add_and_dot(a, V, W)
1045 } -> std::convertible_to<typename VectorType::value_type>;
1046
1047 // Some queries:
1048 {
1049 U.size()
1050 } -> std::convertible_to<typename VectorType::size_type>;
1051
1052 {
1053 U.all_zero()
1054 } -> std::same_as<bool>;
1055
1056 {
1057 U.get_mpi_communicator()
1058 } -> std::same_as<MPI_Comm>;
1059 };
1060
1061
1070 template <typename MatrixType, typename VectorType>
1072 requires(const MatrixType &A, VectorType &dst, const VectorType &src) {
1073 A.vmult(dst, src);
1074 };
1075
1076
1085 template <typename MatrixType, typename VectorType>
1087 requires(const MatrixType &A, VectorType &dst, const VectorType &src) {
1088 A.Tvmult(dst, src);
1089 };
1090
1091#endif
1092
1093} // namespace concepts
1094
1095
1096
1098
1099#endif
#define DEAL_II_DEPRECATED
Definition config.h:205
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:175
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
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