]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Introduce signed_global_dof_index 14811/head
authorDaniel Arndt <arndtd@ornl.gov>
Sat, 18 Feb 2023 04:11:09 +0000 (23:11 -0500)
committerDaniel Arndt <arndtd@ornl.gov>
Mon, 20 Feb 2023 12:58:46 +0000 (07:58 -0500)
cmake/configure/configure_20_trilinos.cmake
include/deal.II/base/index_set.h
include/deal.II/base/types.h
include/deal.II/lac/read_write_vector.h
include/deal.II/lac/read_write_vector.templates.h
include/deal.II/lac/trilinos_tpetra_communication_pattern.h
include/deal.II/lac/trilinos_tpetra_vector.h
include/deal.II/lac/trilinos_tpetra_vector.templates.h
include/deal.II/lac/vector_element_access.h
source/base/index_set.cc
source/lac/trilinos_tpetra_communication_pattern.cc

index 3655c55c9225ea62cd63b7c95627c3ec97c2a5cb..b13a0a72243c94e3ffb7b09b25b74d72cae46473 100644 (file)
@@ -241,9 +241,9 @@ macro(feature_trilinos_find_external var)
       add_flags(CMAKE_REQUIRED_FLAGS "${DEAL_II_CXX_FLAGS} ${DEAL_II_LINKER_FLAGS}")
 
       if(DEAL_II_WITH_64BIT_INDICES)
-        set(_global_index_type "std::uint64_t")
+        set(_global_index_type "std::int64_t")
       else()
