1cd05a4c0SHong Zhang #define PETSC_DLL 2cd05a4c0SHong Zhang /* 3cd05a4c0SHong Zhang Provides utility routines for split MPI communicator. 4cd05a4c0SHong Zhang */ 5d382aafbSBarry Smith #include "petscsys.h" /*I "petscsys.h" I*/ 6cd05a4c0SHong Zhang 7638faf0bSHong Zhang extern PetscErrorCode PetscSubcommCreate_contiguous(MPI_Comm,PetscInt,PetscSubcomm*); 8638faf0bSHong Zhang extern PetscErrorCode PetscSubcommCreate_interlaced(MPI_Comm,PetscInt,PetscSubcomm*); 9638faf0bSHong Zhang 10cd05a4c0SHong Zhang #undef __FUNCT__ 11cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy" 12c540e29cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommDestroy(PetscSubcomm psubcomm) 13cd05a4c0SHong Zhang { 14cd05a4c0SHong Zhang PetscErrorCode ierr; 15cd05a4c0SHong Zhang 16cd05a4c0SHong Zhang PetscFunctionBegin; 1751e7a19cSBarry Smith ierr = MPI_Comm_free(&psubcomm->dupparent);CHKERRQ(ierr); 1851e7a19cSBarry Smith ierr = MPI_Comm_free(&psubcomm->comm);CHKERRQ(ierr); 19cd05a4c0SHong Zhang ierr = PetscFree(psubcomm);CHKERRQ(ierr); 20cd05a4c0SHong Zhang PetscFunctionReturn(0); 21cd05a4c0SHong Zhang } 22cd05a4c0SHong Zhang 23cd05a4c0SHong Zhang #undef __FUNCT__ 24cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate" 25ab8c242fSMatthew Knepley /*@C 26cd05a4c0SHong Zhang PetscSubcommCreate - Create a PetscSubcomm context. 27cd05a4c0SHong Zhang 28cd05a4c0SHong Zhang Collective on MPI_Comm 29cd05a4c0SHong Zhang 30cd05a4c0SHong Zhang Input Parameter: 31cd05a4c0SHong Zhang + comm - MPI communicator 32638faf0bSHong Zhang . nsubcomm - the number of subcommunicators to be created 33638faf0bSHong Zhang - subcommtype - subcommunicator type 34cd05a4c0SHong Zhang 35cd05a4c0SHong Zhang Output Parameter: 36cd05a4c0SHong Zhang . psubcomm - location to store the PetscSubcomm context 37cd05a4c0SHong Zhang 38638faf0bSHong Zhang Level: advanced 39cd05a4c0SHong Zhang 40638faf0bSHong Zhang .keywords: communicator, create 41638faf0bSHong Zhang 42638faf0bSHong Zhang .seealso: PetscSubcommDestroy() 43638faf0bSHong Zhang @*/ 44638faf0bSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommCreate(MPI_Comm comm,PetscInt nsubcomm,PetscSubcommType subcommtype,PetscSubcomm *psubcomm) 45638faf0bSHong Zhang { 46638faf0bSHong Zhang PetscErrorCode ierr; 47638faf0bSHong Zhang PetscMPIInt rank,size; 48638faf0bSHong Zhang 49638faf0bSHong Zhang PetscFunctionBegin; 50638faf0bSHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 51638faf0bSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 52638faf0bSHong Zhang if (nsubcomm < 1 || nsubcomm > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Num of subcommunicators %D cannot be < 1 or > input comm size %D",nsubcomm,size); 53638faf0bSHong Zhang 54638faf0bSHong Zhang if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS){ 55638faf0bSHong Zhang ierr = PetscSubcommCreate_contiguous(comm,nsubcomm,psubcomm);CHKERRQ(ierr); 56638faf0bSHong Zhang } else if (subcommtype == PETSC_SUBCOMM_INTERLACED){ 57638faf0bSHong Zhang ierr = PetscSubcommCreate_interlaced(comm,nsubcomm,psubcomm);CHKERRQ(ierr); 58638faf0bSHong Zhang } else { 59638faf0bSHong Zhang SETERRQ1(comm,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",subcommtype); 60638faf0bSHong Zhang } 61638faf0bSHong Zhang PetscFunctionReturn(0); 62638faf0bSHong Zhang } 63638faf0bSHong Zhang 64638faf0bSHong Zhang #undef __FUNCT__ 65638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced" 66638faf0bSHong Zhang PetscErrorCode PetscSubcommCreate_contiguous(MPI_Comm comm,PetscInt nsubcomm,PetscSubcomm *psubcomm) 67638faf0bSHong Zhang { 68638faf0bSHong Zhang PetscErrorCode ierr; 69638faf0bSHong Zhang PetscMPIInt rank,size,*subsize,duprank,subrank; 70638faf0bSHong Zhang PetscInt np_subcomm,nleftover,i,color,rankstart; 71638faf0bSHong Zhang MPI_Comm subcomm=0,dupcomm=0; 72638faf0bSHong Zhang PetscSubcomm psubcomm_tmp; 73638faf0bSHong Zhang 74638faf0bSHong Zhang PetscFunctionBegin; 75*55e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 76*55e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 77*55e3b8d2SHong Zhang 78638faf0bSHong Zhang /* get size of each subcommunicator */ 79638faf0bSHong Zhang ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr); 80638faf0bSHong Zhang np_subcomm = size/nsubcomm; 81638faf0bSHong Zhang nleftover = size - nsubcomm*np_subcomm; 82638faf0bSHong Zhang for (i=0; i<nsubcomm; i++){ 83638faf0bSHong Zhang subsize[i] = np_subcomm; 84638faf0bSHong Zhang if (i<nleftover) subsize[i]++; 85638faf0bSHong Zhang } 86638faf0bSHong Zhang 87638faf0bSHong Zhang /* get color and subrank of this proc */ 88638faf0bSHong Zhang rankstart = 0; 89638faf0bSHong Zhang for (i=0; i<nsubcomm; i++){ 90638faf0bSHong Zhang if ( rank >= rankstart && rank < rankstart+subsize[i]) { 91638faf0bSHong Zhang color = i; 92638faf0bSHong Zhang subrank = rank - rankstart; 93638faf0bSHong Zhang duprank = rank; 94638faf0bSHong Zhang break; 95638faf0bSHong Zhang } else { 96638faf0bSHong Zhang rankstart += subsize[i]; 97638faf0bSHong Zhang } 98638faf0bSHong Zhang } 99638faf0bSHong Zhang ierr = PetscFree(subsize);CHKERRQ(ierr); 100638faf0bSHong Zhang 101638faf0bSHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 102638faf0bSHong Zhang 103638faf0bSHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 104638faf0bSHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 105638faf0bSHong Zhang 106638faf0bSHong Zhang ierr = PetscNew(struct _n_PetscSubcomm,&psubcomm_tmp);CHKERRQ(ierr); 107638faf0bSHong Zhang psubcomm_tmp->parent = comm; 108638faf0bSHong Zhang psubcomm_tmp->dupparent = dupcomm; 109638faf0bSHong Zhang psubcomm_tmp->comm = subcomm; 110638faf0bSHong Zhang psubcomm_tmp->n = nsubcomm; 111638faf0bSHong Zhang psubcomm_tmp->color = color; 112638faf0bSHong Zhang *psubcomm = psubcomm_tmp; 113638faf0bSHong Zhang PetscFunctionReturn(0); 114638faf0bSHong Zhang } 115638faf0bSHong Zhang 116638faf0bSHong Zhang #undef __FUNCT__ 117638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced" 118638faf0bSHong Zhang /* 119638faf0bSHong Zhang Note: 120638faf0bSHong Zhang In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators 12145fc02eaSBarry Smith by iteratively taking a process into a subcommunicator. 122cd05a4c0SHong Zhang Example: size=4, nsubcomm=(*psubcomm)->n=3 123cd05a4c0SHong Zhang comm=(*psubcomm)->parent: 124cd05a4c0SHong Zhang rank: [0] [1] [2] [3] 125cd05a4c0SHong Zhang color: 0 1 2 0 126cd05a4c0SHong Zhang 127cd05a4c0SHong Zhang subcomm=(*psubcomm)->comm: 128cd05a4c0SHong Zhang subrank: [0] [0] [0] [1] 129cd05a4c0SHong Zhang 130cd05a4c0SHong Zhang dupcomm=(*psubcomm)->dupparent: 131cd05a4c0SHong Zhang duprank: [0] [2] [3] [1] 132cd05a4c0SHong Zhang 133cd05a4c0SHong Zhang Here, subcomm[color = 0] has subsize=2, owns process [0] and [3] 134cd05a4c0SHong Zhang subcomm[color = 1] has subsize=1, owns process [1] 135cd05a4c0SHong Zhang subcomm[color = 2] has subsize=1, owns process [2] 136cd05a4c0SHong Zhang dupcomm has same number of processes as comm, and its duprank maps 137cd05a4c0SHong Zhang processes in subcomm contiguously into a 1d array: 138cd05a4c0SHong Zhang duprank: [0] [1] [2] [3] 139cd05a4c0SHong Zhang rank: [0] [3] [1] [2] 140cd05a4c0SHong Zhang subcomm[0] subcomm[1] subcomm[2] 141638faf0bSHong Zhang */ 142cd05a4c0SHong Zhang 143638faf0bSHong Zhang PetscErrorCode PetscSubcommCreate_interlaced(MPI_Comm comm,PetscInt nsubcomm,PetscSubcomm *psubcomm) 144cd05a4c0SHong Zhang { 145cd05a4c0SHong Zhang PetscErrorCode ierr; 146cd05a4c0SHong Zhang PetscMPIInt rank,size,*subsize,duprank,subrank; 147cd05a4c0SHong Zhang PetscInt np_subcomm,nleftover,i,j,color; 148cd05a4c0SHong Zhang MPI_Comm subcomm=0,dupcomm=0; 149c540e29cSHong Zhang PetscSubcomm psubcomm_tmp; 150cd05a4c0SHong Zhang 151cd05a4c0SHong Zhang PetscFunctionBegin; 152*55e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 153*55e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 154*55e3b8d2SHong Zhang 155cd05a4c0SHong Zhang /* get size of each subcommunicator */ 156cd05a4c0SHong Zhang ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr); 157cd05a4c0SHong Zhang np_subcomm = size/nsubcomm; 158cd05a4c0SHong Zhang nleftover = size - nsubcomm*np_subcomm; 159cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++){ 160cd05a4c0SHong Zhang subsize[i] = np_subcomm; 161cd05a4c0SHong Zhang if (i<nleftover) subsize[i]++; 162cd05a4c0SHong Zhang } 163cd05a4c0SHong Zhang 164cd05a4c0SHong Zhang /* find color for this proc */ 165cd05a4c0SHong Zhang color = rank%nsubcomm; 166cd05a4c0SHong Zhang subrank = rank/nsubcomm; 167cd05a4c0SHong Zhang 168cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 169cd05a4c0SHong Zhang 170cd05a4c0SHong Zhang j = 0; duprank = 0; 171cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++){ 172cd05a4c0SHong Zhang if (j == color){ 173cd05a4c0SHong Zhang duprank += subrank; 174cd05a4c0SHong Zhang break; 175cd05a4c0SHong Zhang } 176cd05a4c0SHong Zhang duprank += subsize[i]; j++; 177cd05a4c0SHong Zhang } 17845fc02eaSBarry Smith ierr = PetscFree(subsize);CHKERRQ(ierr); 179cd05a4c0SHong Zhang 180cd05a4c0SHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 181cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 182cd05a4c0SHong Zhang 18304033c30SSatish Balay ierr = PetscNew(struct _n_PetscSubcomm,&psubcomm_tmp);CHKERRQ(ierr); 184cd05a4c0SHong Zhang psubcomm_tmp->parent = comm; 185cd05a4c0SHong Zhang psubcomm_tmp->dupparent = dupcomm; 186cd05a4c0SHong Zhang psubcomm_tmp->comm = subcomm; 187cd05a4c0SHong Zhang psubcomm_tmp->n = nsubcomm; 188cd05a4c0SHong Zhang psubcomm_tmp->color = color; 189cd05a4c0SHong Zhang *psubcomm = psubcomm_tmp; 190cd05a4c0SHong Zhang PetscFunctionReturn(0); 191cd05a4c0SHong Zhang } 192638faf0bSHong Zhang 193638faf0bSHong Zhang 194