1cd05a4c0SHong Zhang #define PETSC_DLL 2cd05a4c0SHong Zhang /* 3cd05a4c0SHong Zhang Provides utility routines for split MPI communicator. 4cd05a4c0SHong Zhang */ 5cd05a4c0SHong Zhang #include "petsc.h" /*I "petsc.h" I*/ 6*ab8c242fSMatthew Knepley #include "petscsys.h" /*I "petscsys.h" I*/ 7cd05a4c0SHong Zhang 8cd05a4c0SHong Zhang #undef __FUNCT__ 9cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy" 10cd05a4c0SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommDestroy(PetscSubcomm *psubcomm) 11cd05a4c0SHong Zhang { 12cd05a4c0SHong Zhang PetscErrorCode ierr; 13cd05a4c0SHong Zhang 14cd05a4c0SHong Zhang PetscFunctionBegin; 15cd05a4c0SHong Zhang ierr = PetscFree(psubcomm);CHKERRQ(ierr); 16cd05a4c0SHong Zhang PetscFunctionReturn(0); 17cd05a4c0SHong Zhang } 18cd05a4c0SHong Zhang 19cd05a4c0SHong Zhang #undef __FUNCT__ 20cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate" 21*ab8c242fSMatthew Knepley /*@C 22cd05a4c0SHong Zhang PetscSubcommCreate - Create a PetscSubcomm context. 23cd05a4c0SHong Zhang 24cd05a4c0SHong Zhang Collective on MPI_Comm 25cd05a4c0SHong Zhang 26cd05a4c0SHong Zhang Input Parameter: 27cd05a4c0SHong Zhang + comm - MPI communicator 28cd05a4c0SHong Zhang - nsubcomm - the number of subcommunicators to be created 29cd05a4c0SHong Zhang 30cd05a4c0SHong Zhang Output Parameter: 31cd05a4c0SHong Zhang . psubcomm - location to store the PetscSubcomm context 32cd05a4c0SHong Zhang 33cd05a4c0SHong Zhang Options Database Keys: 34cd05a4c0SHong Zhang 35cd05a4c0SHong Zhang Notes: 36cd05a4c0SHong Zhang To avoid data scattering from subcomm back to original comm, we create subcommunicators 37cd05a4c0SHong Zhang by iteratively taking a processe into a subcommunicator. 38cd05a4c0SHong Zhang Example: size=4, nsubcomm=(*psubcomm)->n=3 39cd05a4c0SHong Zhang comm=(*psubcomm)->parent: 40cd05a4c0SHong Zhang rank: [0] [1] [2] [3] 41cd05a4c0SHong Zhang color: 0 1 2 0 42cd05a4c0SHong Zhang 43cd05a4c0SHong Zhang subcomm=(*psubcomm)->comm: 44cd05a4c0SHong Zhang subrank: [0] [0] [0] [1] 45cd05a4c0SHong Zhang 46cd05a4c0SHong Zhang dupcomm=(*psubcomm)->dupparent: 47cd05a4c0SHong Zhang duprank: [0] [2] [3] [1] 48cd05a4c0SHong Zhang 49cd05a4c0SHong Zhang Here, subcomm[color = 0] has subsize=2, owns process [0] and [3] 50cd05a4c0SHong Zhang subcomm[color = 1] has subsize=1, owns process [1] 51cd05a4c0SHong Zhang subcomm[color = 2] has subsize=1, owns process [2] 52cd05a4c0SHong Zhang dupcomm has same number of processes as comm, and its duprank maps 53cd05a4c0SHong Zhang processes in subcomm contiguously into a 1d array: 54cd05a4c0SHong Zhang duprank: [0] [1] [2] [3] 55cd05a4c0SHong Zhang rank: [0] [3] [1] [2] 56cd05a4c0SHong Zhang subcomm[0] subcomm[1] subcomm[2] 57cd05a4c0SHong Zhang 58cd05a4c0SHong Zhang Level: advanced 59cd05a4c0SHong Zhang 60cd05a4c0SHong Zhang .keywords: communicator, create 61cd05a4c0SHong Zhang 62cd05a4c0SHong Zhang .seealso: PetscSubcommDestroy() 63cd05a4c0SHong Zhang @*/ 64cd05a4c0SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommCreate(MPI_Comm comm,PetscInt nsubcomm,PetscSubcomm **psubcomm) 65cd05a4c0SHong Zhang { 66cd05a4c0SHong Zhang PetscErrorCode ierr; 67cd05a4c0SHong Zhang PetscMPIInt rank,size,*subsize,duprank,subrank; 68cd05a4c0SHong Zhang PetscInt np_subcomm,nleftover,i,j,color; 69cd05a4c0SHong Zhang MPI_Comm subcomm=0,dupcomm=0; 70cd05a4c0SHong Zhang PetscSubcomm *psubcomm_tmp; 71cd05a4c0SHong Zhang 72cd05a4c0SHong Zhang PetscFunctionBegin; 73cd05a4c0SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 74cd05a4c0SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 75cd05a4c0SHong Zhang if (nsubcomm < 1 || nsubcomm > size) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Num of subcommunicators %D cannot be < 1 or > input comm size %D",nsubcomm,size); 76cd05a4c0SHong Zhang 77cd05a4c0SHong Zhang /* get size of each subcommunicator */ 78cd05a4c0SHong Zhang ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr); 79cd05a4c0SHong Zhang np_subcomm = size/nsubcomm; 80cd05a4c0SHong Zhang nleftover = size - nsubcomm*np_subcomm; 81cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++){ 82cd05a4c0SHong Zhang subsize[i] = np_subcomm; 83cd05a4c0SHong Zhang if (i<nleftover) subsize[i]++; 84cd05a4c0SHong Zhang } 85cd05a4c0SHong Zhang 86cd05a4c0SHong Zhang /* find color for this proc */ 87cd05a4c0SHong Zhang color = rank%nsubcomm; 88cd05a4c0SHong Zhang subrank = rank/nsubcomm; 89cd05a4c0SHong Zhang 90cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 91cd05a4c0SHong Zhang 92cd05a4c0SHong Zhang j = 0; duprank = 0; 93cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++){ 94cd05a4c0SHong Zhang if (j == color){ 95cd05a4c0SHong Zhang duprank += subrank; 96cd05a4c0SHong Zhang break; 97cd05a4c0SHong Zhang } 98cd05a4c0SHong Zhang duprank += subsize[i]; j++; 99cd05a4c0SHong Zhang } 100cd05a4c0SHong Zhang 101cd05a4c0SHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 102cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 103cd05a4c0SHong Zhang ierr = PetscFree(subsize);CHKERRQ(ierr); 104cd05a4c0SHong Zhang 105cd05a4c0SHong Zhang ierr = PetscNew(PetscSubcomm,&psubcomm_tmp);CHKERRQ(ierr); 106cd05a4c0SHong Zhang psubcomm_tmp->parent = comm; 107cd05a4c0SHong Zhang psubcomm_tmp->dupparent = dupcomm; 108cd05a4c0SHong Zhang psubcomm_tmp->comm = subcomm; 109cd05a4c0SHong Zhang psubcomm_tmp->n = nsubcomm; 110cd05a4c0SHong Zhang psubcomm_tmp->color = color; 111cd05a4c0SHong Zhang *psubcomm = psubcomm_tmp; 112cd05a4c0SHong Zhang PetscFunctionReturn(0); 113cd05a4c0SHong Zhang } 114