-        set(_global_index_type "unsigned int")
+        set(_global_index_type "int")
       endif()
 
       CHECK_CXX_SOURCE_COMPILES(
@@ -253,7 +253,7 @@ macro(feature_trilinos_find_external var)
         int main()
         {
           using LO       = int;
-          using GO       = std::make_signed_t<${_global_index_type}>;
+          using GO       = ${_global_index_type};
           using map_type = Tpetra::Map<LO, GO>;
           Teuchos::RCP<const map_type>   dummy_map = Teuchos::rcp(new map_type());
           Tpetra::Vector<double, LO, GO> dummy_vector(dummy_map);
index 722a3010ad1f2a616b8e123d250e3b8600dc116a..18be4ffe8ac089a3f3ea9dff61c201c70452a561 100644 (file)
@@ -496,7 +496,7 @@ public:
                     const bool      overlapping  = false) const;
 
 #  ifdef DEAL_II_TRILINOS_WITH_TPETRA
-  Tpetra::Map<int, std::make_signed_t<types::global_dof_index>>
+  Tpetra::Map<int, types::signed_global_dof_index>
   make_tpetra_map(const MPI_Comm &communicator = MPI_COMM_WORLD,
                   const bool      overlapping  = false) const;
 #  endif
index 6bb80727087a0f7b723776b0cc832ed6ca99f199..88c2c2b613fb4a8519f1a3d3f11e22472b9fd880 100644 (file)
@@ -20,6 +20,7 @@
 #include <deal.II/base/config.h>
 
 #include <cstdint>
+#include <type_traits> // make_signed_t
 
 
 DEAL_II_NAMESPACE_OPEN
@@ -81,6 +82,13 @@ namespace types
   using global_dof_index  = unsigned int;
 #endif
 
+  /**
+   * Signed version of global_dof_index.
+   * This is useful for interacting with Trilinos' Tpetra that only works well
+   * with a signed global ordinal type.
+   */
+  using signed_global_dof_index = std::make_signed_t<global_dof_index>;
+
   /**
    * An identifier that denotes the MPI type associated with
    * types::global_dof_index.
index 6236f67e0c9ae45ba065dc9ed0ca51230573d05a..62358e09177bb6e0cf8d92aaa59faf4ae2b5b583 100644 (file)
@@ -657,9 +657,8 @@ namespace LinearAlgebra
      * used directly.
      */
     void
-    import(const Tpetra::
-             Vector<Number, int, std::make_signed_t<types::global_dof_index>>
-               &                   tpetra_vector,
+    import(const Tpetra::Vector<Number, int, types::signed_global_dof_index>
+             &                     tpetra_vector,
            const IndexSet &        locally_owned_elements,
            VectorOperation::values operation,
            const MPI_Comm &        mpi_comm,
index a85caa9a18ffc8da84c232c20d9d27d39d648286..66e4f5a276e4b77ff52d1dd38a91a3e543e3c3b0 100644 (file)
@@ -577,8 +577,7 @@ namespace LinearAlgebra
   template <typename Number>
   void
   ReadWriteVector<Number>::import(
-    const Tpetra::
-      Vector<Number, int, std::make_signed_t<types::global_dof_index>> &vector,
+    const Tpetra::Vector<Number, int, types::signed_global_dof_index> &vector,
     const IndexSet &        source_elements,
     VectorOperation::values operation,
     const MPI_Comm &        mpi_comm,
@@ -621,11 +620,11 @@ namespace LinearAlgebra
                       "LinearAlgebra::TpetraWrappers::CommunicationPattern."));
       }
 
-    Tpetra::Export<int, std::make_signed_t<types::global_dof_index>>
-      tpetra_export(tpetra_comm_pattern->get_tpetra_export());
+    Tpetra::Export<int, types::signed_global_dof_index> tpetra_export(
+      tpetra_comm_pattern->get_tpetra_export());
 
-    Tpetra::Vector<Number, int, std::make_signed_t<types::global_dof_index>>
-      target_vector(tpetra_export.getSourceMap());
+    Tpetra::Vector<Number, int, types::signed_global_dof_index> target_vector(
+      tpetra_export.getSourceMap());
     target_vector.doImport(vector, tpetra_export, Tpetra::REPLACE);
 
     const auto *new_values = target_vector.getData().get();
index 862dadc24d05c56521f9760e3ad384f4dfa549c3..927c9353704608e7248df09f05278504ded0c485 100644 (file)
@@ -68,13 +68,13 @@ namespace LinearAlgebra
       /**
        * Return the underlying Tpetra::Import object.
        */
-      const Tpetra::Import<int, std::make_signed_t<types::global_dof_index>> &
+      const Tpetra::Import<int, types::signed_global_dof_index> &
       get_tpetra_import() const;
 
       /**
        * Return the underlying Tpetra::Export object.
        */
-      const Tpetra::Export<int, std::make_signed_t<types::global_dof_index>> &
+      const Tpetra::Export<int, types::signed_global_dof_index> &
       get_tpetra_export() const;
 
     private:
@@ -86,15 +86,13 @@ namespace LinearAlgebra
       /**
        * Shared pointer to the Tpetra::Import object used.
        */
-      std::unique_ptr<
-        Tpetra::Import<int, std::make_signed_t<types::global_dof_index>>>
+      std::unique_ptr<Tpetra::Import<int, types::signed_global_dof_index>>
         tpetra_import;
 
       /**
        * Shared pointer to the Tpetra::Export object used.
        */
-      std::unique_ptr<
-        Tpetra::Export<int, std::make_signed_t<types::global_dof_index>>>
+      std::unique_ptr<Tpetra::Export<int, types::signed_global_dof_index>>
         tpetra_export;
     };
   } // end of namespace TpetraWrappers
index 6e3bdab2f735542f89099d7f90dfb0adab30c6fd..6cf689122b56cfb923a54559c1a11d11841385f5 100644 (file)
@@ -337,15 +337,14 @@ namespace LinearAlgebra
        * Return a const reference to the underlying Trilinos
        * Tpetra::Vector class.
        */
-      const Tpetra::
-        Vector<Number, int, std::make_signed_t<types::global_dof_index>> &
-        trilinos_vector() const;
+      const Tpetra::Vector<Number, int, types::signed_global_dof_index> &
+      trilinos_vector() const;
 
       /**
        * Return a (modifiable) reference to the underlying Trilinos
        * Tpetra::Vector class.
        */
