]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Replace is_detected (type) by is_supported_operation (variable).
authorWolfgang Bangerth <bangerth@colostate.edu>
Tue, 1 Feb 2022 17:31:40 +0000 (10:31 -0700)
committerWolfgang Bangerth <bangerth@colostate.edu>
Tue, 1 Feb 2022 17:31:40 +0000 (10:31 -0700)
Also provide a bunch of doc strings for the involved functions. Then also replace
all of the current uses of the is_detected type by constexpr variables.

include/deal.II/base/template_constraints.h
include/deal.II/lac/la_parallel_vector.h
include/deal.II/lac/precondition.h
include/deal.II/matrix_free/fe_evaluation.h
include/deal.II/matrix_free/matrix_free.h
include/deal.II/matrix_free/type_traits.h
include/deal.II/matrix_free/vector_access_internal.h

index 915732c0842fc087858c923bc3c9933810f3adc6..1d0e3b4ba29d0f6f51d195d656727016119e6bb5 100644 (file)
 
 DEAL_II_NAMESPACE_OPEN
 
-// Detection idiom from Version 2 of the C++ Extensions for Library
+// Detection idiom adapted from Version 2 of the C++ Extensions for Library
 // Fundamentals, ISO/IEC TS 19568:2017
 namespace internal
 {
+  /**
+   * A namespace used to declare the machinery for detecting whether a specific
+   * class supports an operation. This approach simulates C++20-style
+   * concepts with language standards before C++20.
+   */
   namespace SupportsOperation
   {
     template <class...>
     using void_t = void;
 
-    // primary template handles all types not supporting the archetypal Op
+    /**
+     * The primary template class used in detecting operations. If the
+     * compiler does not choose the specialization, then the fall-back
+     * case is this general template, which then declares member variables
+     * and types according to the failed detection.
+     */
     template <class Default,
-              class /*AlwaysVoid*/,
+              class AlwaysVoid,
               template <class...>
               class Op,
-              class... /*Args*/>
+              class... Args>
     struct detector
     {
       using value_t = std::false_type;
       using type    = Default;
     };
 
-    // specialization recognizes and handles only types supporting Op
+    /**
+     * A specialization of the general template.
+     *
+     * The trick this class uses is that, just like the general template,
+     * its second argument is always `void`, but here it is written as
+     * `void_t<Op<Args...>>` and consequently the compiler will only select this
+     * specialization if `Op<Args...>` is in fact a valid type. This means that
+     * the operation we seek to understand is indeed supported.
+     *
+     * This specialization then declares member variables and types according
+     * to the successful detection.
+     */
     template <class Default, template <class...> class Op, class... Args>
     struct detector<Default, void_t<Op<Args...>>, Op, Args...>
     {
@@ -57,10 +78,17 @@ namespace internal
     };
 
 
-    // base class for nonesuch to inherit from so it is not an aggregate
+    /**
+     * A base class for the `nonesuch` type to inherit from so it is not an
+     * aggregate.
+     */
     struct nonesuch_base
     {};
 
+    /**
+     * A type that can not be used in any reasonable way and consequently
+     * can be used to indicate a failed detection in template metaprogramming.
+     */
     struct nonesuch : private nonesuch_base
     {
       ~nonesuch()                = delete;
@@ -69,29 +97,64 @@ namespace internal
       operator=(nonesuch const &) = delete;
     };
 
-  } // namespace SupportsOperation
+    template <class Default, template <class...> class Op, class... Args>
+    using detected_or = detector<Default, void, Op, Args...>;
 
-  template <class Default, template <class...> class Op, class... Args>
-  using detected_or =
-    internal::SupportsOperation::detector<Default, void, Op, Args...>;
+    template <template <class...> class Op, class... Args>
+    using is_detected = typename detected_or<nonesuch, Op, Args...>::value_t;
 
-  template <template <class...> class Op, class... Args>
-  using is_detected =
-    typename detected_or<SupportsOperation::nonesuch, Op, Args...>::value_t;
+    template <template <class...> class Op, class... Args>
+    using detected_t = typename detected_or<nonesuch, Op, Args...>::type;
 
-  template <template <class...> class Op, class... Args>
-  using detected_t =
-    typename detected_or<SupportsOperation::nonesuch, Op, Args...>::type;
+    template <class Default, template <class...> class Op, class... Args>
+    using detected_or_t = typename detected_or<Default, Op, Args...>::type;
 
-  template <class Default, template <class...> class Op, class... Args>
-  using detected_or_t = typename detected_or<Default, Op, Args...>::type;
+    template <class Expected, template <class...> class Op, class... Args>
+    using is_detected_exact = std::is_same<Expected, detected_t<Op, Args...>>;
+
+    template <class To, template <class...> class Op, class... Args>
+    using is_detected_convertible =
+      std::is_convertible<detected_t<Op, Args...>, To>;
+  } // namespace SupportsOperation
 
-  template <class Expected, template <class...> class Op, class... Args>
-  using is_detected_exact = std::is_same<Expected, detected_t<Op, Args...>>;
 
