class Interface
{
public:
+ /**
+ * Default constructor.
+ */
+ Interface();
+
/**
* Constructor. @p process is an object that provides information
* about what processes the current process wants to communicate with,
* and the data to be sent/received. @p comm is the communicator on
* which this communication is to happen.
+ *
+ * @deprecated This constructor stores the Process object and the
+ * communicator so that one can later call the run() function
+ * without arguments. This approach is deprecated. Instead, use
+ * the default constructor of this class along with the run()
+ * function that takes an argument.
*/
+ DEAL_II_DEPRECATED
Interface(Process<T1, T2> &process, const MPI_Comm &comm);
/**
virtual ~Interface() = default;
/**
- * Run the consensus algorithm and return the requesting processes.
+ * Run the consensus algorithm and return a vector of process ranks
+ * that have requested answers from the current process.
+ *
+ * @deprecated This function is deprecated. It can be called
+ * if the Process object and communicator to be used have previously
+ * been provided to the non-default constructor. Use the run()
+ * functions taking arguments instead.
+ */
+ DEAL_II_DEPRECATED
+ std::vector<unsigned int>
+ run();
+
+ /**
+ * Run the consensus algorithm and return a vector of process ranks
+ * that have requested answers from the current process.
+ *
+ * This version of the run() function simply unpacks the functions
+ * packaged in `process` and calls the version of the run() function
+ * that takes a number of `std::function` arguments.
+ */
+ std::vector<unsigned int>
+ run(Process<T1, T2> &process, const MPI_Comm &comm);
+
+ /**
+ * Run the consensus algorithm and return a vector of process ranks
+ * that have requested answers from the current process.
+ *
+ * @param[in] targets A vector that contains the ranks of processes
+ * to which requests should be sent and from which answers need
+ * to be received.
+ * @param[in] create_request A function object that takes the rank
+ * of a target process as argument and returns the message that
+ * forms the request to this target.
+ * @param[in] answer_request A function that takes as arguments the
+ * rank of the process that has sent a request to us, along with
+ * the message of the request, and returns the message that forms
+ * the answer that should be sent back to the requesting process.
+ * @param[in] process_answer A function object that takes as argument
+ * the rank of a process from which we have received an answer
+ * to a previously sent request, along with the message that
+ * forms this answer. This function is used to describe what
+ * the caller of the consensus algorithm wants to do with the
+ * received answer.
+ * @param[in] comm The MPI communicator on which the whole algorithm
+ * is to be performed.
*/
virtual std::vector<unsigned int>
- run() = 0;
+ run(const std::vector<unsigned int> &targets,
+ const std::function<std::vector<T1>(const unsigned int)>
+ &create_request,
+ const std::function<std::vector<T2>(const unsigned int,
+ const std::vector<T1> &)>
+ & answer_request,
+ const std::function<void(const unsigned int,
+ const std::vector<T2> &)> &process_answer,
+ const MPI_Comm & comm) = 0;
- protected:
+ private:
/**
* Reference to the process provided by the user.
+ *
+ * This member variable is only used in the deprecated constructor
+ * and the run() function without argument. It is a `nullptr`
+ * otherwise
*/
- Process<T1, T2> &process;
+ DEAL_II_DEPRECATED
+ Process<T1, T2> *process;
/**
* MPI communicator.
+ *
+ * This member variable is only used in the deprecated constructor
+ * and the run() function without argument.
*/
- const MPI_Comm &comm;
+ DEAL_II_DEPRECATED
+ MPI_Comm comm;
};
class NBX : public Interface<T1, T2>
{
public:
+ /**
+ * Default constructor.
+ */
+ NBX() = default;
+
/**
* Constructor.
*
* @param process Process to be run during consensus algorithm.
* @param comm MPI Communicator
+ *
+ * @deprecated This constructor stores the Process object and the
+ * communicator so that one can later call the run() function
+ * without arguments. This approach is deprecated. Instead, use
+ * the default constructor of this class along with the run()
+ * function that takes an argument.
*/
+ DEAL_II_DEPRECATED
NBX(Process<T1, T2> &process, const MPI_Comm &comm);
/**
*/
virtual ~NBX() = default;
+ // Import the declarations from the base class.
+ using Interface<T1, T2>::run;
+
/**
* @copydoc Interface::run()
*/
virtual std::vector<unsigned int>
- run() override;
+ run(const std::vector<unsigned int> &targets,
+ const std::function<std::vector<T1>(const unsigned int)>
+ &create_request,
+ const std::function<std::vector<T2>(const unsigned int,
+ const std::vector<T1> &)>
+ & answer_request,
+ const std::function<void(const unsigned int,
+ const std::vector<T2> &)> &process_answer,
+ const MPI_Comm & comm) override;
private:
#ifdef DEAL_II_WITH_MPI
- /**
- * List of processes this process wants to send requests to.
- */
- std::vector<unsigned int> targets;
-
/**
* Buffers for sending requests.
*/
* have been satisfied.
*/
bool
- all_locally_originated_receives_are_completed();
+ all_locally_originated_receives_are_completed(
+ const std::function<void(const unsigned int, const std::vector<T2> &)>
+ & process_answer,
+ const MPI_Comm &comm);
/**
* Signal to all other ranks that this rank has received all request
* answers via entering IBarrier.
*/
void
- signal_finish();
+ signal_finish(const MPI_Comm &comm);
/**
* Check whether all of the requests for answers that were created by
* answer.
*/
void
- maybe_answer_one_request();
+ maybe_answer_one_request(
+ const std::function<std::vector<T2>(const unsigned int,
+ const std::vector<T1> &)>
+ & answer_request,
+ const MPI_Comm &comm);
/**
* Start to send all requests via ISend and post IRecvs for the incoming
* answer messages.
*/
void
- start_communication();
+ start_communication(
+ const std::vector<unsigned int> &targets,
+ const std::function<std::vector<T1>(const unsigned int)>
+ & create_request,
+ const MPI_Comm &comm);
/**
* After all rank has received all answers, the MPI data structures can
* be freed and the received answers can be processed.
*/
void
- clean_up_and_end_communication();
+ clean_up_and_end_communication(const MPI_Comm &comm);
};
/**
class PEX : public Interface<T1, T2>
{
public:
+ /**
+ * Default constructor.
+ */
+ PEX() = default;
+
+
/**
* Constructor.
*
* @param process Process to be run during consensus algorithm.
* @param comm MPI Communicator
+ *
+ * @deprecated This constructor stores the Process object and the
+ * communicator so that one can later call the run() function
+ * without arguments. This approach is deprecated. Instead, use
+ * the default constructor of this class along with the run()
+ * function that takes an argument.
*/
+ DEAL_II_DEPRECATED
PEX(Process<T1, T2> &process, const MPI_Comm &comm);
/**
*/
virtual ~PEX() = default;
+ // Import the declarations from the base class.
+ using Interface<T1, T2>::run;
+
/**
* @copydoc Interface::run()
*/
virtual std::vector<unsigned int>
- run() override;
+ run(const std::vector<unsigned int> &targets,
+ const std::function<std::vector<T1>(const unsigned int)>
+ &create_request,
+ const std::function<std::vector<T2>(const unsigned int,
+ const std::vector<T1> &)>
+ & answer_request,
+ const std::function<void(const unsigned int,
+ const std::vector<T2> &)> &process_answer,
+ const MPI_Comm & comm) override;
private:
#ifdef DEAL_II_WITH_MPI
- /**
- * List of ranks of processes this processes wants to send a request to.
- */
- std::vector<unsigned int> targets;
-
- // data structures to send and receive requests
-
/**
* Buffers for sending requests.
*/
* answer messages.
*/
unsigned int
- start_communication();
+ start_communication(
+ const std::vector<unsigned int> &targets,
+ const std::function<std::vector<T1>(const unsigned int)>
+ & create_request,
+ const MPI_Comm &comm);
/**
- * The ith request message from another rank has been received: process
- * the request and send an answer.
+ * The `index`th request message from another rank has been received:
+ * process the request and send an answer.
*/
void
- answer_one_request(const unsigned int index);
+ answer_one_request(
+ const unsigned int index,
+ const std::function<std::vector<T2>(const unsigned int,
+ const std::vector<T1> &)>
+ & answer_request,
+ const MPI_Comm &comm);
/**
* Receive and process all of the incoming responses to the
* requests we sent.
*/
void
- process_incoming_answers();
+ process_incoming_answers(
+ const unsigned int n_targets,
+ const std::function<void(const unsigned int, const std::vector<T2> &)>
+ & process_answer,
+ const MPI_Comm &comm);
/**
* After all answers have been exchanged, the MPI data structures can be
class Serial : public Interface<T1, T2>
{
public:
+ /**
+ * Default constructor.
+ */
+ Serial() = default;
+
/**
* Constructor.
*
* @param process Process to be run during consensus algorithm.
* @param comm MPI Communicator (ignored)
+ *
+ * @deprecated This constructor stores the Process object and the
+ * communicator so that one can later call the run() function
+ * without arguments. This approach is deprecated. Instead, use
+ * the default constructor of this class along with the run()
+ * function that takes an argument.
*/
+ DEAL_II_DEPRECATED
Serial(Process<T1, T2> &process, const MPI_Comm &comm);
+ // Import the declarations from the base class.
+ using Interface<T1, T2>::run;
+
/**
* @copydoc Interface::run()
*/
virtual std::vector<unsigned int>
- run() override;
+ run(const std::vector<unsigned int> &targets,
+ const std::function<std::vector<T1>(const unsigned int)>
+ &create_request,
+ const std::function<std::vector<T2>(const unsigned int,
+ const std::vector<T1> &)>
+ & answer_request,
+ const std::function<void(const unsigned int,
+ const std::vector<T2> &)> &process_answer,
+ const MPI_Comm & comm) override;
};
class Selector : public Interface<T1, T2>
{
public:
+ /**
+ * Default constructor.
+ */
+ Selector() = default;
+
/**
* Constructor.
*
* @param process Process to be run during consensus algorithm.
* @param comm MPI Communicator.
+ *
+ * @deprecated This constructor stores the Process object and the
+ * communicator so that one can later call the run() function
+ * without arguments. This approach is deprecated. Instead, use
+ * the default constructor of this class along with the run()
+ * function that takes an argument.
*/
+ DEAL_II_DEPRECATED
Selector(Process<T1, T2> &process, const MPI_Comm &comm);
/**
*/
virtual ~Selector() = default;
+ // Import the declarations from the base class.
+ using Interface<T1, T2>::run;
+
/**
* @copydoc Interface::run()
*
* @note The function call is delegated to another ConsensusAlgorithms::Interface implementation.
*/
virtual std::vector<unsigned int>
- run() override;
+ run(const std::vector<unsigned int> &targets,
+ const std::function<std::vector<T1>(const unsigned int)>
+ &create_request,
+ const std::function<std::vector<T2>(const unsigned int,
+ const std::vector<T1> &)>
+ & answer_request,
+ const std::function<void(const unsigned int,
+ const std::vector<T2> &)> &process_answer,
+ const MPI_Comm & comm) override;
private:
// Pointer to the actual ConsensusAlgorithms::Interface implementation.
template <typename T1, typename T2>
Interface<T1, T2>::Interface(Process<T1, T2> &process,
const MPI_Comm & comm)
- : process(process)
+ : process(&process)
, comm(comm)
{}
+ template <typename T1, typename T2>
+ Interface<T1, T2>::Interface()
+ : process(nullptr)
+ , comm(MPI_COMM_NULL)
+ {}
+
+
+
+ template <typename T1, typename T2>
+ std::vector<unsigned int>
+ Interface<T1, T2>::run()
+ {
+ Assert(process != nullptr,
+ ExcMessage("This function can only be called if the "
+ "deprecated non-default constructor of this class "
+ "has previously been called to set the Process "
+ "object and a communicator."));
+ return run(*process, comm);
+ }
+
+
+
+ template <typename T1, typename T2>
+ std::vector<unsigned int>
+ Interface<T1, T2>::run(Process<T1, T2> &process, const MPI_Comm &comm)
+ {
+ // Unpack the 'process' object and call the function that takes
+ // function objects for all operations.
+ return run(
+ process.compute_targets(),
+ /* create_request: */
+ [&process](const unsigned int target) {
+ std::vector<T1> request;
+ process.create_request(target, request);
+ return request;
+ },
+ /* answer_request: */
+ [&process](const unsigned int source,
+ const std::vector<T1> &request) {
+ std::vector<T2> answer;
+ process.answer_request(source, request, answer);
+ return answer;
+ },
+ /* process_answer: */
+ [&process](const unsigned int target, const std::vector<T2> &answer) {
+ process.read_answer(target, answer);
+ },
+ comm);
+ }
+
+
+
template <typename T1, typename T2>
NBX<T1, T2>::NBX(Process<T1, T2> &process, const MPI_Comm &comm)
: Interface<T1, T2>(process, comm)
template <typename T1, typename T2>
std::vector<unsigned int>
- NBX<T1, T2>::run()
+ NBX<T1, T2>::run(
+ const std::vector<unsigned int> &targets,
+ const std::function<std::vector<T1>(const unsigned int)>
+ &create_request,
+ const std::function<std::vector<T2>(const unsigned int,
+ const std::vector<T1> &)>
+ &answer_request,
+ const std::function<void(const unsigned int, const std::vector<T2> &)>
+ & process_answer,
+ const MPI_Comm &comm)
{
+ Assert(has_unique_elements(targets),
+ ExcMessage("The consensus algorithms expect that each process "
+ "only sends a single message to another process, "
+ "but the targets provided include duplicates."));
+
static CollectiveMutex mutex;
- CollectiveMutex::ScopedLock lock(mutex, this->comm);
+ CollectiveMutex::ScopedLock lock(mutex, comm);
// 1) Send data to identified targets and start receiving
// the answers from these very same processes.
- start_communication();
+ start_communication(targets, create_request, comm);
// 2) Until all posted receive operations are known to have completed,
// answer requests and keep checking whether all requests of
// That's ok: We will get around to dealing with all remaining
// message later. We just want to move on to the next step
// as early as possible.
- while (all_locally_originated_receives_are_completed() == false)
- maybe_answer_one_request();
+ while (all_locally_originated_receives_are_completed(process_answer,
+ comm) == false)
+ maybe_answer_one_request(answer_request, comm);
// 3) Signal to all other processes that all requests of this process
// have been answered
- signal_finish();
+ signal_finish(comm);
// 4) Nevertheless, this process has to keep on answering (potential)
// incoming requests until all processes have received the
// answer to all requests
while (all_remotely_originated_receives_are_completed() == false)
- maybe_answer_one_request();
+ maybe_answer_one_request(answer_request, comm);
// 5) process the answer to all requests
- clean_up_and_end_communication();
+ clean_up_and_end_communication(comm);
return std::vector<unsigned int>(requesting_processes.begin(),
requesting_processes.end());
template <typename T1, typename T2>
void
- NBX<T1, T2>::start_communication()
+ NBX<T1, T2>::start_communication(
+ const std::vector<unsigned int> &targets,
+ const std::function<std::vector<T1>(const unsigned int)>
+ & create_request,
+ const MPI_Comm &comm)
{
#ifdef DEAL_II_WITH_MPI
// 1)
- targets = this->process.compute_targets();
- Assert(has_unique_elements(targets),
- ExcMessage("The consensus algorithms expect that each process "
- "only sends a single message to another process, "
- "but the targets provided include duplicates."));
-
const auto n_targets = targets.size();
const int tag_request = Utilities::MPI::internal::Tags::
for (unsigned int index = 0; index < n_targets; ++index)
{
const unsigned int rank = targets[index];
- AssertIndexRange(rank,
- Utilities::MPI::n_mpi_processes(this->comm));
+ AssertIndexRange(rank, Utilities::MPI::n_mpi_processes(comm));
auto &send_buffer = send_buffers[index];
- this->process.create_request(rank, send_buffer);
+ send_buffer = create_request(rank);
// Post a request to send data
auto ierr = MPI_Isend(send_buffer.data(),
MPI_BYTE,
rank,
tag_request,
- this->comm,
+ comm,
&send_requests[index]);
AssertThrowMPI(ierr);
}
// a request to:
n_outstanding_answers = n_targets;
}
+#else
+ (void)targets;
+ (void)create_request;
+ (void)comm;
#endif
}
template <typename T1, typename T2>
bool
- NBX<T1, T2>::all_locally_originated_receives_are_completed()
+ NBX<T1, T2>::all_locally_originated_receives_are_completed(
+ const std::function<void(const unsigned int, const std::vector<T2> &)>
+ & process_answer,
+ const MPI_Comm &comm)
{
#ifdef DEAL_II_WITH_MPI
// We know that all requests have come in when we have pending
int request_is_pending;
MPI_Status status;
- const auto ierr = MPI_Iprobe(MPI_ANY_SOURCE,
- tag_deliver,
- this->comm,
- &request_is_pending,
- &status);
+ const auto ierr = MPI_Iprobe(
+ MPI_ANY_SOURCE, tag_deliver, comm, &request_is_pending, &status);
AssertThrowMPI(ierr);
// If there is no pending message with this tag,
MPI_BYTE,
target,
tag_deliver,
- this->comm,
+ comm,
MPI_STATUS_IGNORE);
AssertThrowMPI(ierr);
}
- this->process.read_answer(target, recv_buffer);
+ process_answer(target, recv_buffer);
// Finally, remove this rank from the list of outstanding
// targets:
}
#else
+ (void)process_answer;
+ (void)comm;
+
return true;
#endif
}
template <typename T1, typename T2>
void
- NBX<T1, T2>::maybe_answer_one_request()
+ NBX<T1, T2>::maybe_answer_one_request(
+ const std::function<std::vector<T2>(const unsigned int,
+ const std::vector<T1> &)>
+ & answer_request,
+ const MPI_Comm &comm)
{
#ifdef DEAL_II_WITH_MPI
// only answer one.
MPI_Status status;
int request_is_pending;
- const auto ierr = MPI_Iprobe(MPI_ANY_SOURCE,
- tag_request,
- this->comm,
- &request_is_pending,
- &status);
+ const auto ierr = MPI_Iprobe(
+ MPI_ANY_SOURCE, tag_request, comm, &request_is_pending, &status);
AssertThrowMPI(ierr);
if (request_is_pending != 0)
MPI_BYTE,
other_rank,
tag_request,
- this->comm,
+ comm,
MPI_STATUS_IGNORE);
AssertThrowMPI(ierr);
// and ask the 'process' object to produce an answer:
request_buffers.emplace_back(std::make_unique<std::vector<T2>>());
auto &request_buffer = *request_buffers.back();
- this->process.answer_request(other_rank,
- buffer_recv,
- request_buffer);
+ request_buffer = answer_request(other_rank, buffer_recv);
// Then initiate sending the answer back to the requester.
request_requests.emplace_back(std::make_unique<MPI_Request>());
MPI_BYTE,
other_rank,
tag_deliver,
- this->comm,
+ comm,
request_requests.back().get());
AssertThrowMPI(ierr);
}
+#else
+ (void)answer_request;
+ (void)comm;
#endif
}
template <typename T1, typename T2>
void
- NBX<T1, T2>::signal_finish()
+ NBX<T1, T2>::signal_finish(const MPI_Comm &comm)
{
#ifdef DEAL_II_WITH_MPI
# if DEAL_II_MPI_VERSION_GTE(3, 0)
- const auto ierr = MPI_Ibarrier(this->comm, &barrier_request);
+ const auto ierr = MPI_Ibarrier(comm, &barrier_request);
AssertThrowMPI(ierr);
# else
- AssertThrow(
- false,
- ExcMessage(
- "ConsensusAlgorithms::NBX uses MPI 3.0 features. You should compile with at least MPI 3.0."));
+ AssertThrow(false,
+ ExcMessage(
+ "ConsensusAlgorithms::NBX uses MPI 3.0 features. "
+ "You should compile with at least MPI 3.0."));
# endif
+#else
+ (void)comm;
#endif
}
template <typename T1, typename T2>
void
- NBX<T1, T2>::clean_up_and_end_communication()
+ NBX<T1, T2>::clean_up_and_end_communication(const MPI_Comm &comm)
{
#ifdef DEAL_II_WITH_MPI
// clean up
# ifdef DEBUG
// note: IBarrier seems to make problem during testing, this
// additional Barrier seems to help
- ierr = MPI_Barrier(this->comm);
+ ierr = MPI_Barrier(comm);
AssertThrowMPI(ierr);
# endif
}
+#else
+ (void)comm;
#endif
}
template <typename T1, typename T2>
std::vector<unsigned int>
- PEX<T1, T2>::run()
+ PEX<T1, T2>::run(
+ const std::vector<unsigned int> &targets,
+ const std::function<std::vector<T1>(const unsigned int)>
+ &create_request,
+ const std::function<std::vector<T2>(const unsigned int,
+ const std::vector<T1> &)>
+ &answer_request,
+ const std::function<void(const unsigned int, const std::vector<T2> &)>
+ & process_answer,
+ const MPI_Comm &comm)
{
+ Assert(has_unique_elements(targets),
+ ExcMessage("The consensus algorithms expect that each process "
+ "only sends a single message to another process, "
+ "but the targets provided include duplicates."));
+
static CollectiveMutex mutex;
- CollectiveMutex::ScopedLock lock(mutex, this->comm);
+ CollectiveMutex::ScopedLock lock(mutex, comm);
// 1) Send requests and start receiving the answers.
// In particular, determine how many requests we should expect
// on the current process.
- const unsigned int n_requests = start_communication();
+ const unsigned int n_requests =
+ start_communication(targets, create_request, comm);
// 2) Answer requests:
for (unsigned int request = 0; request < n_requests; ++request)
- answer_one_request(request);
+ answer_one_request(request, answer_request, comm);
// 3) Process answers:
- process_incoming_answers();
+ process_incoming_answers(targets.size(), process_answer, comm);
// 4) Make sure all sends have successfully terminated:
clean_up_and_end_communication();
template <typename T1, typename T2>
unsigned int
- PEX<T1, T2>::start_communication()
+ PEX<T1, T2>::start_communication(
+ const std::vector<unsigned int> &targets,
+ const std::function<std::vector<T1>(const unsigned int)>
+ & create_request,
+ const MPI_Comm &comm)
{
#ifdef DEAL_II_WITH_MPI
const int tag_request = Utilities::MPI::internal::Tags::
// 1) determine with which processes this process wants to communicate
// with
- targets = this->process.compute_targets();
- Assert(has_unique_elements(targets),
- ExcMessage("The consensus algorithms expect that each process "
- "only sends a single message to another process, "
- "but the targets provided include duplicates."));
const unsigned int n_targets = targets.size();
// 2) determine who wants to communicate with this process
const unsigned int n_sources =
- compute_n_point_to_point_communications(this->comm, targets);
+ compute_n_point_to_point_communications(comm, targets);
// 2) allocate memory
recv_buffers.resize(n_targets);
for (unsigned int i = 0; i < n_targets; ++i)
{
const unsigned int rank = targets[i];
- AssertIndexRange(rank, Utilities::MPI::n_mpi_processes(this->comm));
+ AssertIndexRange(rank, Utilities::MPI::n_mpi_processes(comm));
// pack data which should be sent
auto &send_buffer = send_buffers[i];
- this->process.create_request(rank, send_buffer);
+ send_buffer = create_request(rank);
// start to send data
auto ierr = MPI_Isend(send_buffer.data(),
MPI_BYTE,
rank,
tag_request,
- this->comm,
+ comm,
&send_request_requests[i]);
AssertThrowMPI(ierr);
}
return n_sources;
#else
+ (void)targets;
+ (void)create_request;
+ (void)comm;
return 0;
#endif
}
template <typename T1, typename T2>
void
- PEX<T1, T2>::answer_one_request(const unsigned int index)
+ PEX<T1, T2>::answer_one_request(
+ const unsigned int index,
+ const std::function<std::vector<T2>(const unsigned int,
+ const std::vector<T1> &)>
+ & answer_request,
+ const MPI_Comm &comm)
{
#ifdef DEAL_II_WITH_MPI
const int tag_request = Utilities::MPI::internal::Tags::
// Wait until we have a message ready for retrieval, though we don't
// care which process it is from.
MPI_Status status;
- int ierr = MPI_Probe(MPI_ANY_SOURCE, tag_request, this->comm, &status);
+ int ierr = MPI_Probe(MPI_ANY_SOURCE, tag_request, comm, &status);
AssertThrowMPI(ierr);
// Get rank of incoming message and verify that it makes sense
MPI_BYTE,
other_rank,
tag_request,
- this->comm,
+ comm,
&status);
AssertThrowMPI(ierr);
// Process request by asking the user-provided function for
// the answer and post a send for it.
auto &request_buffer = requests_buffers[index];
- this->process.answer_request(other_rank, buffer_recv, request_buffer);
+ request_buffer = answer_request(other_rank, buffer_recv);
ierr = MPI_Isend(request_buffer.data(),
request_buffer.size() * sizeof(T2),
MPI_BYTE,
other_rank,
tag_deliver,
- this->comm,
+ comm,
&send_answer_requests[index]);
AssertThrowMPI(ierr);
#else
+ (void)answer_request;
+ (void)comm;
(void)index;
#endif
}
template <typename T1, typename T2>
void
- PEX<T1, T2>::process_incoming_answers()
+ PEX<T1, T2>::process_incoming_answers(
+ const unsigned int n_targets,
+ const std::function<void(const unsigned int, const std::vector<T2> &)>
+ & process_answer,
+ const MPI_Comm &comm)
{
#ifdef DEAL_II_WITH_MPI
const int tag_deliver = Utilities::MPI::internal::Tags::
// we need not process them in order -- rather, just see what
// comes in and then look at message originators' ranks and
// message sizes
- for (unsigned int i = 0; i < targets.size(); ++i)
+ for (unsigned int i = 0; i < n_targets; ++i)
{
MPI_Status status;
{
const int ierr =
- MPI_Probe(MPI_ANY_SOURCE, tag_deliver, this->comm, &status);
+ MPI_Probe(MPI_ANY_SOURCE, tag_deliver, comm, &status);
AssertThrowMPI(ierr);
}
MPI_BYTE,
other_rank,
tag_deliver,
- this->comm,
+ comm,
MPI_STATUS_IGNORE);
AssertThrowMPI(ierr);
}
- this->process.read_answer(other_rank, recv_buffer);
+ process_answer(other_rank, recv_buffer);
}
+#else
+ (void)n_targets;
+ (void)process_answer;
+ (void)comm;
#endif
}
template <typename T1, typename T2>
std::vector<unsigned int>
- Serial<T1, T2>::run()
+ Serial<T1, T2>::run(
+ const std::vector<unsigned int> &targets,
+ const std::function<std::vector<T1>(const unsigned int)>
+ &create_request,
+ const std::function<std::vector<T2>(const unsigned int,
+ const std::vector<T1> &)>
+ &answer_request,
+ const std::function<void(const unsigned int, const std::vector<T2> &)>
+ & process_answer,
+ const MPI_Comm &comm)
{
- const auto targets = this->process.compute_targets();
+ (void)comm;
+ Assert((Utilities::MPI::job_supports_mpi() == false) ||
+ (Utilities::MPI::n_mpi_processes(comm) == 1),
+ ExcMessage("You shouldn't use the 'Serial' class on "
+ "communicators that have more than one process "
+ "associated with it."));
// The only valid target for a serial program is itself.
if (targets.size() != 0)
// Since the caller indicates that there is a target, and since we
// know that it is the current process, let the process send
// something to itself.
- std::vector<T1> send_buffer;
- std::vector<T2> recv_buffer;
- std::vector<T2> request_buffer;
-
- this->process.create_request(0, send_buffer);
- this->process.answer_request(0, send_buffer, request_buffer);
- recv_buffer = request_buffer;
- this->process.read_answer(0, recv_buffer);
+ const std::vector<T1> request = create_request(0);
+ const std::vector<T2> answer = answer_request(0, request);
+ process_answer(0, answer);
}
return targets; // nothing to do
template <typename T1, typename T2>
Selector<T1, T2>::Selector(Process<T1, T2> &process, const MPI_Comm &comm)
: Interface<T1, T2>(process, comm)
+ {}
+
+
+
+ template <typename T1, typename T2>
+ std::vector<unsigned int>
+ Selector<T1, T2>::run(
+ const std::vector<unsigned int> &targets,
+ const std::function<std::vector<T1>(const unsigned int)>
+ &create_request,
+ const std::function<std::vector<T2>(const unsigned int,
+ const std::vector<T1> &)>
+ &answer_request,
+ const std::function<void(const unsigned int, const std::vector<T2> &)>
+ & process_answer,
+ const MPI_Comm &comm)
{
// Depending on the number of processes we switch between
// implementations. We reduce the threshold for debug mode to be
# else
if (n_procs > 99)
# endif
- consensus_algo.reset(new NBX<T1, T2>(process, comm));
+ consensus_algo.reset(new NBX<T1, T2>());
else
# endif
#endif
if (n_procs > 1)
- consensus_algo.reset(new PEX<T1, T2>(process, comm));
+ consensus_algo.reset(new PEX<T1, T2>());
else
- consensus_algo.reset(new Serial<T1, T2>(process, comm));
- }
-
-
+ consensus_algo.reset(new Serial<T1, T2>());
- template <typename T1, typename T2>
- std::vector<unsigned int>
- Selector<T1, T2>::run()
- {
- return consensus_algo->run();
+ return consensus_algo->run(
+ targets, create_request, answer_request, process_answer, comm);
}