From 40c975893b96c3a538895bcfda935e6431afc7b5 Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Tue, 18 Jan 2022 08:44:11 -0700 Subject: [PATCH] Do not store some things in the CA::Interface base class. The base class currently stores a number of things about the communicator that are either unused in any of the implementations, or only used in one place in the Selector::run() function. Remove these member variables and just query what we need in the one place where it's used. --- include/deal.II/base/mpi_consensus_algorithms.h | 15 --------------- .../base/mpi_consensus_algorithms.templates.h | 13 +++++++------ 2 files changed, 7 insertions(+), 21 deletions(-) diff --git a/include/deal.II/base/mpi_consensus_algorithms.h b/include/deal.II/base/mpi_consensus_algorithms.h index 3c34b591b4..4319fe0b75 100644 --- a/include/deal.II/base/mpi_consensus_algorithms.h +++ b/include/deal.II/base/mpi_consensus_algorithms.h @@ -244,21 +244,6 @@ namespace Utilities * MPI communicator. */ const MPI_Comm &comm; - - /** - * Cache if job supports MPI. - */ - const bool job_supports_mpi; - - /** - * Rank of this process. - */ - const unsigned int my_rank; - - /** - * Number of processes in the communicator. - */ - const unsigned int n_procs; }; diff --git a/include/deal.II/base/mpi_consensus_algorithms.templates.h b/include/deal.II/base/mpi_consensus_algorithms.templates.h index b0c2993bed..e2c7090a2c 100644 --- a/include/deal.II/base/mpi_consensus_algorithms.templates.h +++ b/include/deal.II/base/mpi_consensus_algorithms.templates.h @@ -116,9 +116,6 @@ namespace Utilities const MPI_Comm & comm) : process(process) , comm(comm) - , job_supports_mpi(Utilities::MPI::job_supports_mpi()) - , my_rank(job_supports_mpi ? this_mpi_process(comm) : 0) - , n_procs(job_supports_mpi ? n_mpi_processes(comm) : 1) {} @@ -692,18 +689,22 @@ namespace Utilities // able to test also the non-blocking implementation. This feature // is tested by: // tests/multigrid/transfer_matrix_free_06.with_mpi=true.with_p4est=true.with_trilinos=true.mpirun=10.output + + const unsigned int n_procs = (Utilities::MPI::job_supports_mpi() ? + Utilities::MPI::n_mpi_processes(comm) : + 1); #ifdef DEAL_II_WITH_MPI # if DEAL_II_MPI_VERSION_GTE(3, 0) # ifdef DEBUG - if (this->n_procs > 10) + if (n_procs > 10) # else - if (this->n_procs > 99) + if (n_procs > 99) # endif consensus_algo.reset(new NBX(process, comm)); else # endif #endif - if (this->n_procs > 1) + if (n_procs > 1) consensus_algo.reset(new PEX(process, comm)); else consensus_algo.reset(new Serial(process, comm)); -- 2.39.5