-  template <class To, template <class...> class Op, class... Args>
-  using is_detected_convertible =
-    std::is_convertible<detected_t<Op, Args...>, To>;
+  /**
+   * A `constexpr` variable that describes whether or not `Op<Args...>` is a
+   * valid expression.
+   *
+   * The way this is used is to define an `Op` operation template that
+   * describes the operation we want to perform, and `Args` is a template
+   * pack that describes the arguments to the operation. This variable
+   * then states whether the operation, with these arguments, leads to
+   * a valid C++ expression.
+   *
+   * An example is if one wanted to find out whether a type `T` has
+   * a `get_mpi_communicator()` member function. In that case, one would write
+   * the operation as
+   * @code
+   * template <typename T>
+   * using get_mpi_communicator_op
+   *   = decltype(std::declval<T>().get_mpi_communicator());
+   * @endcode
+   * and could define a variable like
+   * @code
+   * template <typename T>
+   * constexpr bool has_get_mpi_communicator =
+   * is_supported_operation<get_mpi_communicator_op, T>;
+   * @endcode
+   *
+   * The trick used here is that `get_mpi_communicator_op` is a general
+   * template, but when used with a type that does *not* have a
+   * `get_mpi_communicator()` member variable, the `decltype(...)` operation
+   * will fail because its argument does not represent a valid expression for
+   * such a type. In other words, for such types `T` that do not have such a
+   * member function, the general template `get_mpi_communicator_op` represents
+   * a valid declaration, but the instantiation `get_mpi_communicator_op<T>`
+   * is not, and the variable declared here detects and reports this.
+   */
+  template <template <class...> class Op, class... Args>
+  constexpr bool is_supported_operation =
+    SupportsOperation::is_detected<Op, Args...>::value;
 } // namespace internal
 
 
@@ -171,7 +234,8 @@ using begin_and_end_t =
   decltype(std::begin(std::declval<T>()), std::end(std::declval<T>()));
 
 template <typename T>
-using has_begin_and_end = internal::is_detected<begin_and_end_t, T>;
+using has_begin_and_end =
+  internal::SupportsOperation::is_detected<begin_and_end_t, T>;
 
 
 
index c2e1d121c9cfea54796e25295a5012098797b17f..7c37a6588e78cdd2b838014da4fadc29dbc080e9 100644 (file)
@@ -1891,7 +1891,8 @@ namespace internal
         decltype(std::declval<T>().get_mpi_communicator());
 
       template <typename T>
-      using has_get_mpi_communicator = is_detected<get_mpi_communicator_t, T>;
+      static constexpr bool has_get_mpi_communicator =
+        is_supported_operation<get_mpi_communicator_t, T>;
 
       // A helper type-trait that leverage SFINAE to figure out if type T has
       // void T::locally_owned_domain_indices()
@@ -1900,8 +1901,8 @@ namespace internal
         decltype(std::declval<T>().locally_owned_domain_indices());
 
       template <typename T>
-      using has_locally_owned_domain_indices =
-        is_detected<locally_owned_domain_indices_t, T>;
+      static constexpr bool has_locally_owned_domain_indices =
+        is_supported_operation<locally_owned_domain_indices_t, T>;
 
       // A helper type-trait that leverage SFINAE to figure out if type T has
       // void T::locally_owned_range_indices()
@@ -1910,8 +1911,8 @@ namespace internal
         decltype(std::declval<T>().locally_owned_range_indices());
 
       template <typename T>
-      using has_locally_owned_range_indices =
-        is_detected<locally_owned_range_indices_t, T>;
+      static constexpr bool has_locally_owned_range_indices =
+        is_supported_operation<locally_owned_range_indices_t, T>;
 
       // A helper type-trait that leverage SFINAE to figure out if type T has
       // void T::initialize_dof_vector(VectorType v)
@@ -1920,14 +1921,15 @@ namespace internal
         decltype(std::declval<T>().initialize_dof_vector());
 
       template <typename T>
-      using has_initialize_dof_vector = is_detected<initialize_dof_vector_t, T>;
+      static constexpr bool has_initialize_dof_vector =
+        is_supported_operation<initialize_dof_vector_t, T>;
 
       // Used for (Trilinos/PETSc)Wrappers::SparseMatrix
-      template <typename MatrixType,
-                typename std::enable_if<
-                  has_get_mpi_communicator<MatrixType>::value &&
-                    has_locally_owned_domain_indices<MatrixType>::value,
-                  MatrixType>::type * = nullptr>
+      template <
+        typename MatrixType,
+        typename std::enable_if<has_get_mpi_communicator<MatrixType> &&
+                                  has_locally_owned_domain_indices<MatrixType>,
+                                MatrixType>::type * = nullptr>
       static void
       reinit_domain_vector(MatrixType &                                mat,
                            LinearAlgebra::distributed::Vector<Number> &vec,
@@ -1938,10 +1940,9 @@ namespace internal
       }
 
       // Used for MatrixFree and DiagonalMatrix
-      template <
-        typename MatrixType,
-        typename std::enable_if<has_initialize_dof_vector<MatrixType>::value,
-                                MatrixType>::type * = nullptr>
+      template <typename MatrixType,
+                typename std::enable_if<has_initialize_dof_vector<MatrixType>,
+                                        MatrixType>::type * = nullptr>
       static void
       reinit_domain_vector(MatrixType &                                mat,
                            LinearAlgebra::distributed::Vector<Number> &vec,
@@ -1953,11 +1954,11 @@ namespace internal
       }
 
       // Used for (Trilinos/PETSc)Wrappers::SparseMatrix
-      template <typename MatrixType,
-                typename std::enable_if<
-                  has_get_mpi_communicator<MatrixType>::value &&
-                    has_locally_owned_range_indices<MatrixType>::value,
-                  MatrixType>::type * = nullptr>
+      template <
+        typename MatrixType,
+        typename std::enable_if<has_get_mpi_communicator<MatrixType> &&
+                                  has_locally_owned_range_indices<MatrixType>,
+                                MatrixType>::type * = nullptr>
       static void
       reinit_range_vector(MatrixType &                                mat,
                           LinearAlgebra::distributed::Vector<Number> &vec,
@@ -1968,10 +1969,9 @@ namespace internal
       }
 
       // Used for MatrixFree and DiagonalMatrix
