/*###################################################################### Example 5 : MPI_Send and MPI_Recv Description: Examples 5 and 6 demonstrate the difference between blocking and non-blocking point-to-point communication. Example 5: MPI_Send/MPI_Recv (blocking) Example 6: MPI_Isend/MPI_Irecv (non-blocking) sendbuff recvbuff sendbuff recvbuff ######## ######## ######## ######## # # # # # # # # 0 # AA # # # # AA # # EE # # # # # # # # # ######## ######## ######## ######## T # # # # # # # # 1 # BB # # # # BB # # AA # a # # # # # # # # ######## ######## ######## ######## s # # # # # # # # 2 # CC # # # # CC # # BB # k # # # # # # # # ######## ######## ######## ######## s # # # # # # # # 3 # DD # # # # DD # # CC # # # # # # # # # ######## ######## ######## ######## # # # # # # # # 4 # EE # # # # EE # # DD # # # # # # # # # ######## ######## ######## ######## BEFORE AFTER Each task transfers a vector of random numbers (sendbuff) to the next task (taskid+1). The last task transfers it to task 0. Consequently, each task receives a vector from the preceding task and puts it in recvbuff. This example shows that MPI_Send and MPI_Recv are not the most efficient functions to perform this work. Since they work in blocking mode (i.e. waits for the transfer to finish before continuing its execution). In order to receive their vector, each task must post an MPI_Recv corresponding to the MPI_Send of the sender, and so wait for the reception to complete before sending sendbuff to the next task. In order to avoid a deadlock, one of the task must initiate the communication and post its MPI_Recv after the MPI_Send. In this example, the last task initiate this cascading process where each consecutive task is waiting for the complete reception from the preceding task before starting to send "sendbuff" to the next. The whole process completes when the last task have finished its reception. Before the communication, task 0 gather the sum of all the vectors to sent from each tasks, and prints them out. Similarly after the communication, task 0 gathers all the sum of the vectors received by each task and prints them out along with the communication times. Example 6 show how to use non-blocking communication (MPI_Isend and MPI_Irecv) to accomplish the same work is much less time. The size of the vector (buffsize) is given as an argument to the program at run time. Author: Carol Gauthier Centre de Calcul scientifique Universite de Sherbrooke Last revision: September 2005 ######################################################################*/ #include #include #include #include #include "math.h" #include "mpi.h" int main(int argc,char** argv) { int taskid, ntasks; MPI_Status status; MPI_Request req[1024]; int ierr,i,j,itask,recvtaskid; int buffsize; double *sendbuff,*recvbuff; double sendbuffsum,recvbuffsum; double sendbuffsums[1024],recvbuffsums[1024]; double inittime,totaltime,recvtime,recvtimes[1024]; /*===============================================================*/ /* MPI Initialisation. Its important to put this call at the */ /* begining of the program, after variable declarations. */ MPI_Init(&argc, &argv); /*===============================================================*/ /* Get the number of MPI tasks and the taskid of this task. */ MPI_Comm_rank(MPI_COMM_WORLD,&taskid); MPI_Comm_size(MPI_COMM_WORLD,&ntasks); /*===============================================================*/ /* Get buffsize value from program arguments. */ buffsize=atoi(argv[1]); /*===============================================================*/ /* Printing out the description of the example. */ if ( taskid == 0 ){ printf("\n\n\n"); printf("##########################################################\n\n"); printf(" Example 5 \n\n"); printf(" Point-to-point Communication: MPI_Send MPI_Recv \n\n"); printf(" Vector size: %d\n",buffsize); printf(" Number of tasks: %d\n\n",ntasks); printf("##########################################################\n\n"); printf(" --> BEFORE COMMUNICATION <--\n\n"); } /*=============================================================*/ /* Memory allocation. */ sendbuff=(double *)malloc(sizeof(double)*buffsize); recvbuff=(double *)malloc(sizeof(double)*buffsize); /*=============================================================*/ /* Vectors and/or matrices initalisation. */ srand((unsigned)time( NULL ) + taskid); for(i=0;i AFTER COMMUNICATION <-- \n\n"); for(itask=0;itask