static char help[] = "Demonstrates BuildTwoSided functions.\n"; #include typedef struct { PetscInt rank; PetscScalar value; char ok[3]; } Unit; static PetscErrorCode MakeDatatype(MPI_Datatype *dtype) { MPI_Datatype dtypes[3],tmptype; PetscMPIInt lengths[3]; MPI_Aint displs[3]; Unit dummy; PetscFunctionBegin; dtypes[0] = MPIU_INT; dtypes[1] = MPIU_SCALAR; dtypes[2] = MPI_CHAR; lengths[0] = 1; lengths[1] = 1; lengths[2] = 3; /* Curse the evil beings that made std::complex a non-POD type. */ displs[0] = (char*)&dummy.rank - (char*)&dummy; /* offsetof(Unit,rank); */ displs[1] = (char*)&dummy.value - (char*)&dummy; /* offsetof(Unit,value); */ displs[2] = (char*)&dummy.ok - (char*)&dummy; /* offsetof(Unit,ok); */ PetscCallMPI(MPI_Type_create_struct(3,lengths,displs,dtypes,&tmptype)); PetscCallMPI(MPI_Type_commit(&tmptype)); PetscCallMPI(MPI_Type_create_resized(tmptype,0,sizeof(Unit),dtype)); PetscCallMPI(MPI_Type_commit(dtype)); PetscCallMPI(MPI_Type_free(&tmptype)); { MPI_Aint lb,extent; PetscCallMPI(MPI_Type_get_extent(*dtype,&lb,&extent)); PetscCheck(extent == sizeof(Unit),PETSC_COMM_WORLD,PETSC_ERR_LIB,"New type has extent %d != sizeof(Unit) %d",(int)extent,(int)sizeof(Unit)); } PetscFunctionReturn(0); } struct FCtx { PetscMPIInt rank; PetscMPIInt nto; PetscMPIInt *toranks; Unit *todata; PetscSegBuffer seg; }; static PetscErrorCode FSend(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt tonum,PetscMPIInt rank,void *todata,MPI_Request req[],void *ctx) { struct FCtx *fctx = (struct FCtx*)ctx; PetscFunctionBegin; PetscCheck(rank == fctx->toranks[tonum],PETSC_COMM_SELF,PETSC_ERR_PLIB,"Rank %d does not match toranks[%d] %d",rank,tonum,fctx->toranks[tonum]); PetscCheck(fctx->rank == *(PetscMPIInt*)todata,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Todata %d does not match rank %d",*(PetscMPIInt*)todata,fctx->rank); PetscCallMPI(MPI_Isend(&fctx->todata[tonum].rank,1,MPIU_INT,rank,tag[0],comm,&req[0])); PetscCallMPI(MPI_Isend(&fctx->todata[tonum].value,1,MPIU_SCALAR,rank,tag[1],comm,&req[1])); PetscFunctionReturn(0); } static PetscErrorCode FRecv(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt rank,void *fromdata,MPI_Request req[],void *ctx) { struct FCtx *fctx = (struct FCtx*)ctx; Unit *buf; PetscFunctionBegin; PetscCheck(*(PetscMPIInt*)fromdata == rank,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Dummy data %d from rank %d corrupt",*(PetscMPIInt*)fromdata,rank); PetscCall(PetscSegBufferGet(fctx->seg,1,&buf)); PetscCallMPI(MPI_Irecv(&buf->rank,1,MPIU_INT,rank,tag[0],comm,&req[0])); PetscCallMPI(MPI_Irecv(&buf->value,1,MPIU_SCALAR,rank,tag[1],comm,&req[1])); buf->ok[0] = 'o'; buf->ok[1] = 'k'; buf->ok[2] = 0; PetscFunctionReturn(0); } int main(int argc,char **argv) { PetscMPIInt rank,size,*toranks,*fromranks,nto,nfrom; PetscInt i,n; PetscBool verbose,build_twosided_f; Unit *todata,*fromdata; MPI_Datatype dtype; PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); verbose = PETSC_FALSE; PetscCall(PetscOptionsGetBool(NULL,NULL,"-verbose",&verbose,NULL)); build_twosided_f = PETSC_FALSE; PetscCall(PetscOptionsGetBool(NULL,NULL,"-build_twosided_f",&build_twosided_f,NULL)); for (i=1,nto=0; i