]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use standardized detection idiom
authorDaniel Arndt <arndtd@ornl.gov>
Mon, 24 Jan 2022 23:17:02 +0000 (18:17 -0500)
committerDaniel Arndt <arndtd@ornl.gov>
Tue, 25 Jan 2022 15:16:00 +0000 (10:16 -0500)
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/matrix_free.h
include/deal.II/matrix_free/type_traits.h

index 6a253619b089b3c71305e2242f7073866d2d3a65..31c4fab04182719503fc0f898568ec951870dfe5 100644 (file)
 
 DEAL_II_NAMESPACE_OPEN
 
+template <class...>
+using void_t = void;
+
+// Detection idiom from Version 2 of the C++ Extensions for Library
+// Fundamentals, ISO/IEC TS 19568:2017
+
+namespace internal
+{
+  // base class for nonesuch to inherit from so it is not an aggregate
+  struct nonesuch_base
+  {};
+
+  // primary template handles all types not supporting the archetypal Op
+  template <class Default,
+            class /*AlwaysVoid*/,
+            template <class...>
+            class Op,
+            class... /*Args*/>
+  struct detector
+  {
+    using value_t = std::false_type;
+    using type    = Default;
+  };
+
+  // specialization recognizes and handles only types supporting Op
+  template <class Default, template <class...> class Op, class... Args>
+  struct detector<Default, void_t<Op<Args...>>, Op, Args...>
+  {
+    using value_t = std::true_type;
+    using type    = Op<Args...>;
+  };
+} // namespace internal
+
+struct nonesuch : private internal::nonesuch_base
+{
+  ~nonesuch()                = delete;
+  nonesuch(nonesuch const &) = delete;
+  void
+  operator=(nonesuch const &) = delete;
+};
+
+template <class Default, template <class...> class Op, class... Args>
+using detected_or = internal::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 detected_t = typename detected_or<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 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 internal
 {
   namespace TemplateConstraints
@@ -98,23 +160,11 @@ struct enable_if_all
  * the `begin()` and `end()` functions, or is a C-style array.
  */
 template <typename T>
-class has_begin_and_end
-{
-  template <typename C>
-  static std::false_type
-  test(...);
-
-  template <typename C>
-  static auto
-  test(int) -> decltype(std::begin(std::declval<C>()),
-                        std::end(std::declval<C>()),
-                        std::true_type());
+using begin_and_end_t =
+  decltype(std::begin(std::declval<T>()), std::end(std::declval<T>()));
 
-public:
-  using type = decltype(test<T>(0));
-
-  static const bool value = type::value;
-};
+template <typename T>
+using has_begin_and_end = is_detected<begin_and_end_t, T>;
 
 
 
index 44d0161fbf2f24da8f00f26c5c0afae5a22c4042..c2e1d121c9cfea54796e25295a5012098797b17f 100644 (file)
@@ -1887,75 +1887,40 @@ namespace internal
       // A helper type-trait that leverage SFINAE to figure out if type T has
       // void T::get_mpi_communicator()
       template <typename T>
-      struct has_get_mpi_communicator
-      {
-      private:
-        static bool
-        detect(...);
-
-        template <typename U>
-        static decltype(std::declval<U>().get_mpi_communicator())
-        detect(const U &);
+      using get_mpi_communicator_t =
+        decltype(std::declval<T>().get_mpi_communicator());
 
-      public:
-        static const bool value =
-          !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-      };
+      template <typename T>
+      using has_get_mpi_communicator = is_detected<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()
       template <typename T>
-      struct has_locally_owned_domain_indices
-      {
-      private:
-        static bool
-        detect(...);
-
-        template <typename U>
-        static decltype(std::declval<U>().locally_owned_domain_indices())
-        detect(const U &);
+      using locally_owned_domain_indices_t =
+        decltype(std::declval<T>().locally_owned_domain_indices());
 
-      public:
-        static const bool value =
-          !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-      };
+      template <typename T>
+      using has_locally_owned_domain_indices =
+        is_detected<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()
       template <typename T>
-      struct has_locally_owned_range_indices
-      {
-      private:
-        static bool
-        detect(...);
-
-        template <typename U>
-        static decltype(std::declval<U>().locally_owned_range_indices())
-        detect(const U &);
+      using locally_owned_range_indices_t =
+        decltype(std::declval<T>().locally_owned_range_indices());
 
-      public:
-        static const bool value =
-          !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-      };
+      template <typename T>
+      using has_locally_owned_range_indices =
+        is_detected<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)
       template <typename T>
