MPI Subarray Type

We used the subarray type for mpi4py in the two-dimensional exchange example due to pointer issues. We can create subarrays for Fortran and C++ as well, and they can be useful if we wish to define larger blocks than a row or column, such as a submatrix within a computational grid. For example, many algorithms can utilize more than one ghost zone in each dimension. In Fortran we could use array slices, or in C++ we could create an intermediate buffer, but both these solutions would require a copy, which can be slow. With a derived type such as the subarray, MPI can pluck its buffer directly from memory without the programmer needing to set up another array or slice.

C/C++

In this example we pick an array of size 3x4 out of the full local array, which has been declared as nrl,ncl. We start the selection at w[2][3]. Always be sure that the subsize fits within the full size.

We specify MPI_ORDER_C because we are using the normal row-major layout.

    int ndims=2;
    int sizes[]={nrl,ncl};

    int starts[ndims];
    int subsizes[ndims];

    subsizes[0]=3; subsizes[1]=4;

    MPI_Datatype sendtype;

    starts[0]=2; starts[1]=3;
    MPI_Type_create_subarray(ndims,sizes,subsizes,starts,MPI_ORDER_C,MPI_INT,&sendtype);
    MPI_Type_commit(&sendtype);

We will send this selection into a receive buffer of the same subsize. The data will be written into a different section of the local array, one that starts at w[3][2]. Be sure that the sizes of sending and receiving buffers always match.

MPI_Datatype recvtype;
starts[0]=3; starts[1]=2;
MPI_Type_create_subarray(ndims,sizes,subsizes,starts,MPI_ORDER_C,MPI_INT,&recvtype);
MPI_Type_commit(&recvtype);

When we send the subarray, we must provide the pointer to the first element in the send buffer, and similarly for the receive buffer. How we obtain this may depend on how our array was set up. In our examples throughout this course we have set up two-dimensional arrays in C++ as C-style dynamically-sized arrays of pointers. This means that the array name w just represents the first pointer-to-pointer. So we pass the pointer to the first value directly.

  MPI_Send(&w[0][0], 1, sendtype, 1, tag, MPI_COMM_WORLD);
  MPI_Recv(&w[0][0], 1, recvtype, root, tag, MPI_COMM_WORLD, MPI_STATUS_IGNORE);

The full example program sets up a 6x8 local array in each of two processes, then sends the contents of the selection in one process to the different subsection of the other process.

C++

Contents of mpi_subarray.cxx

#include <cstring>
#include <iostream> 
#include <iomanip> 
#include <mpi.h> 

using namespace std; 

// MPI subarray example