-      template <
-        typename MatrixType,
-        typename std::enable_if<has_initialize_dof_vector<MatrixType>::value,
-                                MatrixType>::type * = nullptr>
+      template <typename MatrixType,
+                typename std::enable_if<has_initialize_dof_vector<MatrixType>,
+                                        MatrixType>::type * = nullptr>
       static void
       reinit_range_vector(MatrixType &                                mat,
                           LinearAlgebra::distributed::Vector<Number> &vec,
index 6e0b1ebbb55040e934a52d688ecb0ecd88af0d31..142993d20d2ef68ca7897d1b02fdd2dc8f1afe3e 100644 (file)
@@ -533,7 +533,7 @@ namespace internal
                                      std::declval<const VectorType &>()));
 
     template <typename T, typename VectorType>
-    using has_Tvmult = is_detected<Tvmult_t, T, VectorType>;
+    constexpr bool has_Tvmult = is_supported_operation<Tvmult_t, T, VectorType>;
 
     template <typename T, typename VectorType>
     using step_t = decltype(
@@ -541,7 +541,7 @@ namespace internal
                                    std::declval<const VectorType &>()));
 
     template <typename T, typename VectorType>
-    using has_step = is_detected<step_t, T, VectorType>;
+    constexpr bool has_step = is_supported_operation<step_t, T, VectorType>;
 
     template <typename T, typename VectorType>
     using step_omega_t =
@@ -550,7 +550,8 @@ namespace internal
                                             std::declval<const double>()));
 
     template <typename T, typename VectorType>
-    using has_step_omega = is_detected<step_omega_t, T, VectorType>;
+    constexpr bool has_step_omega =
+      is_supported_operation<step_omega_t, T, VectorType>;
 
     template <typename T, typename VectorType>
     using Tstep_t = decltype(
@@ -558,7 +559,7 @@ namespace internal
                                     std::declval<const VectorType &>()));
 
     template <typename T, typename VectorType>
-    using has_Tstep = is_detected<Tstep_t, T, VectorType>;
+    constexpr bool has_Tstep = is_supported_operation<Tstep_t, T, VectorType>;
 
     template <typename T, typename VectorType>
     using Tstep_omega_t =
@@ -567,7 +568,8 @@ namespace internal
                                              std::declval<const double>()));
 
     template <typename T, typename VectorType>
-    using has_Tstep_omega = is_detected<Tstep_omega_t, T, VectorType>;
+    constexpr bool has_Tstep_omega =
+      is_supported_operation<Tstep_omega_t, T, VectorType>;
 
     template <typename T, typename VectorType>
     using jacobi_step_t = decltype(
@@ -576,7 +578,8 @@ namespace internal
                                           std::declval<const double>()));
 
     template <typename T, typename VectorType>
-    using has_jacobi_step = is_detected<jacobi_step_t, T, VectorType>;
+    constexpr bool has_jacobi_step =
+      is_supported_operation<jacobi_step_t, T, VectorType>;
 
     template <typename T, typename VectorType>
     using SOR_step_t = decltype(
@@ -585,7 +588,8 @@ namespace internal
                                        std::declval<const double>()));
 
     template <typename T, typename VectorType>
-    using has_SOR_step = is_detected<SOR_step_t, T, VectorType>;
+    constexpr bool has_SOR_step =
+      is_supported_operation<SOR_step_t, T, VectorType>;
 
     template <typename T, typename VectorType>
     using SSOR_step_t = decltype(
@@ -594,7 +598,8 @@ namespace internal
                                         std::declval<const double>()));
 
     template <typename T, typename VectorType>
-    using has_SSOR_step = is_detected<SSOR_step_t, T, VectorType>;
+    constexpr bool has_SSOR_step =
+      is_supported_operation<SSOR_step_t, T, VectorType>;
 
     template <typename MatrixType>
     class PreconditionJacobiImpl
@@ -620,10 +625,9 @@ namespace internal
         this->vmult(dst, src);
       }
 
-      template <
-        typename VectorType,
-        typename std::enable_if<has_jacobi_step<MatrixType, VectorType>::value,
-                                MatrixType>::type * = nullptr>
+      template <typename VectorType,
+                typename std::enable_if<has_jacobi_step<MatrixType, VectorType>,
+                                        MatrixType>::type * = nullptr>
       void
       step(VectorType &dst, const VectorType &src) const
       {
@@ -632,7 +636,7 @@ namespace internal
 
       template <
         typename VectorType,
-        typename std::enable_if<!has_jacobi_step<MatrixType, VectorType>::value,
+        typename std::enable_if<!has_jacobi_step<MatrixType, VectorType>,
                                 MatrixType>::type * = nullptr>
       void
       step(VectorType &, const VectorType &) const
@@ -678,20 +682,18 @@ namespace internal
         this->A->precondition_TSOR(dst, src, this->relaxation);
       }
 
-      template <
-        typename VectorType,
-        typename std::enable_if<has_SOR_step<MatrixType, VectorType>::value,
-                                MatrixType>::type * = nullptr>
+      template <typename VectorType,
+                typename std::enable_if<has_SOR_step<MatrixType, VectorType>,
+                                        MatrixType>::type * = nullptr>
       void
       step(VectorType &dst, const VectorType &src) const
       {
         this->A->SOR_step(dst, src, this->relaxation);
       }
 
-      template <
-        typename VectorType,
-        typename std::enable_if<!has_SOR_step<MatrixType, VectorType>::value,
-                                MatrixType>::type * = nullptr>
+      template <typename VectorType,
+                typename std::enable_if<!has_SOR_step<MatrixType, VectorType>,
+                                        MatrixType>::type * = nullptr>
       void
       step(VectorType &, const VectorType &) const
       {
@@ -699,20 +701,18 @@ namespace internal
                ExcMessage("Matrix A does not provide a SOR_step() function!"));
       }
 
