]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix Tpetra detection and unit tests
authorDaniel Arndt <arndtd@ornl.gov>
Fri, 17 Feb 2023 19:10:54 +0000 (14:10 -0500)
committerDaniel Arndt <arndtd@ornl.gov>
Mon, 20 Feb 2023 12:58:08 +0000 (07:58 -0500)
cmake/configure/configure_20_trilinos.cmake
include/deal.II/base/index_set.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 593e02647138a6ef2950638920a9a7036937bc84..3655c55c9225ea62cd63b7c95627c3ec97c2a5cb 100644 (file)
@@ -233,6 +233,7 @@ macro(feature_trilinos_find_external var)
       list(APPEND CMAKE_REQUIRED_INCLUDES ${MPI_CXX_INCLUDE_PATH})
 
       list(APPEND CMAKE_REQUIRED_LIBRARIES ${Trilinos_LIBRARIES} ${MPI_LIBRARIES})
+      list(APPEND CMAKE_REQUIRED_FLAGS ${TRILINOS_CXX_FLAGS})
 
       # For the case of Trilinos being compiled with openmp support the
       # following Tpetra test needs -fopenmp to succeed. Make sure that we
@@ -249,11 +250,10 @@ macro(feature_trilinos_find_external var)
         "
         #include <cstdint>
         #include <Tpetra_Vector.hpp>
