]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add a function to convert a Teuchos::Comm into an MPI_Comm. 16282/head
authorSebastian Kinnewig <kinnewig@ifam.uni-hannover.de>
Mon, 20 Nov 2023 16:22:38 +0000 (17:22 +0100)
committerSebastian Kinnewig <kinnewig@ifam.uni-hannover.de>
Mon, 20 Nov 2023 17:19:24 +0000 (18:19 +0100)
doc/news/changes/minor/20231120SebastianKinnewig [new file with mode: 0644]
include/deal.II/base/trilinos_utilities.h
include/deal.II/lac/trilinos_tpetra_vector.templates.h
source/base/trilinos_utilities.cc

diff --git a/doc/news/changes/minor/20231120SebastianKinnewig b/doc/news/changes/minor/20231120SebastianKinnewig
new file mode 100644 (file)
index 0000000..08c5101
--- /dev/null
@@ -0,0 +1,4 @@
+New: Introduce a new function, Utilities::Trilinos::teuchos_comm_to_mpi_comm(),
+to convert a Teuchos::RCP<Teuchos::Comm<int>> communicator into an MPI_Comm communicator.
+<br>
+(Sebastian Kinnewig, 2023/11/20)
index 6e6ffcb19a280524fcd4706b79bf077891e40df0..d9b29c22b67e7b7486c42582cd74bf841dede10e 100644 (file)
@@ -189,6 +189,16 @@ namespace Utilities
 #ifdef DEAL_II_TRILINOS_WITH_TPETRA
   namespace Trilinos
   {
+    /**
+     * Return the underlying MPI_Comm communicator from the
+     * <a
+     * href="https://docs.trilinos.org/dev/packages/teuchos/doc/html/classTeuchos_1_1Comm.html">Teuchos::Comm</a>
+     * communicator.
+     */
+    MPI_Comm
+    teuchos_comm_to_mpi_comm(
+      const Teuchos::RCP<const Teuchos::Comm<int>> &teuchos_comm);
+
     namespace internal
     {
       /**
index 06f32439f4a759196128925f78e38a54542d582c..bb122eb8908b2ee132dcf35c5f35de5413d8d73d 100644 (file)
@@ -235,12 +235,10 @@ namespace LinearAlgebra
                V.get_stored_elements().size()) ||
               (source_stored_elements != V.get_stored_elements()))
             {
-              const Teuchos::MpiComm<int> *mpi_comm =
-                dynamic_cast<const Teuchos::MpiComm<int> *>(
-                  vector->getMap()->getComm().get());
-              Assert(mpi_comm != nullptr, ExcInternalError());
-              create_tpetra_comm_pattern(V.get_stored_elements(),
-                                         *(mpi_comm->getRawMpiComm())());
+              create_tpetra_comm_pattern(
+                V.get_stored_elements(),
+                Utilities::Trilinos::teuchos_comm_to_mpi_comm(
+                  vector->getMap()->getComm()));
             }
         }
       else
@@ -526,12 +524,10 @@ namespace LinearAlgebra
         }
 
       // Check that the vector is zero on _all_ processors.
-      const Teuchos::MpiComm<int> *mpi_comm =
-        dynamic_cast<const Teuchos::MpiComm<int> *>(
-          vector->getMap()->getComm().get());
-      Assert(mpi_comm != nullptr, ExcInternalError());
       unsigned int num_nonzero =
-        Utilities::MPI::sum(flag, *(mpi_comm->getRawMpiComm())());
+        Utilities::MPI::sum(flag,
+                            Utilities::Trilinos::teuchos_comm_to_mpi_comm(
+                              vector->getMap()->getComm()));
 
       return num_nonzero == 0;
     }
@@ -609,11 +605,8 @@ namespace LinearAlgebra
     MPI_Comm
     Vector<Number>::get_mpi_communicator() const
     {
-      const auto *const tpetra_comm =
-        dynamic_cast<const Teuchos::MpiComm<int> *>(
-          vector->getMap()->getComm().get());
-      Assert(tpetra_comm != nullptr, ExcInternalError());
-      return *(tpetra_comm->getRawMpiComm())();
+      return Utilities::Trilinos::teuchos_comm_to_mpi_comm(
+        vector->getMap()->getComm());
     }
 
 
@@ -739,16 +732,8 @@ namespace LinearAlgebra
     MPI_Comm
     Vector<Number>::mpi_comm() const
     {
-      MPI_Comm out;
-#  ifdef DEAL_II_WITH_MPI
-      const Teuchos::MpiComm<int> *mpi_comm =
-        dynamic_cast<const Teuchos::MpiComm<int> *>(
-          vector->getMap()->getComm().get());
-      out = *mpi_comm->getRawMpiComm();
-#  else
-      out            = MPI_COMM_SELF;
-#  endif
-      return out;
+      return Utilities::Trilinos::teuchos_comm_to_mpi_comm(
+        vector->getMap()->getComm());
     }
 
 
index e9f2b105e5a3bb32f4b6f0de92a390b54ccfbe8b..1c275690292c376f87cc6e65ed4bba115ab2dd77 100644 (file)
@@ -169,6 +169,32 @@ namespace Utilities
     }
   } // namespace Trilinos
 #endif
+
+
+#ifdef DEAL_II_TRILINOS_WITH_TPETRA
+  namespace Trilinos
+  {
+    MPI_Comm
+    teuchos_comm_to_mpi_comm(
+      const Teuchos::RCP<const Teuchos::Comm<int>> &teuchos_comm)
+    {
+      MPI_Comm out;
+#  ifdef DEAL_II_WITH_MPI
+      // Cast from Teuchos::Comm<int> to Teuchos::MpiComm<int>.
+      const Teuchos::MpiComm<int> *mpi_comm =
+        dynamic_cast<const Teuchos::MpiComm<int> *>(teuchos_comm.get());
+      Assert(mpi_comm != nullptr, ExcInternalError());
+      // From the Teuchos::MpiComm<int> object we can extract
+      // the MPI_Comm object via the getRawMpiComm() function.
+      out = *(mpi_comm->getRawMpiComm())();
+#  else
+      out = MPI_COMM_SELF;
+#  endif
+      return out;
+    }
+  }    // namespace Trilinos
+#endif // DEAL_II_TRILINOS_WITH_TPETRA
+
 } // namespace Utilities
 
 DEAL_II_NAMESPACE_CLOSE

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.