1e5c89e4eSSatish Balay 2*c6db04a5SJed Brown #include <petscsys.h> /*I "petscsys.h" I*/ 3e5c89e4eSSatish Balay 4e5c89e4eSSatish Balay 5e5c89e4eSSatish Balay #undef __FUNCT__ 6e5c89e4eSSatish Balay #define __FUNCT__ "PetscGatherNumberOfMessages" 7e5c89e4eSSatish Balay /*@C 8e5c89e4eSSatish Balay PetscGatherNumberOfMessages - Computes the number of messages a node expects to receive 9e5c89e4eSSatish Balay 10e5c89e4eSSatish Balay Collective on MPI_Comm 11e5c89e4eSSatish Balay 12e5c89e4eSSatish Balay Input Parameters: 13e5c89e4eSSatish Balay + comm - Communicator 14e5c89e4eSSatish Balay . iflags - an array of integers of length sizeof(comm). A '1' in ilengths[i] represent a 15e5c89e4eSSatish Balay message from current node to ith node. Optionally PETSC_NULL 16e5c89e4eSSatish Balay - ilengths - Non zero ilengths[i] represent a message to i of length ilengths[i]. 17e5c89e4eSSatish Balay Optionally PETSC_NULL. 18e5c89e4eSSatish Balay 19e5c89e4eSSatish Balay Output Parameters: 20e5c89e4eSSatish Balay . nrecvs - number of messages received 21e5c89e4eSSatish Balay 22e5c89e4eSSatish Balay Level: developer 23e5c89e4eSSatish Balay 24e5c89e4eSSatish Balay Concepts: mpi utility 25e5c89e4eSSatish Balay 26e5c89e4eSSatish Balay Notes: 27e5c89e4eSSatish Balay With this info, the correct message lengths can be determined using 28e5c89e4eSSatish Balay PetscGatherMessageLengths() 29e5c89e4eSSatish Balay 30e5c89e4eSSatish Balay Either iflags or ilengths should be provided. If iflags is not 31e5c89e4eSSatish Balay provided (PETSC_NULL) it can be computed from ilengths. If iflags is 32e5c89e4eSSatish Balay provided, ilengths is not required. 33e5c89e4eSSatish Balay 34e5c89e4eSSatish Balay .seealso: PetscGatherMessageLengths() 35e5c89e4eSSatish Balay @*/ 367087cfbeSBarry Smith PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm comm,const PetscMPIInt iflags[],const PetscMPIInt ilengths[],PetscMPIInt *nrecvs) 37e5c89e4eSSatish Balay { 38910ba992SMatthew Knepley PetscMPIInt size,rank,*recv_buf,i,*iflags_local = PETSC_NULL,*iflags_localm = PETSC_NULL; 39e5c89e4eSSatish Balay PetscErrorCode ierr; 40e5c89e4eSSatish Balay 41e5c89e4eSSatish Balay PetscFunctionBegin; 42e5c89e4eSSatish Balay 43e5c89e4eSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 44e5c89e4eSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 45e5c89e4eSSatish Balay 46e5c89e4eSSatish Balay ierr = PetscMalloc2(size,PetscMPIInt,&recv_buf,size,PetscMPIInt,&iflags_localm);CHKERRQ(ierr); 47e5c89e4eSSatish Balay 48e5c89e4eSSatish Balay /* If iflags not provided, compute iflags from ilengths */ 49e5c89e4eSSatish Balay if (!iflags) { 50e32f2f54SBarry Smith if (!ilengths) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Either iflags or ilengths should be provided"); 51e5c89e4eSSatish Balay iflags_local = iflags_localm; 52e5c89e4eSSatish Balay for (i=0; i<size; i++) { 53e5c89e4eSSatish Balay if (ilengths[i]) iflags_local[i] = 1; 54e5c89e4eSSatish Balay else iflags_local[i] = 0; 55e5c89e4eSSatish Balay } 56e5c89e4eSSatish Balay } else { 57300a7f5bSBarry Smith iflags_local = (PetscMPIInt *) iflags; 58e5c89e4eSSatish Balay } 59e5c89e4eSSatish Balay 60e5c89e4eSSatish Balay /* Post an allreduce to determine the numer of messages the current node will receive */ 61e5c89e4eSSatish Balay ierr = MPI_Allreduce(iflags_local,recv_buf,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); 62e5c89e4eSSatish Balay *nrecvs = recv_buf[rank]; 63e5c89e4eSSatish Balay 64e5c89e4eSSatish Balay ierr = PetscFree2(recv_buf,iflags_localm);CHKERRQ(ierr); 65e5c89e4eSSatish Balay PetscFunctionReturn(0); 66e5c89e4eSSatish Balay } 67e5c89e4eSSatish Balay 68e5c89e4eSSatish Balay 69e5c89e4eSSatish Balay #undef __FUNCT__ 70e5c89e4eSSatish Balay #define __FUNCT__ "PetscGatherMessageLengths" 71e5c89e4eSSatish Balay /*@C 72e5c89e4eSSatish Balay PetscGatherMessageLengths - Computes info about messages that a MPI-node will receive, 73e5c89e4eSSatish Balay including (from-id,length) pairs for each message. 74e5c89e4eSSatish Balay 75e5c89e4eSSatish Balay Collective on MPI_Comm 76e5c89e4eSSatish Balay 77e5c89e4eSSatish Balay Input Parameters: 78e5c89e4eSSatish Balay + comm - Communicator 79e5c89e4eSSatish Balay . nsends - number of messages that are to be sent. 80e5c89e4eSSatish Balay . nrecvs - number of messages being received 81e5c89e4eSSatish Balay - ilengths - an array of integers of length sizeof(comm) 82e5c89e4eSSatish Balay a non zero ilengths[i] represent a message to i of length ilengths[i] 83e5c89e4eSSatish Balay 84e5c89e4eSSatish Balay 85e5c89e4eSSatish Balay Output Parameters: 86e5c89e4eSSatish Balay + onodes - list of node-ids from which messages are expected 87e5c89e4eSSatish Balay - olengths - corresponding message lengths 88e5c89e4eSSatish Balay 89e5c89e4eSSatish Balay Level: developer 90e5c89e4eSSatish Balay 91e5c89e4eSSatish Balay Concepts: mpi utility 92e5c89e4eSSatish Balay 93e5c89e4eSSatish Balay Notes: 94e5c89e4eSSatish Balay With this info, the correct MPI_Irecv() can be posted with the correct 95e5c89e4eSSatish Balay from-id, with a buffer with the right amount of memory required. 96e5c89e4eSSatish Balay 97e5c89e4eSSatish Balay The calling function deallocates the memory in onodes and olengths 98e5c89e4eSSatish Balay 99e5c89e4eSSatish Balay To determine nrecevs, one can use PetscGatherNumberOfMessages() 100e5c89e4eSSatish Balay 101e5c89e4eSSatish Balay .seealso: PetscGatherNumberOfMessages() 102e5c89e4eSSatish Balay @*/ 1037087cfbeSBarry Smith PetscErrorCode PetscGatherMessageLengths(MPI_Comm comm,PetscMPIInt nsends,PetscMPIInt nrecvs,const PetscMPIInt ilengths[],PetscMPIInt **onodes,PetscMPIInt **olengths) 104e5c89e4eSSatish Balay { 105e5c89e4eSSatish Balay PetscErrorCode ierr; 106e5c89e4eSSatish Balay PetscMPIInt size,tag,i,j; 107910ba992SMatthew Knepley MPI_Request *s_waits = PETSC_NULL,*r_waits = PETSC_NULL; 108910ba992SMatthew Knepley MPI_Status *w_status = PETSC_NULL; 109e5c89e4eSSatish Balay 110e5c89e4eSSatish Balay PetscFunctionBegin; 111e5c89e4eSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 112e5c89e4eSSatish Balay ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr); 113e5c89e4eSSatish Balay 114e5c89e4eSSatish Balay /* cannot use PetscMalloc3() here because in the call to MPI_Waitall() they MUST be contiguous */ 115e5c89e4eSSatish Balay ierr = PetscMalloc2(nrecvs+nsends,MPI_Request,&r_waits,nrecvs+nsends,MPI_Status,&w_status);CHKERRQ(ierr); 116e5c89e4eSSatish Balay s_waits = r_waits+nrecvs; 117e5c89e4eSSatish Balay 118e5c89e4eSSatish Balay /* Post the Irecv to get the message length-info */ 119e5c89e4eSSatish Balay ierr = PetscMalloc(nrecvs*sizeof(PetscMPIInt),olengths);CHKERRQ(ierr); 120e5c89e4eSSatish Balay for (i=0; i<nrecvs; i++) { 121e5c89e4eSSatish Balay ierr = MPI_Irecv((*olengths)+i,1,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits+i);CHKERRQ(ierr); 122e5c89e4eSSatish Balay } 123e5c89e4eSSatish Balay 124e5c89e4eSSatish Balay /* Post the Isends with the message length-info */ 125e5c89e4eSSatish Balay for (i=0,j=0; i<size; ++i) { 126e5c89e4eSSatish Balay if (ilengths[i]) { 127300a7f5bSBarry Smith ierr = MPI_Isend((void*)(ilengths+i),1,MPI_INT,i,tag,comm,s_waits+j);CHKERRQ(ierr); 128e5c89e4eSSatish Balay j++; 129e5c89e4eSSatish Balay } 130e5c89e4eSSatish Balay } 131e5c89e4eSSatish Balay 132e5c89e4eSSatish Balay /* Post waits on sends and receivs */ 133e5c89e4eSSatish Balay if (nrecvs+nsends) {ierr = MPI_Waitall(nrecvs+nsends,r_waits,w_status);CHKERRQ(ierr);} 134e5c89e4eSSatish Balay 135e5c89e4eSSatish Balay /* Pack up the received data */ 136e5c89e4eSSatish Balay ierr = PetscMalloc(nrecvs*sizeof(PetscMPIInt),onodes);CHKERRQ(ierr); 137e5c89e4eSSatish Balay for (i=0; i<nrecvs; ++i) { 138e5c89e4eSSatish Balay (*onodes)[i] = w_status[i].MPI_SOURCE; 139e5c89e4eSSatish Balay } 140e5c89e4eSSatish Balay ierr = PetscFree2(r_waits,w_status);CHKERRQ(ierr); 141e5c89e4eSSatish Balay PetscFunctionReturn(0); 142e5c89e4eSSatish Balay } 143dd6ea824SBarry Smith 144e5c89e4eSSatish Balay #undef __FUNCT__ 145e5c89e4eSSatish Balay #define __FUNCT__ "PetscGatherMessageLengths2" 146e5c89e4eSSatish Balay /*@C 147e5c89e4eSSatish Balay PetscGatherMessageLengths2 - Computes info about messages that a MPI-node will receive, 148e5c89e4eSSatish Balay including (from-id,length) pairs for each message. Same functionality as PetscGatherMessageLengths() 149e5c89e4eSSatish Balay except it takes TWO ilenths and output TWO olengths. 150e5c89e4eSSatish Balay 151e5c89e4eSSatish Balay Collective on MPI_Comm 152e5c89e4eSSatish Balay 153e5c89e4eSSatish Balay Input Parameters: 154e5c89e4eSSatish Balay + comm - Communicator 155e5c89e4eSSatish Balay . nsends - number of messages that are to be sent. 156e5c89e4eSSatish Balay . nrecvs - number of messages being received 157e5c89e4eSSatish Balay - ilengths1, ilengths2 - array of integers of length sizeof(comm) 158e5c89e4eSSatish Balay a non zero ilengths[i] represent a message to i of length ilengths[i] 159e5c89e4eSSatish Balay 160e5c89e4eSSatish Balay Output Parameters: 161e5c89e4eSSatish Balay + onodes - list of node-ids from which messages are expected 162e5c89e4eSSatish Balay - olengths1, olengths2 - corresponding message lengths 163e5c89e4eSSatish Balay 164e5c89e4eSSatish Balay Level: developer 165e5c89e4eSSatish Balay 166e5c89e4eSSatish Balay Concepts: mpi utility 167e5c89e4eSSatish Balay 168e5c89e4eSSatish Balay Notes: 169e5c89e4eSSatish Balay With this info, the correct MPI_Irecv() can be posted with the correct 170e5c89e4eSSatish Balay from-id, with a buffer with the right amount of memory required. 171e5c89e4eSSatish Balay 172e5c89e4eSSatish Balay The calling function deallocates the memory in onodes and olengths 173e5c89e4eSSatish Balay 174e5c89e4eSSatish Balay To determine nrecevs, one can use PetscGatherNumberOfMessages() 175e5c89e4eSSatish Balay 176e5c89e4eSSatish Balay .seealso: PetscGatherMessageLengths() and PetscGatherNumberOfMessages() 177e5c89e4eSSatish Balay @*/ 1787087cfbeSBarry Smith PetscErrorCode PetscGatherMessageLengths2(MPI_Comm comm,PetscMPIInt nsends,PetscMPIInt nrecvs,const PetscMPIInt ilengths1[],const PetscMPIInt ilengths2[],PetscMPIInt **onodes,PetscMPIInt **olengths1,PetscMPIInt **olengths2) 179e5c89e4eSSatish Balay { 180e5c89e4eSSatish Balay PetscErrorCode ierr; 181910ba992SMatthew Knepley PetscMPIInt size,tag,i,j,*buf_s = PETSC_NULL,*buf_r = PETSC_NULL,*buf_j = PETSC_NULL; 182910ba992SMatthew Knepley MPI_Request *s_waits = PETSC_NULL,*r_waits = PETSC_NULL; 183910ba992SMatthew Knepley MPI_Status *w_status = PETSC_NULL; 184e5c89e4eSSatish Balay 185e5c89e4eSSatish Balay PetscFunctionBegin; 186e5c89e4eSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 187e5c89e4eSSatish Balay ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr); 188e5c89e4eSSatish Balay 189e5c89e4eSSatish Balay /* cannot use PetscMalloc5() because r_waits and s_waits must be contiquous for the call to MPI_Waitall() */ 190e5c89e4eSSatish Balay ierr = PetscMalloc4(nrecvs+nsends,MPI_Request,&r_waits,2*nrecvs,PetscMPIInt,&buf_r,2*nsends,PetscMPIInt,&buf_s,nrecvs+nsends,MPI_Status,&w_status);CHKERRQ(ierr); 191e5c89e4eSSatish Balay s_waits = r_waits + nrecvs; 192e5c89e4eSSatish Balay 193e5c89e4eSSatish Balay /* Post the Irecv to get the message length-info */ 194e5c89e4eSSatish Balay ierr = PetscMalloc((nrecvs+1)*sizeof(PetscMPIInt),olengths1);CHKERRQ(ierr); 195e5c89e4eSSatish Balay ierr = PetscMalloc((nrecvs+1)*sizeof(PetscMPIInt),olengths2);CHKERRQ(ierr); 196e5c89e4eSSatish Balay for (i=0; i<nrecvs; i++) { 197e5c89e4eSSatish Balay buf_j = buf_r + (2*i); 198e5c89e4eSSatish Balay ierr = MPI_Irecv(buf_j,2,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits+i);CHKERRQ(ierr); 199e5c89e4eSSatish Balay } 200e5c89e4eSSatish Balay 201e5c89e4eSSatish Balay /* Post the Isends with the message length-info */ 202e5c89e4eSSatish Balay for (i=0,j=0; i<size; ++i) { 203e5c89e4eSSatish Balay if (ilengths1[i]) { 204e5c89e4eSSatish Balay buf_j = buf_s + (2*j); 205e5c89e4eSSatish Balay buf_j[0] = *(ilengths1+i); 206e5c89e4eSSatish Balay buf_j[1] = *(ilengths2+i); 207e5c89e4eSSatish Balay ierr = MPI_Isend(buf_j,2,MPI_INT,i,tag,comm,s_waits+j);CHKERRQ(ierr); 208e5c89e4eSSatish Balay j++; 209e5c89e4eSSatish Balay } 210e5c89e4eSSatish Balay } 211e5c89e4eSSatish Balay 212e5c89e4eSSatish Balay /* Post waits on sends and receivs */ 213e5c89e4eSSatish Balay if (nrecvs+nsends) {ierr = MPI_Waitall(nrecvs+nsends,r_waits,w_status);CHKERRQ(ierr);} 214e5c89e4eSSatish Balay 215e5c89e4eSSatish Balay 216e5c89e4eSSatish Balay /* Pack up the received data */ 217e5c89e4eSSatish Balay ierr = PetscMalloc((nrecvs+1)*sizeof(PetscMPIInt),onodes);CHKERRQ(ierr); 218e5c89e4eSSatish Balay for (i=0; i<nrecvs; ++i) { 219e5c89e4eSSatish Balay (*onodes)[i] = w_status[i].MPI_SOURCE; 220e5c89e4eSSatish Balay buf_j = buf_r + (2*i); 221e5c89e4eSSatish Balay (*olengths1)[i] = buf_j[0]; 222e5c89e4eSSatish Balay (*olengths2)[i] = buf_j[1]; 223e5c89e4eSSatish Balay } 224e5c89e4eSSatish Balay 225e5c89e4eSSatish Balay ierr = PetscFree4(r_waits,buf_r,buf_s,w_status);CHKERRQ(ierr); 226e5c89e4eSSatish Balay PetscFunctionReturn(0); 227e5c89e4eSSatish Balay } 228e5c89e4eSSatish Balay 229e5c89e4eSSatish Balay /* 230e5c89e4eSSatish Balay 231e5c89e4eSSatish Balay Allocate a bufffer sufficient to hold messages of size specified in olengths. 232e5c89e4eSSatish Balay And post Irecvs on these buffers using node info from onodes 233e5c89e4eSSatish Balay 234e5c89e4eSSatish Balay */ 235e5c89e4eSSatish Balay #undef __FUNCT__ 236e5c89e4eSSatish Balay #define __FUNCT__ "PetscPostIrecvInt" 2377087cfbeSBarry Smith PetscErrorCode PetscPostIrecvInt(MPI_Comm comm,PetscMPIInt tag,PetscMPIInt nrecvs,const PetscMPIInt onodes[],const PetscMPIInt olengths[],PetscInt ***rbuf,MPI_Request **r_waits) 238e5c89e4eSSatish Balay { 239e5c89e4eSSatish Balay PetscErrorCode ierr; 240c05d87d6SBarry Smith PetscInt **rbuf_t,i,len = 0; 241e5c89e4eSSatish Balay MPI_Request *r_waits_t; 242e5c89e4eSSatish Balay 243e5c89e4eSSatish Balay PetscFunctionBegin; 244e5c89e4eSSatish Balay /* compute memory required for recv buffers */ 245e5c89e4eSSatish Balay for (i=0; i<nrecvs; i++) len += olengths[i]; /* each message length */ 246e5c89e4eSSatish Balay 247e5c89e4eSSatish Balay /* allocate memory for recv buffers */ 248c05d87d6SBarry Smith ierr = PetscMalloc((nrecvs+1)*sizeof(PetscInt*),&rbuf_t);CHKERRQ(ierr); 249c05d87d6SBarry Smith ierr = PetscMalloc(len*sizeof(PetscInt),&rbuf_t[0]);CHKERRQ(ierr); 250e5c89e4eSSatish Balay for (i=1; i<nrecvs; ++i) rbuf_t[i] = rbuf_t[i-1] + olengths[i-1]; 251e5c89e4eSSatish Balay 252e5c89e4eSSatish Balay /* Post the receives */ 253533163c2SBarry Smith ierr = PetscMalloc(nrecvs*sizeof(MPI_Request),&r_waits_t);CHKERRQ(ierr); 254e5c89e4eSSatish Balay for (i=0; i<nrecvs; ++i) { 255e5c89e4eSSatish Balay ierr = MPI_Irecv(rbuf_t[i],olengths[i],MPIU_INT,onodes[i],tag,comm,r_waits_t+i);CHKERRQ(ierr); 256e5c89e4eSSatish Balay } 257e5c89e4eSSatish Balay 258e5c89e4eSSatish Balay *rbuf = rbuf_t; 259e5c89e4eSSatish Balay *r_waits = r_waits_t; 260e5c89e4eSSatish Balay PetscFunctionReturn(0); 261e5c89e4eSSatish Balay } 262e5c89e4eSSatish Balay 263e5c89e4eSSatish Balay #undef __FUNCT__ 264e5c89e4eSSatish Balay #define __FUNCT__ "PetscPostIrecvScalar" 2657087cfbeSBarry Smith PetscErrorCode PetscPostIrecvScalar(MPI_Comm comm,PetscMPIInt tag,PetscMPIInt nrecvs,const PetscMPIInt onodes[],const PetscMPIInt olengths[],PetscScalar ***rbuf,MPI_Request **r_waits) 266e5c89e4eSSatish Balay { 267e5c89e4eSSatish Balay PetscErrorCode ierr; 268052f0c41SBarry Smith PetscMPIInt i; 269e5c89e4eSSatish Balay PetscScalar **rbuf_t; 270e5c89e4eSSatish Balay MPI_Request *r_waits_t; 271c05d87d6SBarry Smith PetscInt len = 0; 272e5c89e4eSSatish Balay 273fe28d99cSBarry Smith PetscFunctionBegin; 274e5c89e4eSSatish Balay /* compute memory required for recv buffers */ 275e5c89e4eSSatish Balay for (i=0; i<nrecvs; i++) len += olengths[i]; /* each message length */ 276e5c89e4eSSatish Balay 277e5c89e4eSSatish Balay /* allocate memory for recv buffers */ 278533163c2SBarry Smith ierr = PetscMalloc((nrecvs+1)*sizeof(PetscScalar*),&rbuf_t);CHKERRQ(ierr); 279533163c2SBarry Smith ierr = PetscMalloc(len*sizeof(PetscScalar),&rbuf_t[0]);CHKERRQ(ierr); 280e5c89e4eSSatish Balay for (i=1; i<nrecvs; ++i) rbuf_t[i] = rbuf_t[i-1] + olengths[i-1]; 281e5c89e4eSSatish Balay 282e5c89e4eSSatish Balay /* Post the receives */ 283533163c2SBarry Smith ierr = PetscMalloc(nrecvs*sizeof(MPI_Request),&r_waits_t);CHKERRQ(ierr); 284e5c89e4eSSatish Balay for (i=0; i<nrecvs; ++i) { 285e5c89e4eSSatish Balay ierr = MPI_Irecv(rbuf_t[i],olengths[i],MPIU_SCALAR,onodes[i],tag,comm,r_waits_t+i);CHKERRQ(ierr); 286e5c89e4eSSatish Balay } 287e5c89e4eSSatish Balay 288e5c89e4eSSatish Balay *rbuf = rbuf_t; 289e5c89e4eSSatish Balay *r_waits = r_waits_t; 290e5c89e4eSSatish Balay PetscFunctionReturn(0); 291e5c89e4eSSatish Balay } 292