-      Tpetra::Vector<Number, int, std::make_signed_t<types::global_dof_index>> &
+      Tpetra::Vector<Number, int, types::signed_global_dof_index> &
       trilinos_vector();
 
       /**
@@ -400,8 +399,7 @@ namespace LinearAlgebra
        * Pointer to the actual Tpetra vector object.
        */
       std::unique_ptr<
-        Tpetra::
-          Vector<Number, int, std::make_signed_t<types::global_dof_index>>>
+        Tpetra::Vector<Number, int, types::signed_global_dof_index>>
         vector;
 
       /**
index 54db805e3866275a66923d4a2c48e072fd1d7df0..c937f9292ad511de05d9b7852131477e42519cab 100644 (file)
@@ -47,12 +47,9 @@ namespace LinearAlgebra
     template <typename Number>
     Vector<Number>::Vector()
       : Subscriptor()
-      , vector(new Tpetra::Vector<Number,
-                                  int,
-                                  std::make_signed_t<types::global_dof_index>>(
-          Teuchos::RCP<
-            Tpetra::Map<int, std::make_signed_t<types::global_dof_index>>>(
-            new Tpetra::Map<int, std::make_signed_t<types::global_dof_index>>(
+      , vector(new Tpetra::Vector<Number, int, types::signed_global_dof_index>(
+          Teuchos::RCP<Tpetra::Map<int, types::signed_global_dof_index>>(
+            new Tpetra::Map<int, types::signed_global_dof_index>(
               0,
               0,
               Utilities::Trilinos::tpetra_comm_self()))))
@@ -63,9 +60,7 @@ namespace LinearAlgebra
     template <typename Number>
     Vector<Number>::Vector(const Vector<Number> &V)
       : Subscriptor()
-      , vector(new Tpetra::Vector<Number,
-                                  int,
-                                  std::make_signed_t<types::global_dof_index>>(
+      , vector(new Tpetra::Vector<Number, int, types::signed_global_dof_index>(
           V.trilinos_vector(),
           Teuchos::Copy))
     {}
@@ -76,12 +71,9 @@ namespace LinearAlgebra
     Vector<Number>::Vector(const IndexSet &parallel_partitioner,
                            const MPI_Comm &communicator)
       : Subscriptor()
-      , vector(new Tpetra::Vector<Number,
-                                  int,
-                                  std::make_signed_t<types::global_dof_index>>(
-          Teuchos::rcp(
-            new Tpetra::Map<int, std::make_signed_t<types::global_dof_index>>(
-              parallel_partitioner.make_tpetra_map(communicator, false)))))
+      , vector(new Tpetra::Vector<Number, int, types::signed_global_dof_index>(
+          Teuchos::rcp(new Tpetra::Map<int, types::signed_global_dof_index>(
+            parallel_partitioner.make_tpetra_map(communicator, false)))))
     {}
 
 
@@ -92,15 +84,13 @@ namespace LinearAlgebra
                            const MPI_Comm &communicator,
                            const bool      omit_zeroing_entries)
     {
-      Tpetra::Map<int, std::make_signed_t<types::global_dof_index>> input_map =
+      Tpetra::Map<int, types::signed_global_dof_index> input_map =
         parallel_partitioner.make_tpetra_map(communicator, false);
       if (vector->getMap()->isSameAs(input_map) == false)
         vector = std::make_unique<
-          Tpetra::
-            Vector<Number, int, std::make_signed_t<types::global_dof_index>>>(
+          Tpetra::Vector<Number, int, types::signed_global_dof_index>>(
           Teuchos::rcp(
-            new Tpetra::Map<int, std::make_signed_t<types::global_dof_index>>(
-              input_map)));
+            new Tpetra::Map<int, types::signed_global_dof_index>(input_map)));
       else if (omit_zeroing_entries == false)
         {
           vector->putScalar(0.);
@@ -142,8 +132,8 @@ namespace LinearAlgebra
         {
           if (size() == V.size())
             {
-              Tpetra::Import<int, std::make_signed_t<types::global_dof_index>>
-                data_exchange(vector->getMap(), V.trilinos_vector().getMap());
+              Tpetra::Import<int, types::signed_global_dof_index> data_exchange(
+                vector->getMap(), V.trilinos_vector().getMap());
 
               vector->doImport(V.trilinos_vector(),
                                data_exchange,
@@ -151,9 +141,7 @@ namespace LinearAlgebra
             }
           else
             vector = std::make_unique<
-              Tpetra::Vector<Number,
-                             int,
-                             std::make_signed_t<types::global_dof_index>>>(
+              Tpetra::Vector<Number, int, types::signed_global_dof_index>>(
               V.trilinos_vector());
         }
 
@@ -214,10 +202,10 @@ namespace LinearAlgebra
               "LinearAlgebra::TpetraWrappers::CommunicationPattern."));
         }
 
-      Tpetra::Export<int, std::make_signed_t<types::global_dof_index>>
-        tpetra_export(tpetra_comm_pattern->get_tpetra_export());
-      Tpetra::Vector<Number, int, std::make_signed_t<types::global_dof_index>>
-        source_vector(tpetra_export.getSourceMap());
+      Tpetra::Export<int, types::signed_global_dof_index> tpetra_export(
+        tpetra_comm_pattern->get_tpetra_export());
+      Tpetra::Vector<Number, int, types::signed_global_dof_index> source_vector(
+        tpetra_export.getSourceMap());
 
       {
         source_vector.template sync<Kokkos::HostSpace>();
@@ -229,9 +217,8 @@ namespace LinearAlgebra
         for (size_t k = 0; k < localLength; ++k)
           x_1d(k) = *values_it++;
         source_vector.template sync<
-          typename Tpetra::
-            Vector<Number, int, std::make_signed_t<types::global_dof_index>>::
-              device_type::memory_space>();
+          typename Tpetra::Vector<Number, int, types::signed_global_dof_index>::
+            device_type::memory_space>();
       }
       if (operation == VectorOperation::insert)
         vector->doExport(source_vector, tpetra_export, Tpetra::REPLACE);
@@ -290,11 +277,10 @@ namespace LinearAlgebra
 
           // TODO: Tpetra doesn't have a combine mode that also updates local
           // elements, maybe there is a better workaround.
-          Tpetra::
-            Vector<Number, int, std::make_signed_t<types::global_dof_index>>
-              dummy(vector->getMap(), false);
-          Tpetra::Import<int, std::make_signed_t<types::global_dof_index>>
-            data_exchange(down_V.trilinos_vector().getMap(), dummy.getMap());
+          Tpetra::Vector<Number, int, types::signed_global_dof_index> dummy(
+            vector->getMap(), false);
+          Tpetra::Import<int, types::signed_global_dof_index> data_exchange(
+            down_V.trilinos_vector().getMap(), dummy.getMap());
 
           dummy.doImport(down_V.trilinos_vector(),
                          data_exchange,
@@ -355,9 +341,8 @@ namespace LinearAlgebra
           vector_1d(k) += a;
         }
       vector->template sync<
-        typename Tpetra::
-          Vector<Number, int, std::make_signed_t<types::global_dof_index>>::
-            device_type::memory_space>();
+        typename Tpetra::Vector<Number, int, types::signed_global_dof_index>::
+          device_type::memory_space>();
     }
 
 
@@ -621,9 +606,8 @@ namespace LinearAlgebra
 
 
     template <typename Number>
-    const Tpetra::
-      Vector<Number, int, std::make_signed_t<types::global_dof_index>> &
-      Vector<Number>::trilinos_vector() const
+    const Tpetra::Vector<Number, int, types::signed_global_dof_index> &
+    Vector<Number>::trilinos_vector() const
     {
       return *vector;
     }
@@ -631,7 +615,7 @@ namespace LinearAlgebra
 
 
     template <typename Number>
-    Tpetra::Vector<Number, int, std::make_signed_t<types::global_dof_index>> &
+    Tpetra::Vector<Number, int, types::signed_global_dof_index> &
     Vector<Number>::trilinos_vector()
     {
       return *vector;
index 2d2307769a942b02f8cda50082feb7c96a9b16cc..62332e8a04d9771d18010206a0c432ddadd4bd3a 100644 (file)
@@ -157,8 +157,8 @@ namespace internal
     LinearAlgebra::TpetraWrappers::Vector<NumberType> &V)
   {
     // Extract local indices in the vector.
-    Tpetra::Vector<NumberType, int, std::make_signed_t<types::global_dof_index>>
-                                      vector = V.trilinos_vector();
+    Tpetra::Vector<NumberType, int, types::signed_global_dof_index> vector =
+      V.trilinos_vector();
     TrilinosWrappers::types::int_type trilinos_i =
       vector.getMap()->getLocalElement(
         static_cast<TrilinosWrappers::types::int_type>(i));
@@ -170,9 +170,8 @@ namespace internal
     vector.template modify<Kokkos::HostSpace>();
     vector_1d(trilinos_i) += value;
     vector.template sync<
-      typename Tpetra::
-        Vector<NumberType, int, std::make_signed_t<types::global_dof_index>>::
-          device_type::memory_space>();
+      typename Tpetra::Vector<NumberType, int, types::signed_global_dof_index>::
+        device_type::memory_space>();
   }
 
 
@@ -185,8 +184,8 @@ namespace internal
     LinearAlgebra::TpetraWrappers::Vector<NumberType> &V)
   {
     // Extract local indices in the vector.
-    Tpetra::Vector<NumberType, int, std::make_signed_t<types::global_dof_index>>
-                                      vector = V.trilinos_vector();
+    Tpetra::Vector<NumberType, int, types::signed_global_dof_index> vector =
+      V.trilinos_vector();
     TrilinosWrappers::types::int_type trilinos_i =
       vector.getMap()->getLocalElement(
         static_cast<TrilinosWrappers::types::int_type>(i));
@@ -198,9 +197,8 @@ namespace internal
     vector.template modify<Kokkos::HostSpace>();
     vector_1d(trilinos_i) = value;
     vector.template sync<
-      typename Tpetra::
-        Vector<NumberType, int, std::make_signed_t<types::global_dof_index>>::
-          device_type::memory_space>();
+      typename Tpetra::Vector<NumberType, int, types::signed_global_dof_index>::
+        device_type::memory_space>();
   }
 
 
@@ -212,8 +210,8 @@ namespace internal
     const types::global_dof_index                            i)
   {
     // Extract local indices in the vector.
-    Tpetra::Vector<NumberType, int, std::make_signed_t<types::global_dof_index>>
-                                      vector = V.trilinos_vector();
+    Tpetra::Vector<NumberType, int, types::signed_global_dof_index> vector =
+      V.trilinos_vector();
     TrilinosWrappers::types::int_type trilinos_i =
       vector.getMap()->getLocalElement(
         static_cast<TrilinosWrappers::types::int_type>(i));
index 81e38f538f8e39843881cb8bc83b0e48b0de1e8a..630d8ba445757cea7b6124d0495452c618a0b0b3 100644 (file)
@@ -716,7 +716,7 @@ IndexSet::fill_index_vector(std::vector<size_type> &indices) const
 #ifdef DEAL_II_WITH_TRILINOS
 #  ifdef DEAL_II_TRILINOS_WITH_TPETRA
 
-Tpetra::Map<int, std::make_signed_t<types::global_dof_index>>
+Tpetra::Map<int, types::signed_global_dof_index>
 IndexSet::make_tpetra_map(const MPI_Comm &communicator,
                           const bool      overlapping) const
 {
@@ -752,7 +752,7 @@ IndexSet::make_tpetra_map(const MPI_Comm &communicator,
   const bool linear =
     overlapping ? false : is_ascending_and_one_to_one(communicator);
   if (linear)
-    return Tpetra::Map<int, std::make_signed_t<types::global_dof_index>>(
+    return Tpetra::Map<int, types::signed_global_dof_index>(
       size(),
       n_elements(),
       0,
@@ -764,13 +764,12 @@ IndexSet::make_tpetra_map(const MPI_Comm &communicator,
     );
   else
     {
-      const std::vector<size_type> indices = get_index_vector();
-      std::vector<std::make_signed_t<types::global_dof_index>> int_indices(
-        indices.size());
+      const std::vector<size_type>                indices = get_index_vector();
+      std::vector<types::signed_global_dof_index> int_indices(indices.size());
       std::copy(indices.begin(), indices.end(), int_indices.begin());
-      const Teuchos::ArrayView<std::make_signed_t<types::global_dof_index>>
-        arr_view(int_indices);
-      return Tpetra::Map<int, std::make_signed_t<types::global_dof_index>>(
+      const Teuchos::ArrayView<types::signed_global_dof_index> arr_view(
+        int_indices);
+      return Tpetra::Map<int, types::signed_global_dof_index>(
         size(),
         arr_view,
         0,
index 34c08d633279147c7759eed6978f0d2e5a4d1b3a..795ce6818527a83c57ff59734c79894055ad2830 100644 (file)
@@ -51,22 +51,22 @@ namespace LinearAlgebra
     {
       comm = std::make_shared<const MPI_Comm>(communicator);
 
-      auto vector_space_vector_map = Teuchos::rcp(
-        new Tpetra::Map<int, std::make_signed_t<types::global_dof_index>>(
+      auto vector_space_vector_map =
+        Teuchos::rcp(new Tpetra::Map<int, types::signed_global_dof_index>(
           vector_space_vector_index_set.make_tpetra_map(*comm, false)));
-      auto read_write_vector_map = Teuchos::rcp(
-        new Tpetra::Map<int, std::make_signed_t<types::global_dof_index>>(
+      auto read_write_vector_map =
+        Teuchos::rcp(new Tpetra::Map<int, types::signed_global_dof_index>(
           read_write_vector_index_set.make_tpetra_map(*comm, true)));
 
       // Target map is read_write_vector_map
       // Source map is vector_space_vector_map. This map must have uniquely
       // owned GID.
-      tpetra_import = std::make_unique<
-        Tpetra::Import<int, std::make_signed_t<types::global_dof_index>>>(
-        read_write_vector_map, vector_space_vector_map);
-      tpetra_export = std::make_unique<
-        Tpetra::Export<int, std::make_signed_t<types::global_dof_index>>>(
-        read_write_vector_map, vector_space_vector_map);
+      tpetra_import =
+        std::make_unique<Tpetra::Import<int, types::signed_global_dof_index>>(
+          read_write_vector_map, vector_space_vector_map);
+      tpetra_export =
+        std::make_unique<Tpetra::Export<int, types::signed_global_dof_index>>(
+          read_write_vector_map, vector_space_vector_map);
     }
 
 
@@ -79,7 +79,7 @@ namespace LinearAlgebra
 
 
 
-    const Tpetra::Import<int, std::make_signed_t<types::global_dof_index>> &
+    const Tpetra::Import<int, types::signed_global_dof_index> &
     CommunicationPattern::get_tpetra_import() const
     {
       return *tpetra_import;
@@ -87,7 +87,7 @@ namespace LinearAlgebra
 
 
 
-    const Tpetra::Export<int, std::make_signed_t<types::global_dof_index>> &
+    const Tpetra::Export<int, types::signed_global_dof_index> &
     CommunicationPattern::get_tpetra_export() const
     {
       return *tpetra_export;

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.