-        int
-        main()
+        int main()
         {
           using LO       = int;
-          using GO       = ${_global_index_type};
+          using GO       = std::make_signed_t<${_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);
@@ -283,6 +283,7 @@ macro(feature_trilinos_find_external var)
       list(APPEND CMAKE_REQUIRED_INCLUDES ${MPI_CXX_INCLUDE_PATH})
 
       list(APPEND CMAKE_REQUIRED_LIBRARIES ${Trilinos_LIBRARIES} ${MPI_LIBRARIES})
+      list(APPEND CMAKE_REQUIRED_FLAGS ${TRILINOS_CXX_FLAGS})
 
       CHECK_CXX_SOURCE_COMPILES(
         "
index b17c18282ce191948743af8cd1e35294ca449d90..722a3010ad1f2a616b8e123d250e3b8600dc116a 100644 (file)
@@ -496,7 +496,7 @@ public:
                     const bool      overlapping  = false) const;
 
 #  ifdef DEAL_II_TRILINOS_WITH_TPETRA
-  Tpetra::Map<int, types::global_dof_index>
+  Tpetra::Map<int, std::make_signed_t<types::global_dof_index>>
   make_tpetra_map(const MPI_Comm &communicator = MPI_COMM_WORLD,
                   const bool      overlapping  = false) const;
 #  endif
index 14b6f445ab70e642af975606e2915e81b8d57f1b..6236f67e0c9ae45ba065dc9ed0ca51230573d05a 100644 (file)
@@ -657,13 +657,14 @@ namespace LinearAlgebra
      * used directly.
      */
     void
-    import(
-      const Tpetra::Vector<Number, int, types::global_dof_index> &tpetra_vector,
-      const IndexSet &        locally_owned_elements,
-      VectorOperation::values operation,
-      const MPI_Comm &        mpi_comm,
-      const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
-        &communication_pattern);
+    import(const Tpetra::
+             Vector<Number, int, std::make_signed_t<types::global_dof_index>>
+               &                   tpetra_vector,
+           const IndexSet &        locally_owned_elements,
+           VectorOperation::values operation,
+           const MPI_Comm &        mpi_comm,
+           const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
+             &communication_pattern);
 #  endif
 
     /**
index c97afd5e3dab5046575616afc71dd7b374ac0322..a85caa9a18ffc8da84c232c20d9d27d39d648286 100644 (file)
@@ -577,10 +577,11 @@ namespace LinearAlgebra
   template <typename Number>
   void
   ReadWriteVector<Number>::import(
-    const Tpetra::Vector<Number, int, types::global_dof_index> &vector,
-    const IndexSet &                                            source_elements,
-    VectorOperation::values                                     operation,
-    const MPI_Comm &                                            mpi_comm,
+    const Tpetra::
+      Vector<Number, int, std::make_signed_t<types::global_dof_index>> &vector,
+    const IndexSet &        source_elements,
+    VectorOperation::values operation,
+    const MPI_Comm &        mpi_comm,
     const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
       &communication_pattern)
   {
@@ -620,11 +621,11 @@ namespace LinearAlgebra
                       "LinearAlgebra::TpetraWrappers::CommunicationPattern."));
       }
 
-    Tpetra::Export<int, types::global_dof_index> tpetra_export(
-      tpetra_comm_pattern->get_tpetra_export());
+    Tpetra::Export<int, std::make_signed_t<types::global_dof_index>>
+      tpetra_export(tpetra_comm_pattern->get_tpetra_export());
 
-    Tpetra::Vector<Number, int, types::global_dof_index> target_vector(
-      tpetra_export.getSourceMap());
+    Tpetra::Vector<Number, int, std::make_signed_t<types::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 a746580583f34eed84d83d7488f737baa59717df..862dadc24d05c56521f9760e3ad384f4dfa549c3 100644 (file)
@@ -68,13 +68,13 @@ namespace LinearAlgebra
       /**
        * Return the underlying Tpetra::Import object.
        */
-      const Tpetra::Import<int, types::global_dof_index> &
+      const Tpetra::Import<int, std::make_signed_t<types::global_dof_index>> &
       get_tpetra_import() const;
 
       /**
        * Return the underlying Tpetra::Export object.
        */
-      const Tpetra::Export<int, types::global_dof_index> &
+      const Tpetra::Export<int, std::make_signed_t<types::global_dof_index>> &
       get_tpetra_export() const;
 
     private:
@@ -86,13 +86,15 @@ namespace LinearAlgebra
       /**
        * Shared pointer to the Tpetra::Import object used.
        */
-      std::unique_ptr<Tpetra::Import<int, types::global_dof_index>>
+      std::unique_ptr<
+        Tpetra::Import<int, std::make_signed_t<types::global_dof_index>>>
         tpetra_import;
 
       /**
        * Shared pointer to the Tpetra::Export object used.
        */
-      std::unique_ptr<Tpetra::Export<int, types::global_dof_index>>
+      std::unique_ptr<
+        Tpetra::Export<int, std::make_signed_t<types::global_dof_index>>>
         tpetra_export;
     };
   } // end of namespace TpetraWrappers
index 26c1f1675a1c072342a45b72977cf1db3513c52c..6e3bdab2f735542f89099d7f90dfb0adab30c6fd 100644 (file)
@@ -337,14 +337,15 @@ namespace LinearAlgebra
        * Return a const reference to the underlying Trilinos
        * Tpetra::Vector class.
        */
-      const Tpetra::Vector<Number, int, types::global_dof_index> &
-      trilinos_vector() const;
+      const Tpetra::
+        Vector<Number, int, std::make_signed_t<types::global_dof_index>> &
+        trilinos_vector() const;
 
       /**
        * Return a (modifiable) reference to the underlying Trilinos
        * Tpetra::Vector class.
        */
-      Tpetra::Vector<Number, int, types::global_dof_index> &
+      Tpetra::Vector<Number, int, std::make_signed_t<types::global_dof_index>> &
       trilinos_vector();
 
       /**
@@ -398,7 +399,9 @@ namespace LinearAlgebra
       /**
        * Pointer to the actual Tpetra vector object.
        */