int main (int argc, char *argv[]) {

    int nprocs, rank;
    int tag=0;
    int root=0;

    //Initialize MPI
    MPI_Init(&argc, &argv);
    MPI_Comm_size(MPI_COMM_WORLD,&nprocs);
    MPI_Comm_rank(MPI_COMM_WORLD,&rank);

    //Run with two processes only
    if ( nprocs != 2 ) {
        if ( rank==root) { 
            cout<<"Please run with two processes"<<endl; 
            MPI_Finalize();
            exit(1);
        }
    }

    int nrl=6;
    int ncl=8;

    int **w=new int*[nrl];
    int *wptr=new int[(nrl)*(ncl)];

    for (int i=0;i<nrl;++i,wptr+=ncl) {
       w[i] = wptr;
    }

    if (rank==root) {
        for ( int i = 0; i < nrl; i++ ) {
             for (int j = 0; j < ncl; j++ ) {
                 w[i][j] = i*100+j;
             }
        }
    }
    else {
        for ( int i = 0; i < nrl; i++ ) {
             for (int j = 0; j < ncl; j++ ) {
                 w[i][j] = 0;
             }
        }
    }

    //Set up subarrays
    //
    int ndims=2;
    int sizes[]={nrl,ncl};

    int starts[ndims]; 
    int subsizes[ndims];

    subsizes[0]=3; subsizes[1]=4;

    MPI_Datatype sendtype, recvtype;

    starts[0]=2; starts[1]=3;
    MPI_Type_create_subarray(ndims,sizes,subsizes,starts,MPI_ORDER_C,MPI_INT,&sendtype);
    MPI_Type_commit(&sendtype);

    starts[0]=3; starts[1]=2;
    MPI_Type_create_subarray(ndims,sizes,subsizes,starts,MPI_ORDER_C,MPI_INT,&recvtype);
    MPI_Type_commit(&recvtype);

    //Buffers set up

    // Send type from root to rank 1
    if ( rank==root ) {
        cout<<"Rank "<<rank<<" sending array:"<<endl;
        for (int i=0;i<nrl;i++) {
            for (int j=0;j<ncl;j++) {
                cout<<setfill('0')<<setw(3)<<w[i][j]<<" ";
            }
            cout<<endl;
        }
        MPI_Send(&w[0][0], 1, sendtype, 1, tag, MPI_COMM_WORLD);
    }
    else {
        MPI_Recv(&w[0][0], 1, recvtype, root, tag, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
        cout<<" Rank "<<rank<<" received subarray: "<<endl;
        for (int i=0;i<nrl;i++) {
            for (int j=0;j<ncl;j++) {
                cout<<setfill('0')<<setw(3)<<w[i][j]<<" ";
            }
            cout<<endl;
        }
    }

    //Type_free not really necessary but good practice
    MPI_Type_free(&sendtype);
    MPI_Type_free(&recvtype);

    MPI_Finalize();

    return 0;

}

Download mpi_subarray.cxx file

Fortran

Fortran will be very similar to C++ in this case, since row major versus column major is specified when we define the subarray type. When determining the starting points, we must use 0-based array indexing regardless of how our array was declared. MPI will use the starting position and ordering to calculate the elements to be packed into the buffer.

integer :: ndims 
integer, allocatable, dimension(:) :: wsizes,subsizes, starts
type(MPI_Datatype)::sendtype,recvtype

   ! Declare arrays
   ! w to nrl, ncl
   ! wsizes, subsizes, starts to ndims (2 in our examples)

    ! Set up subarrays
    wsizes = [nrl, ncl]                 ! Full array size
    subsizes = [ sub_rows, sub_cols ]   ! Subarray size

    ! Send array
    starts   = [ 2, 3 ]                 ! Starting indices

    ! Subarray type requires counting from 0
    starts=starts-lb

    call MPI_Type_create_subarray(2, wsizes, subsizes, starts, &
                                  MPI_ORDER_FORTRAN, MPI_INTEGER, sendtype)
    call MPI_Type_commit(sendtype)

    ! Receive array
    starts=[ 3,2 ]-lb

    call MPI_Type_create_subarray(2, wsizes, subsizes, starts, &
                                  MPI_ORDER_FORTRAN, MPI_INTEGER, recvtype, ierr)
    call MPI_Type_commit(recvtype)

Fortran passes all arrays by reference, i.e. by passing the pointer, so no special treatment is required in the send or receive.

   call MPI_Send(w, 1, sendtype, 1, 0, MPI_COMM_WORLD)
   call MPI_Recv(w, 1 , recvtype, root, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE)

The full program is similar to the C++ example; it sets up a 6x8 local array in each of two processes, then sends the contents of the selection in one process to the different subsection of the other process.

Fortran

Contents of mpi_subarray.f90

program mpi_subarray
    use mpi_f08
    implicit none

    integer, parameter :: root = 0
    integer :: ierr, rank, size
    integer :: nrl, ncl
    integer :: sub_rows, sub_cols
    integer :: starts(2), subsizes(2), wsizes(2)
    integer :: i, j
    integer :: lb, ubr, ubc
    integer, allocatable :: w(:,:)
    type(MPI_Datatype)::sendtype,recvtype

    call MPI_Init()
    call MPI_Comm_rank(MPI_COMM_WORLD, rank)
    call MPI_Comm_size(MPI_COMM_WORLD, size)

    ! For simplicity, run with 2 processes
    if (size /= 2) then
        if (rank == root) print *, "Please run with 2 processes."
        call MPI_Finalize()
        stop
    end if

    ! Define the full array dimensions
    nrl = 6
    ncl = 8

    ! Define the subarray dimensions
    sub_rows = 3
    sub_cols = 4

    lb=0
    ubr=nrl-1
    ubc=ncl-1

    allocate(w(lb:ubr, lb:ubc))

    ! Fill the arrays.  All zeros on rank 1.
    if (rank == root) then
        do i = lb, ubr
            do j = lb, ubc
                w(i, j) = i * 100 + j
            end do
        end do
    else
        w=0
    end if

    ! Set up subarrays
    wsizes = [nrl, ncl]                 ! Full array size
    subsizes = [ sub_rows, sub_cols ]   ! Subarray size

    ! Send array
    starts   = [ 2, 3 ]                 ! Starting indices

    ! Subarray type requires counting from 0
    starts=starts-lb

    call MPI_Type_create_subarray(2, wsizes, subsizes, starts, &
                                  MPI_ORDER_FORTRAN, MPI_INTEGER, sendtype)
    call MPI_Type_commit(sendtype)

    ! Receive array
    starts=[ 3,2 ]-lb

    call MPI_Type_create_subarray(2, wsizes, subsizes, starts, &
                                  MPI_ORDER_FORTRAN, MPI_INTEGER, recvtype, ierr)
    call MPI_Type_commit(recvtype)

    ! Send subarray from root to rank 1
    if (rank == root) then
        print *, "Rank", rank, "sending array:"
        do i = lb, ubr
            write(*,'(*(i6.3))') w(i,:)
        end do
        call MPI_Send(w, 1, sendtype, 1, 0, MPI_COMM_WORLD)
    else
        call MPI_Recv(w, 1 , recvtype, root, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE)
        ! Print received values
        print *, "Rank", rank, "received subarray:"
        do i = lb, ubr
            write(*,'(*(i6.3))') w(i,:)
        end do
    end if

    call MPI_Type_free(sendtype)
    call MPI_Type_free(recvtype)
    call MPI_Finalize()

end program mpi_subarray

Download mpi_subarray.f90 file

Python

The general syntax for Python subarrays is

sizes=[nrl,ncl]
subsizes=[sub_rows,sub_cols]
starts=[start_row,start_col]
newtype=MPI.INTEGER.Create_subarray(sizes,subsizes,starts)
newtype.Commit()

The complete example is similar to those for C++ and Fortran.

Python

Contents of mpi_subarray.py

import sys
import numpy as np
from mpi4py import MPI

comm = MPI.COMM_WORLD
rank = comm.Get_rank()
nprocs = comm.Get_size()
root = 0

# For simplicity, run with 2 processes
if nprocs != 2:
    print("Please run with exactly 2 processes")
    sys.exit()

nrl=6
ncl=8

# Be sure to match with MPI Type
w = np.zeros((nrl, ncl), dtype='int32')

if rank==root:
    for i in range(nrl):
        for j in range(ncl):
            w[i,j]=i*100+j

# Define the subarrays
sizes = [nrl,ncl]
subsizes=[3,4]

# Set up MPI type for send buffer
starts=[2,3]
sendtype=MPI.INTEGER.Create_subarray(sizes,subsizes,starts)
sendtype.Commit()

# Set up MPI type for receive buffer
starts=[3,2]
recvtype=MPI.INTEGER.Create_subarray(sizes,subsizes,starts)
recvtype.Commit()

# Send subarray from root to rank 1
if rank==root:
    print("Rank ",rank," sending array")
    print(w)
    comm.Send([w,1,sendtype], dest=1, tag=0)
else:
    comm.Recv([w,1,recvtype], source=root, tag=0)
    # Print received values
    print("Rank ",rank," received values")
    print(w)

sendtype.Free()
recvtype.Free()


Download mpi_subarray.py file

Previous
RC Logo RC Logo © 2026 The Rector and Visitors of the University of Virginia