-      struct has_initialize_dof_vector
-      {
-      private:
-        static bool
-        detect(...);
-
-        template <typename U>
-        static decltype(std::declval<U>().initialize_dof_vector(
-          std::declval<LinearAlgebra::distributed::Vector<Number> &>()))
-        detect(const U &);
-
-      public:
-        static const bool value =
-          !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-      };
+      using initialize_dof_vector_t =
+        decltype(std::declval<T>().initialize_dof_vector());
+
+      template <typename T>
+      using has_initialize_dof_vector = is_detected<initialize_dof_vector_t, T>;
 
       // Used for (Trilinos/PETSc)Wrappers::SparseMatrix
       template <typename MatrixType,
index ef2383fecd2ee09150159477b007b63d9b78f8e1..6e0b1ebbb55040e934a52d688ecb0ecd88af0d31 100644 (file)
@@ -528,177 +528,73 @@ namespace internal
   namespace PreconditionRelaxation
   {
     template <typename T, typename VectorType>
-    struct has_Tvmult
-    {
-    private:
-      static bool
-      detect(...);
-
-      template <typename U>
-      static decltype(
-        std::declval<const U>().Tvmult(std::declval<VectorType &>(),
-                                       std::declval<const VectorType &>()))
-      detect(const U &);
-
-    public:
-      static const bool value =
-        !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-    };
+    using Tvmult_t = decltype(
+      std::declval<T const>().Tvmult(std::declval<VectorType &>(),
+                                     std::declval<const VectorType &>()));
 
     template <typename T, typename VectorType>
-    const bool has_Tvmult<T, VectorType>::value;
+    using has_Tvmult = is_detected<Tvmult_t, T, VectorType>;
 
     template <typename T, typename VectorType>
-    struct has_step
-    {
-    private:
-      static bool
-      detect(...);
-
-      template <typename U>
-      static decltype(
-        std::declval<const U>().step(std::declval<VectorType &>(),
-                                     std::declval<const VectorType &>()))
-      detect(const U &);
-
-    public:
-      static const bool value =
-        !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-    };
+    using step_t = decltype(
+      std::declval<T const>().step(std::declval<VectorType &>(),
+                                   std::declval<const VectorType &>()));
 
     template <typename T, typename VectorType>
-    const bool has_step<T, VectorType>::value;
+    using has_step = is_detected<step_t, T, VectorType>;
 
     template <typename T, typename VectorType>
-    struct has_step_omega
-    {
-    private:
-      static bool
-      detect(...);
-
-      template <typename U>
-      static decltype(
-        std::declval<const U>().step(std::declval<VectorType &>(),
-                                     std::declval<const VectorType &>(),
-                                     std::declval<const double>()))
-      detect(const U &);
-
-    public:
-      static const bool value =
-        !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-    };
+    using step_omega_t =
+      decltype(std::declval<T const>().step(std::declval<VectorType &>(),
+                                            std::declval<const VectorType &>(),
+                                            std::declval<const double>()));
 
     template <typename T, typename VectorType>
-    const bool has_step_omega<T, VectorType>::value;
+    using has_step_omega = is_detected<step_omega_t, T, VectorType>;
 
     template <typename T, typename VectorType>
-    struct has_Tstep
-    {
-    private:
-      static bool
-      detect(...);
-
-      template <typename U>
-      static decltype(
-        std::declval<const U>().Tstep(std::declval<VectorType &>(),
-                                      std::declval<const VectorType &>()))
-      detect(const U &);
-
-    public:
-      static const bool value =
-        !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-    };
+    using Tstep_t = decltype(
+      std::declval<T const>().Tstep(std::declval<VectorType &>(),
+                                    std::declval<const VectorType &>()));
 
     template <typename T, typename VectorType>
-    const bool has_Tstep<T, VectorType>::value;
+    using has_Tstep = is_detected<Tstep_t, T, VectorType>;
 
     template <typename T, typename VectorType>
-    struct has_Tstep_omega
-    {
-    private:
-      static bool
-      detect(...);
-
-      template <typename U>
-      static decltype(
-        std::declval<const U>().Tstep(std::declval<VectorType &>(),
-                                      std::declval<const VectorType &>(),
-                                      std::declval<const double>()))
-      detect(const U &);
-
-    public:
-      static const bool value =
-        !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-    };
+    using Tstep_omega_t =
+      decltype(std::declval<T const>().Tstep(std::declval<VectorType &>(),
+                                             std::declval<const VectorType &>(),
+                                             std::declval<const double>()));
 
     template <typename T, typename VectorType>
-    const bool has_Tstep_omega<T, VectorType>::value;
+    using has_Tstep_omega = is_detected<Tstep_omega_t, T, VectorType>;
 
     template <typename T, typename VectorType>
-    struct has_jacobi_step
-    {
-    private:
-      static bool
-      detect(...);
-
-      template <typename U>
-      static decltype(
-        std::declval<const U>().Jacobi_step(std::declval<VectorType &>(),
-                                            std::declval<const VectorType &>(),
-                                            1.0))
-      detect(const U &);
-
-    public:
-      static const bool value =
-        !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-    };
+    using jacobi_step_t = decltype(
+      std::declval<T const>().Jacobi_step(std::declval<VectorType &>(),
+                                          std::declval<const VectorType &>(),
+                                          std::declval<const double>()));
 
     template <typename T, typename VectorType>
-    const bool has_jacobi_step<T, VectorType>::value;
+    using has_jacobi_step = is_detected<jacobi_step_t, T, VectorType>;
 
     template <typename T, typename VectorType>
-    struct has_SOR_step
-    {
-    private:
-      static bool
-      detect(...);
-
-      template <typename U>
-      static decltype(
-        std::declval<const U>().SOR_step(std::declval<VectorType &>(),
-                                         std::declval<const VectorType &>(),
-                                         1.0))
-      detect(const U &);
-
-    public:
-      static const bool value =
-        !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-    };
+    using SOR_step_t = decltype(
+      std::declval<T const>().SOR_step(std::declval<VectorType &>(),
+                                       std::declval<const VectorType &>(),
+                                       std::declval<const double>()));
 
     template <typename T, typename VectorType>
-    const bool has_SOR_step<T, VectorType>::value;
+    using has_SOR_step = is_detected<SOR_step_t, T, VectorType>;
 
     template <typename T, typename VectorType>
-    struct has_SSOR_step
-    {
-    private:
-      static bool
-      detect(...);
-
-      template <typename U>
-      static decltype(
-        std::declval<const U>().SSOR_step(std::declval<VectorType &>(),
-                                          std::declval<const VectorType &>(),
-                                          1.0))
-      detect(const U &);
-
-    public:
-      static const bool value =
-        !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-    };
+    using SSOR_step_t = decltype(
+      std::declval<T const>().SSOR_step(std::declval<VectorType &>(),
+                                        std::declval<const VectorType &>(),
+                                        std::declval<const double>()));
 
     template <typename T, typename VectorType>
-    const bool has_SSOR_step<T, VectorType>::value;
+    using has_SSOR_step = is_detected<SSOR_step_t, T, VectorType>;
 
     template <typename MatrixType>
     class PreconditionJacobiImpl
index d7240ef755699c8e8ad99f896cdce9c31221f325..8731d802b7f11c71362958cd9922a382a4a3259e 100644 (file)
@@ -3359,7 +3359,7 @@ namespace internal
      * Start update_ghost_value for serial vectors
      */
     template <typename VectorType,
-              typename std::enable_if<is_serial_or_dummy<VectorType>::value,
+              typename std::enable_if<is_not_parallel_vector<VectorType>::value,
                                       VectorType>::type * = nullptr>
     void
     update_ghost_values_start(const unsigned int /*component_in_block_vector*/,
@@ -3374,7 +3374,7 @@ namespace internal
     template <typename VectorType,
               typename std::enable_if<
                 !has_update_ghost_values_start<VectorType>::value &&
-                  !is_serial_or_dummy<VectorType>::value,
+                  !is_not_parallel_vector<VectorType>::value,
                 VectorType>::type * = nullptr>
     void
     update_ghost_values_start(const unsigned int component_in_block_vector,
@@ -3580,7 +3580,7 @@ namespace internal
      * Start compress for serial vectors
      */
     template <typename VectorType,
-              typename std::enable_if<is_serial_or_dummy<VectorType>::value,
+              typename std::enable_if<is_not_parallel_vector<VectorType>::value,
                                       VectorType>::type * = nullptr>
     void
     compress_start(const unsigned int /*component_in_block_vector*/,
@@ -3593,10 +3593,11 @@ namespace internal
      * Start compress for vectors that do not support
      * the split into _start() and finish() stages
      */
-    template <typename VectorType,
-              typename std::enable_if<!has_compress_start<VectorType>::value &&
-                                        !is_serial_or_dummy<VectorType>::value,
-                                      VectorType>::type * = nullptr>
+    template <
+      typename VectorType,
+      typename std::enable_if<!has_compress_start<VectorType>::value &&
+                                !is_not_parallel_vector<VectorType>::value,
+                              VectorType>::type * = nullptr>
     void
     compress_start(const unsigned int component_in_block_vector,
                    VectorType &       vec)
@@ -3783,7 +3784,7 @@ namespace internal
      * Reset all ghost values for serial vectors
      */
     template <typename VectorType,
-              typename std::enable_if<is_serial_or_dummy<VectorType>::value,
+              typename std::enable_if<is_not_parallel_vector<VectorType>::value,
                                       VectorType>::type * = nullptr>
     void
     reset_ghost_values(const VectorType & /*vec*/) const
@@ -3798,7 +3799,7 @@ namespace internal
     template <
       typename VectorType,
       typename std::enable_if<!has_exchange_on_subset<VectorType>::value &&
-                                !is_serial_or_dummy<VectorType>::value,
+                                !is_not_parallel_vector<VectorType>::value,
                               VectorType>::type * = nullptr>
     void
     reset_ghost_values(const VectorType &vec) const
index cdb9e3f72882f1733e47e5e64d1728685c81fbde..a03fd80f1b72ea86faff80a26858c22fb1909c53 100644 (file)
@@ -38,83 +38,32 @@ namespace internal
   // a helper type-trait that leverage SFINAE to figure out if type T has
   // ... T::local_element() const
   template <typename T>
-  struct has_local_element
-  {
-  private:
-    // this will work always.
-    // we let it be void as we know T::local_element() (if exists) should
-    // certainly return something
-    static void
-    detect(...);
-
-    // this detecter will work only if we have "... T::local_element() const"
-    // and its return type will be the same as local_element(),
-    // that we expect to be T::value_type
-    template <typename U>
-    static decltype(std::declval<const U>().local_element(0))
-    detect(const U &);
-
-  public:
-    // finally here we check if our detector has non-void return type
-    // T::value_type. This will happen if compiler can use second detector,
-    // otherwise SFINAE let it work with the more general first one that is void
-    static const bool value =
-      !std::is_same<void, decltype(detect(std::declval<T>()))>::value;
-  };
+  using local_element_t = decltype(std::declval<T const>().local_element(0));
 
-  // We need to have a separate declaration for static const members
   template <typename T>
-  const bool has_local_element<T>::value;
+  using has_local_element = is_detected<local_element_t, T>;
 
 
 
   // a helper type-trait that leverage SFINAE to figure out if type T has
   // void T::add_local_element(const uint, const typename T::value_type)
   template <typename T>
-  struct has_add_local_element
-  {
-  private:
-    static int
-    detect(...);
-
-    template <typename U>
-    static decltype(
-      std::declval<U>().add_local_element(0, typename T::value_type()))
-    detect(const U &);
-
-  public:
-    static const bool value =
-      !std::is_same<int, decltype(detect(std::declval<T>()))>::value;
-  };
+  using add_local_element_t =
+    decltype(std::declval<T>().add_local_element(0, typename T::value_type()));
 
-  // We need to have a separate declaration for static const members
   template <typename T>
-  const bool has_add_local_element<T>::value;
+  using has_add_local_element = is_detected<add_local_element_t, T>;
 
 
 
   // a helper type-trait that leverage SFINAE to figure out if type T has
   // void T::set_local_element(const uint, const typename T::value_type)
   template <typename T>
-  struct has_set_local_element
-  {
-  private:
-    static int
-    detect(...);
-
-    template <typename U>
-    static decltype(
-      std::declval<U>().set_local_element(0, typename T::value_type()))
-    detect(const U &);
+  using set_local_element_t =
+    decltype(std::declval<T>().set_local_element(0, typename T::value_type()));
 
-  public:
-    static const bool value =
-      !std::is_same<int, decltype(detect(std::declval<T>()))>::value;
-  };
-
-  // We need to have a separate declaration for static const members
   template <typename T>
-  const bool has_set_local_element<T>::value;
+  using has_set_local_element = is_detected<set_local_element_t, T>;
 
 
 
@@ -122,71 +71,35 @@ namespace internal
   // bool T::partitioners_are_compatible(const Utilities::MPI::Partitioner &)
   // const
   template <typename T>
-  struct has_partitioners_are_compatible
-  {
-  private:
-    static void
-    detect(...);
-
-    template <typename U>
-    static decltype(std::declval<const U>().partitioners_are_compatible(
-      std::declval<Utilities::MPI::Partitioner>()))
-    detect(const U &);
-
-  public:
-    static const bool value =
-      std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-  };
+  using partitioners_are_compatible_t =
+    decltype(std::declval<T const>().partitioners_are_compatible(
+      std::declval<Utilities::MPI::Partitioner>()));
 
-  // We need to have a separate declaration for static const members
   template <typename T>
-  const bool has_partitioners_are_compatible<T>::value;
+  using has_partitioners_are_compatible =
+    is_detected<partitioners_are_compatible_t, T>;
+
 
 
   // same as above to check
   // ... T::begin() const
   template <typename T>
-  struct has_begin
-  {
-  private:
-    static void
-    detect(...);
+  using begin_t = decltype(std::declval<T const>().begin());
 
-    template <typename U>
-    static decltype(std::declval<const U>().begin())
-    detect(const U &);
-
-  public:
-    static const bool value =
-      !std::is_same<void, decltype(detect(std::declval<T>()))>::value;
-  };
-
-  // We need to have a separate declaration for static const members
   template <typename T>
-  const bool has_begin<T>::value;
+  using has_begin = is_detected<begin_t, T>;
+
 
 
   // same as above to check
   // ... T::shared_vector_data() const
   template <typename T>
-  struct has_shared_vector_data
-  {
-  private:
-    static void
-    detect(...);
-
-    template <typename U>
-    static decltype(std::declval<const U>().shared_vector_data())
-    detect(const U &);
+  using shared_vector_data_t =
+    decltype(std::declval<T const>().shared_vector_data());
 
-  public:
-    static const bool value =
-      !std::is_same<void, decltype(detect(std::declval<T>()))>::value;
-  };
-
-  // We need to have a separate declaration for static const members
   template <typename T>
-  const bool has_shared_vector_data<T>::value;
+  using has_shared_vector_data = is_detected<shared_vector_data_t, T>;
+
 
 
   // type trait for vector T and Number to see if
@@ -221,48 +134,23 @@ namespace internal
   // a helper type-trait that leverage SFINAE to figure out if type T has
   // void T::update_ghost_values_start(const uint) const
   template <typename T>
-  struct has_update_ghost_values_start
-  {
-  private:
-    static bool
-    detect(...);
-
-    template <typename U>
-    static decltype(std::declval<const U>().update_ghost_values_start(0))
-    detect(const U &);
-
-  public:
-    static const bool value =
-      !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-  };
+  using update_ghost_values_start_t =
+    decltype(std::declval<T const>().update_ghost_values_start(0));
 
-  // We need to have a separate declaration for static const members
   template <typename T>
-  const bool has_update_ghost_values_start<T>::value;
+  using has_update_ghost_values_start =
+    is_detected<update_ghost_values_start_t, T>;
 
 
 
   // a helper type-trait that leverage SFINAE to figure out if type T has
   // void T::  compress_start(const uint, VectorOperation::values)
   template <typename T>
-  struct has_compress_start
-  {
-  private:
-    static bool
-    detect(...);
-
-    template <typename U>
-    static decltype(std::declval<U>().compress_start(0, VectorOperation::add))
-    detect(const U &);
-
-  public:
-    static const bool value =
-      !std::is_same<decltype(detect(std::declval<T>())), bool>::value;
-  };
+  using compress_start_t =
+    decltype(std::declval<T>().compress_start(0, VectorOperation::add));
 
-  // We need to have a separate declaration for static const members
   template <typename T>
-  const bool has_compress_start<T>::value;
+  using has_compress_start = is_detected<compress_start_t, T>;
 
 
 
@@ -287,24 +175,11 @@ namespace internal
   // a helper type-trait that leverage SFINAE to figure out if type T has
   // T::communication_block_size
   template <typename T>
-  struct has_communication_block_size
-  {
-  private:
-    static void
-    detect(...);
-
-    template <typename U>
-    static decltype(U::communication_block_size)
-    detect(const U &);
-
-  public:
-    static const bool value =
-      !std::is_same<void, decltype(detect(std::declval<T>()))>::value;
-  };
+  using communication_block_size_t = decltype(T::communication_block_size);
 
-  // We need to have a separate declaration for static const members
   template <typename T>
-  const bool has_communication_block_size<T>::value;
+  using has_communication_block_size =
+    is_detected<communication_block_size_t, T>;
 
 
 
@@ -314,38 +189,13 @@ namespace internal
   // (like calculation of diagonals for matrix-free operators)
   // a dummy InVector == unsigned int is provided.
   // Thus we have to treat this case as well.
-  template <typename T>
-  struct is_serial_or_dummy
-  {
-  private:
-    // catches all cases including unsigned int
-    static void
-    detect(...);
-
-    // catches serial vectors
-    template <
-      typename U,
-      typename std::enable_if<is_serial_vector<U>::value, U>::type * = nullptr>
-    static void
-    detect(const U &);
-
-    // catches parallel vectors
-    template <
-      typename U,
-      typename std::enable_if<!is_serial_vector<U>::value, U>::type * = nullptr>
-    static bool
-    detect(const U &);
-
-  public:
-    static const bool value =
-      std::is_same<void, decltype(detect(std::declval<T>()))>::value;
-  };
-
-  // We need to have a separate declaration for static const members
-  template <typename T>
-  const bool is_serial_or_dummy<T>::value;
-
+  template <class T, class IsSerialVectorNotSpecialized = void>
+  using not_parallel_vector_t =
+    std::integral_constant<bool, is_serial_vector<T>::value>;
 
+  template <class T>
+  using is_not_parallel_vector =
+    detected_or_t<std::true_type, not_parallel_vector_t, T>;
 } // namespace internal
 #endif
 

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.