Rookie HPC

Collectives

Non-blocking

Reduction

C | FORTRAN-legacy | FORTRAN-2008

MPI_Ireduce

Definition

MPI_Ireduce is the non-blocking version of MPI_Reduce; it is the means by which MPI processes can apply a reduction calculation. The values sent by the MPI processes will be combined using the reduction operation given and the result will be stored on the MPI process specified as root. Unlike MPI_Reduce however, MPI_Ireduce returns immediately, before the reduction is guaranteed to be complete. The user must therefore explicitly wait (MPI_Wait) or test (MPI_Test) for the completion of MPI_Ireduce before safely reusing the buffers passed. Also, MPI_Ireduce is a collective operation; it must be called by every MPI process in the communicator given. Predefined operations are: MPI_MIN, MPI_MAX, MPI_BOR, MPI_BXOR, MPI_LOR, MPI_LXOR, MPI_BAND, MPI_LAND, MPI_SUM and MPI_PROD. Other variants of MPI_Ireduce are MPI_Reduce, MPI_Allreduce, MPI_Iallreduce. Refer to MPI_Reduce to see the blocking counterpart of MPI_Ireduce.

Copy

Feedback

int MPI_Ireduce(const void* send_buffer,
                void* receive_buffer,
                int count,
                MPI_Datatype datatype,
                MPI_Op operation,
                int root,
                MPI_Comm communicator,
                MPI_Request* request);

Parameters

send_buffer
A pointer on the buffer to send for reduction.
receive_buffer
A pointer on the buffer in which store the result of the reduction. Only the MPI process with the specified as root will receive the reduction result.
count
The number of elements in the send buffer, which is identical to that in the receive buffer as well.
datatype
The type of a buffer element.
operation
The operation to apply to combine messages received in the reduction. This operation must be associative, and commutative for predefined operations while user-defined operations may be non-commutative.
root
The rank of the MPI process that will collect the reduction result.
communicator
The communicator in which the reduction takes place.
request
The variable in which store the handler on the non-blocking operation.

Returned value

The error code returned from the non-blocking reduction.

MPI_SUCCESS
The routine successfully completed.

Example

Copy

Feedback

#include <stdio.h>
#include <stdlib.h>
#include <mpi.h>

/**
 * @brief Illustrates how to use a non-blocking reduce.
 * @details This application consists of a sum reduction; every MPI process
 * sends its rank for reduction before the sum of these ranks is stored in the
 * root MPI process. It can be visualised as follows, with MPI process 0 as
 * root:
 *
 * +-----------+ +-----------+ +-----------+ +-----------+
 * | Process 0 | | Process 1 | | Process 2 | | Process 3 |
 * +-+-------+-+ +-+-------+-+ +-+-------+-+ +-+-------+-+
 *   | Value |     | Value |     | Value |     | Value |
 *   |   0   |     |   1   |     |   2   |     |   3   |
 *   +-------+     +-------+     +-------+     +-------+
 *            \         |           |         /
 *             \        |           |        /
 *              \       |           |       /
 *               \      |           |      /
 *                +-----+-----+-----+-----+
 *                            |
 *                        +---+---+
 *                        |  SUM  |
 *                        +---+---+
 *                            |
 *                        +---+---+
 *                        |   6   |
 *                      +-+-------+-+
 *                      | Process 0 |
 *                      +-----------+
 **/
int main(int argc, char* argv[])
{
    MPI_Init(&argc, &argv);

    // Determine root's rank
    int root_rank = 0;

    // Get my rank
    int my_rank;
    MPI_Comm_rank(MPI_COMM_WORLD, &my_rank);

    // Each MPI process sends its rank to reduction, root MPI process collects the result
    int reduction_result = 0;
    MPI_Request request;
    MPI_Ireduce(&my_rank, &reduction_result, 1, MPI_INT, MPI_SUM, root_rank, MPI_COMM_WORLD, &request);

    // Do some other job
    printf("Process %d issued the MPI_Ireduce and has moved on, printing this message.\n", my_rank);

    // Wait for the MPI_Ireduce to complete
    MPI_Wait(&request, MPI_STATUS_IGNORE);

    if(my_rank == root_rank)
    {
        printf("The sum of all ranks is %d.\n", reduction_result);
    }

    MPI_Finalize();

    return EXIT_SUCCESS;
}