]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add internal type-traits to be used with VectorDataExchange of matrix-free framework
authorDenis Davydov <davydden@gmail.com>
Wed, 6 Mar 2019 16:52:18 +0000 (17:52 +0100)
committerDenis Davydov <davydden@gmail.com>
Fri, 8 Mar 2019 21:04:59 +0000 (22:04 +0100)
include/deal.II/matrix_free/fe_evaluation.h
include/deal.II/matrix_free/matrix_free.h
include/deal.II/matrix_free/type_traits.h [new file with mode: 0644]
tests/matrix_free/matrix_free_type_traits.cc [new file with mode: 0644]
tests/matrix_free/matrix_free_type_traits.with_trilinos=true.output [new file with mode: 0644]

index 5ae7b1905cd032b304d93db3480dc630759eaa6b..0a8e2c29e2599f21ab1ae06f655310edd20bdd1a 100644 (file)
@@ -35,6 +35,7 @@
 #include <deal.II/matrix_free/matrix_free.h>
 #include <deal.II/matrix_free/shape_info.h>
 #include <deal.II/matrix_free/tensor_product_kernels.h>
+#include <deal.II/matrix_free/type_traits.h>
 
 
 DEAL_II_NAMESPACE_OPEN
@@ -3401,155 +3402,7 @@ FEEvaluationBase<dim, n_components_, Number, is_face>::read_cell_data(
 
 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<U const>().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;
-  };
-
-  // We need to have a separate declaration for static const members
-  template <typename T>
-  const bool has_local_element<T>::value;
-
-
-
-  // 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;
-  };
-
-  // We need to have a separate declaration for static const members
-  template <typename T>
-  const bool has_add_local_element<T>::value;
-
-
-
-  // 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 &);
-
-  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;
-
-
-
-  // same as above to check
-  // 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<U const>().partitioners_are_compatible(
-      std::declval<Utilities::MPI::Partitioner>()))
-    detect(const U &);
-
-  public:
-    static const bool value =
-      std::is_same<bool, decltype(detect(std::declval<T>()))>::value;
-  };
-
-  // We need to have a separate declaration for static const members
-  template <typename T>
-  const bool has_partitioners_are_compatible<T>::value;
-
-
-  // same as above to check
-  // ... T::begin() const
-  template <typename T>
-  struct has_begin
-  {
-  private:
-    static void
-    detect(...);
-
-    template <typename U>
-    static decltype(std::declval<U const>().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;
-
-
-  // type trait for vector T and Number to see if
-  // we can do vectorized load/save.
-  // for VectorReader and VectorDistributorLocalToGlobal we assume that
-  // if both begin() and local_element()
-  // exist, then begin() + offset == local_element(offset)
-  template <typename T, typename Number>
-  struct is_vectorizable
-  {
-    static const bool value =
-      has_begin<T>::value &&
-      (has_local_element<T>::value || is_serial_vector<T>::value) &&
-      std::is_same<typename T::value_type, Number>::value;
-  };
-
-  // We need to have a separate declaration for static const members
-  template <typename T, typename Number>
-  const bool is_vectorizable<T, Number>::value;
-
+  // below we use type-traits from matrix-free/type_traits.h
 
   // access to generic const vectors that have operator ().
   // FIXME: this is wrong for Trilinos/Petsc MPI vectors
index e1d18fc083a71ddecda6d58f70a475060e3a9f2d..6397dec747d81c068783a5eebdd73d72e41c99aa 100644 (file)
@@ -44,6 +44,7 @@
 #include <deal.II/matrix_free/mapping_info.h>
 #include <deal.II/matrix_free/shape_info.h>
 #include <deal.II/matrix_free/task_info.h>
+#include <deal.II/matrix_free/type_traits.h>
 
 #include <cstdlib>
 #include <limits>
diff --git a/include/deal.II/matrix_free/type_traits.h b/include/deal.II/matrix_free/type_traits.h
new file mode 100644 (file)
index 0000000..70d784d
--- /dev/null
@@ -0,0 +1,328 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+#ifndef dealii_matrix_free_type_traits_h
+#define dealii_matrix_free_type_traits_h
+
+// various type-traits used exclusively within the matrix-free framework
+
+#include <deal.II/base/config.h>
+
+#include <deal.II/base/partitioner.h>
+
+#include <deal.II/lac/vector_type_traits.h>
+
+
+DEAL_II_NAMESPACE_OPEN
+
+
+namespace internal
+{
+  //
+  // type traits for FEEvaluation
+  //
+
+  // 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<U const>().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;
+  };
+
+  // We need to have a separate declaration for static const members
+  template <typename T>
+  const bool has_local_element<T>::value;
+
+
+
+  // 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;
+  };
+
+  // We need to have a separate declaration for static const members
+  template <typename T>
+  const bool has_add_local_element<T>::value;
+
+
+
+  // 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 &);
+
+  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;
+
+
+
+  // same as above to check
+  // 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<U const>().partitioners_are_compatible(
+      std::declval<Utilities::MPI::Partitioner>()))
+    detect(const U &);
+
+  public:
+    static const bool value =
+      std::is_same<bool, decltype(detect(std::declval<T>()))>::value;
+  };
+
+  // We need to have a separate declaration for static const members
+  template <typename T>
+  const bool has_partitioners_are_compatible<T>::value;
+
+
+  // same as above to check
+  // ... T::begin() const
+  template <typename T>
+  struct has_begin
+  {
+  private:
+    static void
+    detect(...);
+
+    template <typename U>
+    static decltype(std::declval<U const>().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;
+
+
+  // type trait for vector T and Number to see if
+  // we can do vectorized load/save.
+  // for VectorReader and VectorDistributorLocalToGlobal we assume that
+  // if both begin() and local_element()
+  // exist, then begin() + offset == local_element(offset)
+  template <typename T, typename Number>
+  struct is_vectorizable
+  {
+    static const bool value =
+      has_begin<T>::value &&
+      (has_local_element<T>::value || is_serial_vector<T>::value) &&
+      std::is_same<typename T::value_type, Number>::value;
+  };
+
+  // We need to have a separate declaration for static const members
+  template <typename T, typename Number>
+  const bool is_vectorizable<T, Number>::value;
+
+
+  //
+  // type-traits for Matrix-Free
+  //
+  // similar to type traits in FEEvaluation, below we add type-traits
+  // to distinguish between vectors that provide different interface for
+  // operations like update_ghost_values(), compress(), etc.
+  // see internal::has_local_element in fe_evaluation.h that documents
+  // how those type traits work.
+
+  // 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<U const>().update_ghost_values_start(0))
+    detect(const U &);
+
+  public:
+    static const bool value =
+      !std::is_same<bool, decltype(detect(std::declval<T>()))>::value;
+  };
+
+  // We need to have a separate declaration for static const members
+  template <typename T>
+  const bool has_update_ghost_values_start<T>::value;
+
+
+
+  // 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<bool, decltype(detect(std::declval<T>()))>::value;
+  };
+
+  // We need to have a separate declaration for static const members
+  template <typename T>
+  const bool has_compress_start<T>::value;
+
+
+
+  // type trait for vector T to see if
+  // we do a custom data exchange route.
+  // 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;
+  };
+
+  // We need to have a separate declaration for static const members
+  template <typename T>
+  const bool has_exchange_on_subset<T>::value;
+
+
+
+  // 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;
+  };
+
+  // We need to have a separate declaration for static const members
+  template <typename T>
+  const bool has_communication_block_size<T>::value;
+
+
+
+  // type trait for vector T to see if
+  // we need to do any data exchange for this vector type at all.
+  // is_serial_vector<> would have been enough, but in some circumstances
+  // (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;
+
+
+} // namespace internal
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
diff --git a/tests/matrix_free/matrix_free_type_traits.cc b/tests/matrix_free/matrix_free_type_traits.cc
new file mode 100644 (file)
index 0000000..fb0f597
--- /dev/null
@@ -0,0 +1,94 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+
+// test internal typetraits used in matrix_free.h
+
+#include <deal.II/lac/la_parallel_block_vector.h>
+#include <deal.II/lac/la_parallel_vector.h>
+#include <deal.II/lac/trilinos_vector.h>
+
+#include <deal.II/matrix_free/matrix_free.h>
+
+#include "../tests.h"
+
+int
+main()
+{
+  initlog();
+
+  deallog << "has_update_ghost_values_start:" << std::endl
+          << "LinearAlgebra::distributed::Vector = "
+          << internal::has_update_ghost_values_start<
+               LinearAlgebra::distributed::Vector<double>>::value
+          << std::endl
+          << "TrilinosWrappers::MPI::Vector = "
+          << internal::has_update_ghost_values_start<
+               TrilinosWrappers::MPI::Vector>::value
+          << std::endl
+          << "Vector = "
+          << internal::has_update_ghost_values_start<Vector<double>>::value
+          << std::endl;
+
+  deallog << "has_compress_start:" << std::endl
+          << "LinearAlgebra::distributed::Vector = "
+          << internal::has_compress_start<
+               LinearAlgebra::distributed::Vector<double>>::value
+          << std::endl
+          << "TrilinosWrappers::MPI::Vector = "
+          << internal::has_compress_start<TrilinosWrappers::MPI::Vector>::value
+          << std::endl
+          << "Vector = " << internal::has_compress_start<Vector<double>>::value
+          << std::endl;
+
+  deallog
+    << "has_exchange_on_subset:" << std::endl
+    << "LinearAlgebra::distributed::Vector = "
+    << internal::has_exchange_on_subset<
+         LinearAlgebra::distributed::Vector<double>>::value
+    << std::endl
+    << "TrilinosWrappers::MPI::Vector = "
+    << internal::has_exchange_on_subset<TrilinosWrappers::MPI::Vector>::value
+    << std::endl
+    << "Vector = " << internal::has_exchange_on_subset<Vector<double>>::value
+    << std::endl;
+
+
+  deallog << "has_communication_block_size:" << std::endl
+          << "LinearAlgebra::distributed::Vector = "
+          << internal::has_communication_block_size<
+               LinearAlgebra::distributed::Vector<double>>::value
+          << std::endl
+          << "LinearAlgebra::distributed::BlockVector = "
+          << internal::has_communication_block_size<
+               LinearAlgebra::distributed::BlockVector<double>>::value
+          << std::endl;
+
+  deallog << "is_serial_or_dummy:" << std::endl
+          << "LinearAlgebra::distributed::Vector = "
+          << internal::is_serial_or_dummy<
+               LinearAlgebra::distributed::Vector<double>>::value
+          << std::endl
+          << "TrilinosWrappers::MPI::Vector = "
+          << internal::is_serial_or_dummy<TrilinosWrappers::MPI::Vector>::value
+          << std::endl
+          << "Vector = " << internal::is_serial_or_dummy<Vector<double>>::value
+          << std::endl
+          << "unsigned int = "
+          << internal::is_serial_or_dummy<unsigned int>::value << std::endl;
+
+  deallog << "OK" << std::endl;
+}
diff --git a/tests/matrix_free/matrix_free_type_traits.with_trilinos=true.output b/tests/matrix_free/matrix_free_type_traits.with_trilinos=true.output
new file mode 100644 (file)
index 0000000..0382073
--- /dev/null
@@ -0,0 +1,22 @@
+
+DEAL::has_update_ghost_values_start:
+DEAL::LinearAlgebra::distributed::Vector = 1
+DEAL::TrilinosWrappers::MPI::Vector = 0
+DEAL::Vector = 0
+DEAL::has_compress_start:
+DEAL::LinearAlgebra::distributed::Vector = 1
+DEAL::TrilinosWrappers::MPI::Vector = 0
+DEAL::Vector = 0
+DEAL::has_exchange_on_subset:
+DEAL::LinearAlgebra::distributed::Vector = 1
+DEAL::TrilinosWrappers::MPI::Vector = 0
+DEAL::Vector = 0
+DEAL::has_communication_block_size:
+DEAL::LinearAlgebra::distributed::Vector = 0
+DEAL::LinearAlgebra::distributed::BlockVector = 1
+DEAL::is_serial_or_dummy:
+DEAL::LinearAlgebra::distributed::Vector = 0
+DEAL::TrilinosWrappers::MPI::Vector = 0
+DEAL::Vector = 1
+DEAL::unsigned int = 1
+DEAL::OK

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.