-      std::unique_ptr<Tpetra::Vector<Number, int, types::global_dof_index>>
+      std::unique_ptr<
+        Tpetra::
+          Vector<Number, int, std::make_signed_t<types::global_dof_index>>>
         vector;
 
       /**
index 154e45ddf0b3784b3883dad1d5a36a8fbb915975..54db805e3866275a66923d4a2c48e072fd1d7df0 100644 (file)
@@ -47,9 +47,12 @@ namespace LinearAlgebra
     template <typename Number>
     Vector<Number>::Vector()
       : Subscriptor()
-      , vector(new Tpetra::Vector<Number, int, types::global_dof_index>(
-          Teuchos::RCP<Tpetra::Map<int, types::global_dof_index>>(
-            new Tpetra::Map<int, types::global_dof_index>(
+      , 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>>(
               0,
               0,
               Utilities::Trilinos::tpetra_comm_self()))))
@@ -60,7 +63,9 @@ namespace LinearAlgebra
     template <typename Number>
     Vector<Number>::Vector(const Vector<Number> &V)
       : Subscriptor()
-      , vector(new Tpetra::Vector<Number, int, types::global_dof_index>(
+      , vector(new Tpetra::Vector<Number,
+                                  int,
+                                  std::make_signed_t<types::global_dof_index>>(
           V.trilinos_vector(),
           Teuchos::Copy))
     {}
@@ -71,9 +76,12 @@ namespace LinearAlgebra
     Vector<Number>::Vector(const IndexSet &parallel_partitioner,
                            const MPI_Comm &communicator)
       : Subscriptor()
-      , vector(new Tpetra::Vector<Number, int, types::global_dof_index>(
-          Teuchos::rcp(new Tpetra::Map<int, types::global_dof_index>(
-            parallel_partitioner.make_tpetra_map(communicator, false)))))
+      , 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)))))
     {}
 
 
@@ -84,12 +92,15 @@ namespace LinearAlgebra
                            const MPI_Comm &communicator,
                            const bool      omit_zeroing_entries)
     {
-      Tpetra::Map<int, types::global_dof_index> input_map =
+      Tpetra::Map<int, std::make_signed_t<types::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, types::global_dof_index>>(Teuchos::rcp(
-          new Tpetra::Map<int, types::global_dof_index>(input_map)));
+          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>>(
+              input_map)));
       else if (omit_zeroing_entries == false)
         {
           vector->putScalar(0.);
@@ -131,8 +142,8 @@ namespace LinearAlgebra
         {
           if (size() == V.size())
             {
-              Tpetra::Import<int, types::global_dof_index> data_exchange(
-                vector->getMap(), V.trilinos_vector().getMap());
+              Tpetra::Import<int, std::make_signed_t<types::global_dof_index>>
+                data_exchange(vector->getMap(), V.trilinos_vector().getMap());
 
               vector->doImport(V.trilinos_vector(),
                                data_exchange,
@@ -140,7 +151,9 @@ namespace LinearAlgebra
             }
           else
             vector = std::make_unique<
-              Tpetra::Vector<Number, int, types::global_dof_index>>(
+              Tpetra::Vector<Number,
+                             int,
+                             std::make_signed_t<types::global_dof_index>>>(
               V.trilinos_vector());
         }
 
@@ -201,22 +214,25 @@ namespace LinearAlgebra
               "LinearAlgebra::TpetraWrappers::CommunicationPattern."));
         }
 
-      Tpetra::Export<int, types::global_dof_index> tpetra_export(
-        tpetra_comm_pattern->get_tpetra_export());
-      Tpetra::Vector<Number, int, types::global_dof_index> source_vector(
-        tpetra_export.getSourceMap());
-
-      source_vector.template sync<Kokkos::HostSpace>();
-      auto x_2d = source_vector.template getLocalView<Kokkos::HostSpace>();
-      auto x_1d = Kokkos::subview(x_2d, Kokkos::ALL(), 0);
-      source_vector.template modify<Kokkos::HostSpace>();
-      const size_t localLength = source_vector.getLocalLength();
-      auto         values_it   = V.begin();
-      for (size_t k = 0; k < localLength; ++k)
-        x_1d(k) = *values_it++;
-      source_vector.template sync<
-        typename Tpetra::Vector<Number, int, types::global_dof_index>::
-          device_type::memory_space>();
+      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());
+
+      {
+        source_vector.template sync<Kokkos::HostSpace>();
+        auto x_2d = source_vector.template getLocalView<Kokkos::HostSpace>();
+        auto x_1d = Kokkos::subview(x_2d, Kokkos::ALL(), 0);
+        source_vector.template modify<Kokkos::HostSpace>();
+        const size_t localLength = source_vector.getLocalLength();
+        auto         values_it   = V.begin();
+        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>();
+      }
       if (operation == VectorOperation::insert)
         vector->doExport(source_vector, tpetra_export, Tpetra::REPLACE);
       else if (operation == VectorOperation::add)
@@ -274,10 +290,11 @@ 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, types::global_dof_index> dummy(
-            vector->getMap(), false);
-          Tpetra::Import<int, types::global_dof_index> data_exchange(
-            down_V.trilinos_vector().getMap(), dummy.getMap());
+          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());
 
           dummy.doImport(down_V.trilinos_vector(),
                          data_exchange,
@@ -338,8 +355,9 @@ namespace LinearAlgebra
           vector_1d(k) += a;
         }
       vector->template sync<
-        typename Tpetra::Vector<Number, int, types::global_dof_index>::
-          device_type::memory_space>();
+        typename Tpetra::
+          Vector<Number, int, std::make_signed_t<types::global_dof_index>>::
+            device_type::memory_space>();
     }
 
 
@@ -603,8 +621,9 @@ namespace LinearAlgebra
 
 
     template <typename Number>
-    const Tpetra::Vector<Number, int, types::global_dof_index> &
-    Vector<Number>::trilinos_vector() const
+    const Tpetra::
+      Vector<Number, int, std::make_signed_t<types::global_dof_index>> &
+      Vector<Number>::trilinos_vector() const
     {
       return *vector;
     }
@@ -612,7 +631,7 @@ namespace LinearAlgebra
 
 
     template <typename Number>
-    Tpetra::Vector<Number, int, types::global_dof_index> &
+    Tpetra::Vector<Number, int, std::make_signed_t<types::global_dof_index>> &
     Vector<Number>::trilinos_vector()
     {
       return *vector;
index 5893aba997e1d8ca741345beebfe9f3cc8082518..2d2307769a942b02f8cda50082feb7c96a9b16cc 100644 (file)
@@ -157,8 +157,8 @@ namespace internal
     LinearAlgebra::TpetraWrappers::Vector<NumberType> &V)
   {
     // Extract local indices in the vector.
-    Tpetra::Vector<NumberType, int, types::global_dof_index> vector =
-      V.trilinos_vector();
+    Tpetra::Vector<NumberType, int, std::make_signed_t<types::global_dof_index>>
+                                      vector = V.trilinos_vector();
     TrilinosWrappers::types::int_type trilinos_i =
       vector.getMap()->getLocalElement(
         static_cast<TrilinosWrappers::types::int_type>(i));
@@ -170,8 +170,9 @@ namespace internal
     vector.template modify<Kokkos::HostSpace>();
     vector_1d(trilinos_i) += value;
     vector.template sync<
-      typename Tpetra::Vector<NumberType, int, types::global_dof_index>::
-        device_type::memory_space>();
+      typename Tpetra::
+        Vector<NumberType, int, std::make_signed_t<types::global_dof_index>>::
+          device_type::memory_space>();
   }
 
 
@@ -184,8 +185,8 @@ namespace internal
     LinearAlgebra::TpetraWrappers::Vector<NumberType> &V)
   {
     // Extract local indices in the vector.
-    Tpetra::Vector<NumberType, int, types::global_dof_index> vector =
-      V.trilinos_vector();
+    Tpetra::Vector<NumberType, int, std::make_signed_t<types::global_dof_index>>
+                                      vector = V.trilinos_vector();
     TrilinosWrappers::types::int_type trilinos_i =
       vector.getMap()->getLocalElement(
         static_cast<TrilinosWrappers::types::int_type>(i));
@@ -197,8 +198,9 @@ namespace internal
     vector.template modify<Kokkos::HostSpace>();
     vector_1d(trilinos_i) = value;
     vector.template sync<
-      typename Tpetra::Vector<NumberType, int, types::global_dof_index>::
-        device_type::memory_space>();
+      typename Tpetra::
+        Vector<NumberType, int, std::make_signed_t<types::global_dof_index>>::
+          device_type::memory_space>();
   }
 
 
@@ -210,8 +212,8 @@ namespace internal
     const types::global_dof_index                            i)
   {
     // Extract local indices in the vector.
-    Tpetra::Vector<NumberType, int, types::global_dof_index> vector =
-      V.trilinos_vector();
+    Tpetra::Vector<NumberType, int, std::make_signed_t<types::global_dof_index>>
+                                      vector = V.trilinos_vector();
     TrilinosWrappers::types::int_type trilinos_i =
       vector.getMap()->getLocalElement(
         static_cast<TrilinosWrappers::types::int_type>(i));
index bc68b5fefcdf9f25fa9fd8e2b0c8deb74dd3ad98..81e38f538f8e39843881cb8bc83b0e48b0de1e8a 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, types::global_dof_index>
+Tpetra::Map<int, std::make_signed_t<types::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, types::global_dof_index>(
+    return Tpetra::Map<int, std::make_signed_t<types::global_dof_index>>(
       size(),
       n_elements(),
       0,
@@ -764,11 +764,13 @@ IndexSet::make_tpetra_map(const MPI_Comm &communicator,
     );
   else
     {
-      const std::vector<size_type>         indices = get_index_vector();
-      std::vector<types::global_dof_index> int_indices(indices.size());
+      const std::vector<size_type> indices = get_index_vector();
+      std::vector<std::make_signed_t<types::global_dof_index>> int_indices(
+        indices.size());
       std::copy(indices.begin(), indices.end(), int_indices.begin());
-      const Teuchos::ArrayView<types::global_dof_index> arr_view(int_indices);
-      return Tpetra::Map<int, types::global_dof_index>(
+      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>>(
         size(),
         arr_view,
         0,
index afddac8413ffc6a67e69a8d213d2718f7535f711..34c08d633279147c7759eed6978f0d2e5a4d1b3a 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, types::global_dof_index>(
+      auto vector_space_vector_map = Teuchos::rcp(
+        new Tpetra::Map<int, std::make_signed_t<types::global_dof_index>>(
           vector_space_vector_index_set.make_tpetra_map(*comm, false)));
-      auto read_write_vector_map =
-        Teuchos::rcp(new Tpetra::Map<int, types::global_dof_index>(
+      auto read_write_vector_map = Teuchos::rcp(
+        new Tpetra::Map<int, std::make_signed_t<types::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, types::global_dof_index>>(
-          read_write_vector_map, vector_space_vector_map);
-      tpetra_export =
-        std::make_unique<Tpetra::Export<int, types::global_dof_index>>(
-          read_write_vector_map, vector_space_vector_map);
+      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);
     }
 
 
@@ -79,7 +79,7 @@ namespace LinearAlgebra
 
 
 
-    const Tpetra::Import<int, types::global_dof_index> &
+    const Tpetra::Import<int, std::make_signed_t<types::global_dof_index>> &
     CommunicationPattern::get_tpetra_import() const
     {
       return *tpetra_import;
@@ -87,7 +87,7 @@ namespace LinearAlgebra
 
 
 
-    const Tpetra::Export<int, types::global_dof_index> &
+    const Tpetra::Export<int, std::make_signed_t<types::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.