-      template <
-        typename VectorType,
-        typename std::enable_if<has_SOR_step<MatrixType, VectorType>::value,
-                                MatrixType>::type * = nullptr>
+      template <typename VectorType,
+                typename std::enable_if<has_SOR_step<MatrixType, VectorType>,
+                                        MatrixType>::type * = nullptr>
       void
       Tstep(VectorType &dst, const VectorType &src) const
       {
         this->A->TSOR_step(dst, src, this->relaxation);
       }
 
-      template <
-        typename VectorType,
-        typename std::enable_if<!has_SOR_step<MatrixType, VectorType>::value,
-                                MatrixType>::type * = nullptr>
+      template <typename VectorType,
+                typename std::enable_if<!has_SOR_step<MatrixType, VectorType>,
+                                        MatrixType>::type * = nullptr>
       void
       Tstep(VectorType &, const VectorType &) const
       {
@@ -783,20 +783,18 @@ namespace internal
                                    pos_right_of_diagonal);
       }
 
-      template <
-        typename VectorType,
-        typename std::enable_if<has_SSOR_step<MatrixType, VectorType>::value,
-                                MatrixType>::type * = nullptr>
+      template <typename VectorType,
+                typename std::enable_if<has_SSOR_step<MatrixType, VectorType>,
+                                        MatrixType>::type * = nullptr>
       void
       step(VectorType &dst, const VectorType &src) const
       {
         this->A->SSOR_step(dst, src, this->relaxation);
       }
 
-      template <
-        typename VectorType,
-        typename std::enable_if<!has_SSOR_step<MatrixType, VectorType>::value,
-                                MatrixType>::type * = nullptr>
+      template <typename VectorType,
+                typename std::enable_if<!has_SSOR_step<MatrixType, VectorType>,
+                                        MatrixType>::type * = nullptr>
       void
       step(VectorType &, const VectorType &) const
       {
@@ -863,12 +861,12 @@ namespace internal
       const std::vector<size_type> &inverse_permutation;
     };
 
-    template <typename MatrixType,
-              typename PreconditionerType,
-              typename VectorType,
-              typename std::enable_if<
-                has_step_omega<PreconditionerType, VectorType>::value,
-                PreconditionerType>::type * = nullptr>
+    template <
+      typename MatrixType,
+      typename PreconditionerType,
+      typename VectorType,
+      typename std::enable_if<has_step_omega<PreconditionerType, VectorType>,
+                              PreconditionerType>::type * = nullptr>
     void
     step(const MatrixType &,
          const PreconditionerType &preconditioner,
@@ -881,13 +879,13 @@ namespace internal
       preconditioner.step(dst, src, relaxation);
     }
 
-    template <typename MatrixType,
-              typename PreconditionerType,
-              typename VectorType,
-              typename std::enable_if<
-                !has_step_omega<PreconditionerType, VectorType>::value &&
-                  has_step<PreconditionerType, VectorType>::value,
-                PreconditionerType>::type * = nullptr>
+    template <
+      typename MatrixType,
+      typename PreconditionerType,
+      typename VectorType,
+      typename std::enable_if<!has_step_omega<PreconditionerType, VectorType> &&
+                                has_step<PreconditionerType, VectorType>,
+                              PreconditionerType>::type * = nullptr>
     void
     step(const MatrixType &,
          const PreconditionerType &preconditioner,
@@ -904,13 +902,13 @@ namespace internal
       preconditioner.step(dst, src);
     }
 
-    template <typename MatrixType,
-              typename PreconditionerType,
-              typename VectorType,
-              typename std::enable_if<
-                !has_step_omega<PreconditionerType, VectorType>::value &&
-                  !has_step<PreconditionerType, VectorType>::value,
-                PreconditionerType>::type * = nullptr>
+    template <
+      typename MatrixType,
+      typename PreconditionerType,
+      typename VectorType,
+      typename std::enable_if<!has_step_omega<PreconditionerType, VectorType> &&
+                                !has_step<PreconditionerType, VectorType>,
+                              PreconditionerType>::type * = nullptr>
     void
     step(const MatrixType &        A,
          const PreconditionerType &preconditioner,
@@ -930,12 +928,12 @@ namespace internal
       dst.add(relaxation, tmp);
     }
 
