From: Sebastian Kinnewig Date: Mon, 20 Nov 2023 16:22:38 +0000 (+0100) Subject: Add a function to convert a Teuchos::Comm into an MPI_Comm. X-Git-Tag: relicensing~288^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=63846902b8da739320e09166501612f5965bc3fe;p=dealii.git Add a function to convert a Teuchos::Comm into an MPI_Comm. --- diff --git a/doc/news/changes/minor/20231120SebastianKinnewig b/doc/news/changes/minor/20231120SebastianKinnewig new file mode 100644 index 0000000000..08c510180d --- /dev/null +++ b/doc/news/changes/minor/20231120SebastianKinnewig @@ -0,0 +1,4 @@ +New: Introduce a new function, Utilities::Trilinos::teuchos_comm_to_mpi_comm(), +to convert a Teuchos::RCP> communicator into an MPI_Comm communicator. +
+(Sebastian Kinnewig, 2023/11/20) diff --git a/include/deal.II/base/trilinos_utilities.h b/include/deal.II/base/trilinos_utilities.h index 6e6ffcb19a..d9b29c22b6 100644 --- a/include/deal.II/base/trilinos_utilities.h +++ b/include/deal.II/base/trilinos_utilities.h @@ -189,6 +189,16 @@ namespace Utilities #ifdef DEAL_II_TRILINOS_WITH_TPETRA namespace Trilinos { + /** + * Return the underlying MPI_Comm communicator from the + * Teuchos::Comm + * communicator. + */ + MPI_Comm + teuchos_comm_to_mpi_comm( + const Teuchos::RCP> &teuchos_comm); + namespace internal { /** diff --git a/include/deal.II/lac/trilinos_tpetra_vector.templates.h b/include/deal.II/lac/trilinos_tpetra_vector.templates.h index 06f32439f4..bb122eb890 100644 --- a/include/deal.II/lac/trilinos_tpetra_vector.templates.h +++ b/include/deal.II/lac/trilinos_tpetra_vector.templates.h @@ -235,12 +235,10 @@ namespace LinearAlgebra V.get_stored_elements().size()) || (source_stored_elements != V.get_stored_elements())) { - const Teuchos::MpiComm *mpi_comm = - dynamic_cast *>( - 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 *mpi_comm = - dynamic_cast *>( - 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::get_mpi_communicator() const { - const auto *const tpetra_comm = - dynamic_cast *>( - 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::mpi_comm() const { - MPI_Comm out; -# ifdef DEAL_II_WITH_MPI - const Teuchos::MpiComm *mpi_comm = - dynamic_cast *>( - 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()); } diff --git a/source/base/trilinos_utilities.cc b/source/base/trilinos_utilities.cc index e9f2b105e5..1c27569029 100644 --- a/source/base/trilinos_utilities.cc +++ b/source/base/trilinos_utilities.cc @@ -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> &teuchos_comm) + { + MPI_Comm out; +# ifdef DEAL_II_WITH_MPI + // Cast from Teuchos::Comm to Teuchos::MpiComm. + const Teuchos::MpiComm *mpi_comm = + dynamic_cast *>(teuchos_comm.get()); + Assert(mpi_comm != nullptr, ExcInternalError()); + // From the Teuchos::MpiComm 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