-    template <typename MatrixType,
-              typename PreconditionerType,
-              typename VectorType,
-              typename std::enable_if<
-                has_Tstep_omega<PreconditionerType, VectorType>::value,
-                PreconditionerType>::type * = nullptr>
+    template <
+      typename MatrixType,
+      typename PreconditionerType,
+      typename VectorType,
+      typename std::enable_if<has_Tstep_omega<PreconditionerType, VectorType>,
+                              PreconditionerType>::type * = nullptr>
     void
     Tstep(const MatrixType &,
           const PreconditionerType &preconditioner,
@@ -952,8 +950,8 @@ namespace internal
               typename PreconditionerType,
               typename VectorType,
               typename std::enable_if<
-                !has_Tstep_omega<PreconditionerType, VectorType>::value &&
-                  has_Tstep<PreconditionerType, VectorType>::value,
+                !has_Tstep_omega<PreconditionerType, VectorType> &&
+                  has_Tstep<PreconditionerType, VectorType>,
                 PreconditionerType>::type * = nullptr>
     void
     Tstep(const MatrixType &,
@@ -973,7 +971,7 @@ namespace internal
 
     template <typename MatrixType,
               typename VectorType,
-              typename std::enable_if<has_Tvmult<MatrixType, VectorType>::value,
+              typename std::enable_if<has_Tvmult<MatrixType, VectorType>,
                                       MatrixType>::type * = nullptr>
     void
     Tvmult(const MatrixType &A, VectorType &dst, const VectorType &src)
@@ -981,11 +979,10 @@ namespace internal
       A.Tvmult(dst, src);
     }
 
-    template <
-      typename MatrixType,
-      typename VectorType,
-      typename std::enable_if<!has_Tvmult<MatrixType, VectorType>::value,
-                              MatrixType>::type * = nullptr>
+    template <typename MatrixType,
+              typename VectorType,
+              typename std::enable_if<!has_Tvmult<MatrixType, VectorType>,
+                                      MatrixType>::type * = nullptr>
     void
     Tvmult(const MatrixType &, VectorType &, const VectorType &)
     {
@@ -997,8 +994,8 @@ namespace internal
               typename PreconditionerType,
               typename VectorType,
               typename std::enable_if<
-                !has_Tstep_omega<PreconditionerType, VectorType>::value &&
-                  !has_Tstep<PreconditionerType, VectorType>::value,
+                !has_Tstep_omega<PreconditionerType, VectorType> &&
+                  !has_Tstep<PreconditionerType, VectorType>,
                 PreconditionerType>::type * = nullptr>
     void
     Tstep(const MatrixType &        A,
index 8627949db37fe7ef466daaa148e8fa23ce714d5a..449743975f5bc06c49d09527a24636c5818bcc95 100644 (file)
@@ -4080,7 +4080,7 @@ namespace internal
   }
 
   template <typename VectorType,
-            typename std::enable_if<has_shared_vector_data<VectorType>::value,
+            typename std::enable_if<has_shared_vector_data<VectorType>,
                                     VectorType>::type * = nullptr>
   const std::vector<ArrayView<const typename VectorType::value_type>> *
   get_shared_vector_data(VectorType &       vec,
@@ -4099,7 +4099,7 @@ namespace internal
   }
 
   template <typename VectorType,
-            typename std::enable_if<!has_shared_vector_data<VectorType>::value,
+            typename std::enable_if<!has_shared_vector_data<VectorType>,
                                     VectorType>::type * = nullptr>
   const std::vector<ArrayView<const typename VectorType::value_type>> *
   get_shared_vector_data(VectorType &,
@@ -7258,7 +7258,7 @@ namespace internal
             typename VectorType,
             typename EvaluatorType,
             typename std::enable_if<
-              internal::has_begin<VectorType>::value &&
+              internal::has_begin<VectorType> &&
                 std::is_same<decltype(std::declval<VectorType>().begin()),
                              Number *>::value,
               VectorType>::type * = nullptr>
@@ -7308,7 +7308,7 @@ namespace internal
             typename VectorType,
             typename EvaluatorType,
             typename std::enable_if<
-              !internal::has_begin<VectorType>::value ||
+              !internal::has_begin<VectorType> ||
                 !std::is_same<decltype(std::declval<VectorType>().begin()),
                               Number *>::value,
               VectorType>::type * = nullptr>
index 8731d802b7f11c71362958cd9922a382a4a3259e..8f41bb7fe806dfe403e615f769b7993fbd5802f6 100644 (file)
@@ -3371,11 +3371,11 @@ namespace internal
      * Start update_ghost_value for vectors that do not support
      * the split into _start() and finish() stages
      */
-    template <typename VectorType,
-              typename std::enable_if<
-                !has_update_ghost_values_start<VectorType>::value &&
-                  !is_not_parallel_vector<VectorType>::value,
-                VectorType>::type * = nullptr>
+    template <
+      typename VectorType,
+      typename std::enable_if<!has_update_ghost_values_start<VectorType> &&
+                                !is_not_parallel_vector<VectorType>::value,
+                              VectorType>::type * = nullptr>
     void
     update_ghost_values_start(const unsigned int component_in_block_vector,
                               const VectorType & vec)
@@ -3400,11 +3400,11 @@ namespace internal
      * the split into _start() and finish() stages, but don't support
      * exchange on a subset of DoFs
      */
-    template <typename VectorType,
-              typename std::enable_if<
-                has_update_ghost_values_start<VectorType>::value &&
-                  !has_exchange_on_subset<VectorType>::value,
-                VectorType>::type * = nullptr>
+    template <
+      typename VectorType,
+      typename std::enable_if<has_update_ghost_values_start<VectorType> &&
+                                !has_exchange_on_subset<VectorType>,
+                              VectorType>::type * = nullptr>
     void
     update_ghost_values_start(const unsigned int component_in_block_vector,
                               const VectorType & vec)
@@ -3430,11 +3430,11 @@ namespace internal
      * exchange on a subset of DoFs,
      * i.e. LinearAlgebra::distributed::Vector
      */
-    template <typename VectorType,
-              typename std::enable_if<
-                has_update_ghost_values_start<VectorType>::value &&
-                  has_exchange_on_subset<VectorType>::value,
-                VectorType>::type * = nullptr>
+    template <
+      typename VectorType,
+      typename std::enable_if<has_update_ghost_values_start<VectorType> &&
+                                has_exchange_on_subset<VectorType>,
+                              VectorType>::type * = nullptr>
     void
     update_ghost_values_start(const unsigned int component_in_block_vector,
                               const VectorType & vec)
@@ -3492,7 +3492,7 @@ namespace internal
      */
     template <
       typename VectorType,
-      typename std::enable_if<!has_update_ghost_values_start<VectorType>::value,
+      typename std::enable_if<!has_update_ghost_values_start<VectorType>,
                               VectorType>::type * = nullptr>
     void
     update_ghost_values_finish(const unsigned int /*component_in_block_vector*/,
@@ -3506,11 +3506,11 @@ namespace internal
      * the split into _start() and finish() stages, but don't support
      * exchange on a subset of DoFs
      */
-    template <typename VectorType,
-              typename std::enable_if<
-                has_update_ghost_values_start<VectorType>::value &&
-                  !has_exchange_on_subset<VectorType>::value,
-                VectorType>::type * = nullptr>
+    template <
+      typename VectorType,
+      typename std::enable_if<has_update_ghost_values_start<VectorType> &&
+                                !has_exchange_on_subset<VectorType>,
+                              VectorType>::type * = nullptr>
     void
     update_ghost_values_finish(const unsigned int component_in_block_vector,
                                const VectorType & vec)
@@ -3527,11 +3527,11 @@ namespace internal
      * exchange on a subset of DoFs,
      * i.e. LinearAlgebra::distributed::Vector
      */
-    template <typename VectorType,
-              typename std::enable_if<
-                has_update_ghost_values_start<VectorType>::value &&
-                  has_exchange_on_subset<VectorType>::value,
-                VectorType>::type * = nullptr>
+    template <
+      typename VectorType,
+      typename std::enable_if<has_update_ghost_values_start<VectorType> &&
+                                has_exchange_on_subset<VectorType>,
+                              VectorType>::type * = nullptr>
     void
     update_ghost_values_finish(const unsigned int component_in_block_vector,
                                const VectorType & vec)
@@ -3595,7 +3595,7 @@ namespace internal
      */
     template <
       typename VectorType,
-      typename std::enable_if<!has_compress_start<VectorType>::value &&
+      typename std::enable_if<!has_compress_start<VectorType> &&
                                 !is_not_parallel_vector<VectorType>::value,
                               VectorType>::type * = nullptr>
     void
@@ -3614,11 +3614,10 @@ namespace internal
      * the split into _start() and finish() stages, but don't support
      * exchange on a subset of DoFs
      */
-    template <
-      typename VectorType,
-      typename std::enable_if<has_compress_start<VectorType>::value &&
-                                !has_exchange_on_subset<VectorType>::value,
-                              VectorType>::type * = nullptr>
+    template <typename VectorType,
+              typename std::enable_if<has_compress_start<VectorType> &&
+                                        !has_exchange_on_subset<VectorType>,
+                                      VectorType>::type * = nullptr>
     void
     compress_start(const unsigned int component_in_block_vector,
                    VectorType &       vec)
@@ -3636,11 +3635,10 @@ namespace internal
      * exchange on a subset of DoFs,
      * i.e. LinearAlgebra::distributed::Vector
      */
-    template <
-      typename VectorType,
-      typename std::enable_if<has_compress_start<VectorType>::value &&
-                                has_exchange_on_subset<VectorType>::value,
-                              VectorType>::type * = nullptr>
+    template <typename VectorType,
+              typename std::enable_if<has_compress_start<VectorType> &&
+                                        has_exchange_on_subset<VectorType>,
+                                      VectorType>::type * = nullptr>
     void
     compress_start(const unsigned int component_in_block_vector,
                    VectorType &       vec)
@@ -3690,7 +3688,7 @@ namespace internal
      * the split into _start() and finish() stages and serial vectors
      */
     template <typename VectorType,
-              typename std::enable_if<!has_compress_start<VectorType>::value,
+              typename std::enable_if<!has_compress_start<VectorType>,
                                       VectorType>::type * = nullptr>
     void
     compress_finish(const unsigned int /*component_in_block_vector*/,
@@ -3704,11 +3702,10 @@ namespace internal
      * the split into _start() and finish() stages, but don't support
      * exchange on a subset of DoFs
      */
-    template <
-      typename VectorType,
-      typename std::enable_if<has_compress_start<VectorType>::value &&
-                                !has_exchange_on_subset<VectorType>::value,
-                              VectorType>::type * = nullptr>
+    template <typename VectorType,
+              typename std::enable_if<has_compress_start<VectorType> &&
+                                        !has_exchange_on_subset<VectorType>,
+                                      VectorType>::type * = nullptr>
     void
     compress_finish(const unsigned int component_in_block_vector,
                     VectorType &       vec)
@@ -3725,11 +3722,10 @@ namespace internal
      * exchange on a subset of DoFs,
      * i.e. LinearAlgebra::distributed::Vector
      */
-    template <
-      typename VectorType,
-      typename std::enable_if<has_compress_start<VectorType>::value &&
-                                has_exchange_on_subset<VectorType>::value,
-                              VectorType>::type * = nullptr>
+    template <typename VectorType,
+              typename std::enable_if<has_compress_start<VectorType> &&
+                                        has_exchange_on_subset<VectorType>,
+                                      VectorType>::type * = nullptr>
     void
     compress_finish(const unsigned int component_in_block_vector,
                     VectorType &       vec)
@@ -3798,7 +3794,7 @@ namespace internal
      */
     template <
       typename VectorType,
-      typename std::enable_if<!has_exchange_on_subset<VectorType>::value &&
+      typename std::enable_if<!has_exchange_on_subset<VectorType> &&
                                 !is_not_parallel_vector<VectorType>::value,
                               VectorType>::type * = nullptr>
     void
@@ -3818,7 +3814,7 @@ namespace internal
      * LinearAlgebra::distributed::Vector
      */
     template <typename VectorType,
-              typename std::enable_if<has_exchange_on_subset<VectorType>::value,
+              typename std::enable_if<has_exchange_on_subset<VectorType>,
                                       VectorType>::type * = nullptr>
     void
     reset_ghost_values(const VectorType &vec) const
@@ -3862,7 +3858,7 @@ namespace internal
      * i.e. LinearAlgebra::distributed::Vector
      */
     template <typename VectorType,
-              typename std::enable_if<has_exchange_on_subset<VectorType>::value,
+              typename std::enable_if<has_exchange_on_subset<VectorType>,
                                       VectorType>::type * = nullptr>
     void
     zero_vector_region(const unsigned int range_index, VectorType &vec) const
@@ -3903,11 +3899,10 @@ namespace internal
      * subset of DoFs <==> begin() + ind == local_element(ind) but are still a
      * vector type
      */
-    template <
-      typename VectorType,
-      typename std::enable_if<!has_exchange_on_subset<VectorType>::value,
-                              VectorType>::type * = nullptr,
-      typename VectorType::value_type *           = nullptr>
+    template <typename VectorType,
+              typename std::enable_if<!has_exchange_on_subset<VectorType>,
+                                      VectorType>::type * = nullptr,
+              typename VectorType::value_type *           = nullptr>
     void
     zero_vector_region(const unsigned int range_index, VectorType &vec) const
     {
@@ -4002,10 +3997,9 @@ namespace internal
   // would be too many outstanding communication requests.
 
   // default value for vectors that do not have communication_block_size
-  template <
-    typename VectorStruct,
-    typename std::enable_if<!has_communication_block_size<VectorStruct>::value,
-                            VectorStruct>::type * = nullptr>
+  template <typename VectorStruct,
+            typename std::enable_if<!has_communication_block_size<VectorStruct>,
+                                    VectorStruct>::type * = nullptr>
   constexpr unsigned int
   get_communication_block_size(const VectorStruct &)
   {
@@ -4014,10 +4008,9 @@ namespace internal
 
 
 
-  template <
-    typename VectorStruct,
-    typename std::enable_if<has_communication_block_size<VectorStruct>::value,
-                            VectorStruct>::type * = nullptr>
+  template <typename VectorStruct,
+            typename std::enable_if<has_communication_block_size<VectorStruct>,
+                                    VectorStruct>::type * = nullptr>
   constexpr unsigned int
   get_communication_block_size(const VectorStruct &)
   {
index a03fd80f1b72ea86faff80a26858c22fb1909c53..f4c312f3c20c5ff1fbc258548e91f291429be9e2 100644 (file)
@@ -41,7 +41,7 @@ namespace internal
   using local_element_t = decltype(std::declval<T const>().local_element(0));
 
   template <typename T>
-  using has_local_element = is_detected<local_element_t, T>;
+  constexpr bool has_local_element = is_supported_operation<local_element_t, T>;
 
 
 
@@ -52,7 +52,8 @@ namespace internal
     decltype(std::declval<T>().add_local_element(0, typename T::value_type()));
 
   template <typename T>
-  using has_add_local_element = is_detected<add_local_element_t, T>;
+  constexpr bool has_add_local_element =
+    is_supported_operation<add_local_element_t, T>;
 
 
 
@@ -63,7 +64,8 @@ namespace internal
     decltype(std::declval<T>().set_local_element(0, typename T::value_type()));
 
   template <typename T>
-  using has_set_local_element = is_detected<set_local_element_t, T>;
+  constexpr bool has_set_local_element =
+    is_supported_operation<set_local_element_t, T>;
 
 
 
@@ -76,8 +78,8 @@ namespace internal
       std::declval<Utilities::MPI::Partitioner>()));
 
   template <typename T>
-  using has_partitioners_are_compatible =
-    is_detected<partitioners_are_compatible_t, T>;
+  constexpr bool has_partitioners_are_compatible =
+    is_supported_operation<partitioners_are_compatible_t, T>;
 
 
 
@@ -87,7 +89,7 @@ namespace internal
   using begin_t = decltype(std::declval<T const>().begin());
 
   template <typename T>
-  using has_begin = is_detected<begin_t, T>;
+  constexpr bool has_begin = is_supported_operation<begin_t, T>;
 
 
 
@@ -98,7 +100,8 @@ namespace internal
     decltype(std::declval<T const>().shared_vector_data());
 
   template <typename T>
-  using has_shared_vector_data = is_detected<shared_vector_data_t, T>;
+  constexpr bool has_shared_vector_data =
+    is_supported_operation<shared_vector_data_t, T>;
 
 
 
@@ -111,8 +114,8 @@ namespace internal
   struct is_vectorizable
   {
     static const bool value =
-      has_begin<T>::value &&
-      (has_local_element<T>::value ||
+      has_begin<T> &&
+      (has_local_element<T> ||
        is_serial_vector<typename std::remove_const<T>::type>::value) &&
       std::is_same<typename T::value_type, Number>::value;
   };
@@ -138,8 +141,8 @@ namespace internal
     decltype(std::declval<T const>().update_ghost_values_start(0));
 
   template <typename T>
-  using has_update_ghost_values_start =
-    is_detected<update_ghost_values_start_t, T>;
+  constexpr bool has_update_ghost_values_start =
+    is_supported_operation<update_ghost_values_start_t, T>;
 
 
 
@@ -150,7 +153,8 @@ namespace internal
     decltype(std::declval<T>().compress_start(0, VectorOperation::add));
 
   template <typename T>
-  using has_compress_start = is_detected<compress_start_t, T>;
+  constexpr bool has_compress_start =
+    is_supported_operation<compress_start_t, T>;
 
 
 
@@ -159,16 +163,8 @@ namespace internal
   // We assume that if both begin() and local_element()
   // exist, then begin() + offset == local_element(offset)
   template <typename T>
-  struct has_exchange_on_subset
-  {
-    static const bool value = has_begin<T>::value &&
-                              has_local_element<T>::value &&
-                              has_partitioners_are_compatible<T>::value;
-  };
-
-  // We need to have a separate declaration for static const members
-  template <typename T>
-  const bool has_exchange_on_subset<T>::value;
+  constexpr bool   has_exchange_on_subset =
+    has_begin<T> &&has_local_element<T> &&has_partitioners_are_compatible<T>;
 
 
 
@@ -178,8 +174,8 @@ namespace internal
   using communication_block_size_t = decltype(T::communication_block_size);
 
   template <typename T>
-  using has_communication_block_size =
-    is_detected<communication_block_size_t, T>;
+  constexpr bool has_communication_block_size =
+    is_supported_operation<communication_block_size_t, T>;
 
 
 
@@ -195,7 +191,7 @@ namespace internal
 
   template <class T>
   using is_not_parallel_vector =
-    detected_or_t<std::true_type, not_parallel_vector_t, T>;
+    SupportsOperation::detected_or_t<std::true_type, not_parallel_vector_t, T>;
 } // namespace internal
 #endif
 
index 7222b3321bdf155ea5eb2f142165c0f09c358d99..54346511b39190ac6531f9339a8ed1fc8af69be4 100644 (file)
@@ -36,7 +36,7 @@ namespace internal
   // FIXME: this is wrong for Trilinos/Petsc MPI vectors
   // where we should first do Partitioner::local_to_global()
   template <typename VectorType,
-            typename std::enable_if<!has_local_element<VectorType>::value,
+            typename std::enable_if<!has_local_element<VectorType>,
                                     VectorType>::type * = nullptr>
   inline typename VectorType::value_type
   vector_access(const VectorType &vec, const unsigned int entry)
@@ -50,7 +50,7 @@ namespace internal
   // FIXME: this is wrong for Trilinos/Petsc MPI vectors
   // where we should first do Partitioner::local_to_global()
   template <typename VectorType,
-            typename std::enable_if<!has_local_element<VectorType>::value,
+            typename std::enable_if<!has_local_element<VectorType>,
                                     VectorType>::type * = nullptr>
   inline typename VectorType::value_type &
   vector_access(VectorType &vec, const unsigned int entry)
@@ -64,7 +64,7 @@ namespace internal
   // method to access data in local index space, which is what we use in
   // DoFInfo and hence in read_dof_values etc.
   template <typename VectorType,
-            typename std::enable_if<has_local_element<VectorType>::value,
+            typename std::enable_if<has_local_element<VectorType>,
                                     VectorType>::type * = nullptr>
   inline typename VectorType::value_type &
   vector_access(VectorType &vec, const unsigned int entry)
@@ -76,7 +76,7 @@ namespace internal
 
   // same for const access
   template <typename VectorType,
-            typename std::enable_if<has_local_element<VectorType>::value,
+            typename std::enable_if<has_local_element<VectorType>,
                                     VectorType>::type * = nullptr>
   inline typename VectorType::value_type
   vector_access(const VectorType &vec, const unsigned int entry)
@@ -87,7 +87,7 @@ namespace internal
 
 
   template <typename VectorType,
-            typename std::enable_if<has_add_local_element<VectorType>::value,
+            typename std::enable_if<has_add_local_element<VectorType>,
                                     VectorType>::type * = nullptr>
   inline void
   vector_access_add(VectorType &                           vec,
@@ -100,7 +100,7 @@ namespace internal
 
 
   template <typename VectorType,
-            typename std::enable_if<!has_add_local_element<VectorType>::value,
+            typename std::enable_if<!has_add_local_element<VectorType>,
                                     VectorType>::type * = nullptr>
   inline void
   vector_access_add(VectorType &                           vec,
@@ -113,7 +113,7 @@ namespace internal
 
 
   template <typename VectorType,
-            typename std::enable_if<has_add_local_element<VectorType>::value,
+            typename std::enable_if<has_add_local_element<VectorType>,
                                     VectorType>::type * = nullptr>
   inline void
   vector_access_add_global(VectorType &                           vec,
@@ -126,7 +126,7 @@ namespace internal
 
 
   template <typename VectorType,
-            typename std::enable_if<!has_add_local_element<VectorType>::value,
+            typename std::enable_if<!has_add_local_element<VectorType>,
                                     VectorType>::type * = nullptr>
   inline void
   vector_access_add_global(VectorType &                           vec,
@@ -139,7 +139,7 @@ namespace internal
 
 
   template <typename VectorType,
-            typename std::enable_if<has_set_local_element<VectorType>::value,
+            typename std::enable_if<has_set_local_element<VectorType>,
                                     VectorType>::type * = nullptr>
   inline void
   vector_access_set(VectorType &                           vec,
@@ -152,7 +152,7 @@ namespace internal
 
 
   template <typename VectorType,
-            typename std::enable_if<!has_set_local_element<VectorType>::value,
+            typename std::enable_if<!has_set_local_element<VectorType>,
                                     VectorType>::type * = nullptr>
   inline void
   vector_access_set(VectorType &                           vec,
@@ -170,7 +170,7 @@ namespace internal
   // FIXME: this is incorrect for PETSc/Trilinos MPI vectors
   template <
     typename VectorType,
-    typename std::enable_if<!has_partitioners_are_compatible<VectorType>::value,
+    typename std::enable_if<!has_partitioners_are_compatible<VectorType>,
                             VectorType>::type * = nullptr>
   inline void
   check_vector_compatibility(
@@ -186,10 +186,9 @@ namespace internal
 
 
   // same as above for has_partitioners_are_compatible == true
-  template <
-    typename VectorType,
-    typename std::enable_if<has_partitioners_are_compatible<VectorType>::value,
-                            VectorType>::type * = nullptr>
+  template <typename VectorType,
+            typename std::enable_if<has_partitioners_are_compatible<VectorType>,
+                                    VectorType>::type * = nullptr>
   inline void
   check_vector_compatibility(
     const VectorType &                            vec,

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.