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 7*d8a68f86SHong Zhang extern PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm); 8*d8a68f86SHong Zhang extern PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm); 9*d8a68f86SHong Zhang 10*d8a68f86SHong Zhang #undef __FUNCT__ 11*d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetNumber" 12*d8a68f86SHong Zhang /*@C 13*d8a68f86SHong Zhang PetscSubcommSetNumber - Set total number of subcommunicators. 14*d8a68f86SHong Zhang 15*d8a68f86SHong Zhang Collective on MPI_Comm 16*d8a68f86SHong Zhang 17*d8a68f86SHong Zhang Input Parameter: 18*d8a68f86SHong Zhang + psubcomm - PetscSubcomm context 19*d8a68f86SHong Zhang - nsubcomm - the total number of subcommunicators in psubcomm 20*d8a68f86SHong Zhang 21*d8a68f86SHong Zhang Level: advanced 22*d8a68f86SHong Zhang 23*d8a68f86SHong Zhang .keywords: communicator 24*d8a68f86SHong Zhang 25*d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetType(),PetscSubcommSetTypeGeneral() 26*d8a68f86SHong Zhang @*/ 27*d8a68f86SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommSetNumber(PetscSubcomm psubcomm,PetscInt nsubcomm) 28*d8a68f86SHong Zhang { 29*d8a68f86SHong Zhang PetscErrorCode ierr; 30*d8a68f86SHong Zhang MPI_Comm comm=psubcomm->parent; 31*d8a68f86SHong Zhang PetscMPIInt rank,size; 32*d8a68f86SHong Zhang 33*d8a68f86SHong Zhang PetscFunctionBegin; 34*d8a68f86SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate() first"); 35*d8a68f86SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 36*d8a68f86SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 37*d8a68f86SHong 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); 38*d8a68f86SHong Zhang 39*d8a68f86SHong Zhang psubcomm->n = nsubcomm; 40*d8a68f86SHong Zhang PetscFunctionReturn(0); 41*d8a68f86SHong Zhang } 42*d8a68f86SHong Zhang 43*d8a68f86SHong Zhang #undef __FUNCT__ 44*d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetType" 45*d8a68f86SHong Zhang /*@C 46*d8a68f86SHong Zhang PetscSubcommSetType - Set type of subcommunicators. 47*d8a68f86SHong Zhang 48*d8a68f86SHong Zhang Collective on MPI_Comm 49*d8a68f86SHong Zhang 50*d8a68f86SHong Zhang Input Parameter: 51*d8a68f86SHong Zhang + psubcomm - PetscSubcomm context 52*d8a68f86SHong Zhang - subcommtype - subcommunicator type 53*d8a68f86SHong Zhang 54*d8a68f86SHong Zhang Level: advanced 55*d8a68f86SHong Zhang 56*d8a68f86SHong Zhang .keywords: communicator 57*d8a68f86SHong Zhang 58*d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral() 59*d8a68f86SHong Zhang @*/ 60*d8a68f86SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommSetType(PetscSubcomm psubcomm,const PetscSubcommType subcommtype) 61*d8a68f86SHong Zhang { 62*d8a68f86SHong Zhang PetscErrorCode ierr; 63*d8a68f86SHong Zhang 64*d8a68f86SHong Zhang PetscFunctionBegin; 65*d8a68f86SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 66*d8a68f86SHong Zhang if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n); 67*d8a68f86SHong Zhang 68*d8a68f86SHong Zhang if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS){ 69*d8a68f86SHong Zhang ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr); 70*d8a68f86SHong Zhang } else if (subcommtype == PETSC_SUBCOMM_INTERLACED){ 71*d8a68f86SHong Zhang ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr); 72*d8a68f86SHong Zhang } else { 73*d8a68f86SHong Zhang SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",subcommtype); 74*d8a68f86SHong Zhang } 75*d8a68f86SHong Zhang PetscFunctionReturn(0); 76*d8a68f86SHong Zhang } 77*d8a68f86SHong Zhang 78*d8a68f86SHong Zhang #undef __FUNCT__ 79*d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetTypeGeneral" 80*d8a68f86SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm) 81*d8a68f86SHong Zhang { 82*d8a68f86SHong Zhang PetscFunctionBegin; 83*d8a68f86SHong Zhang PetscFunctionReturn(0); 84*d8a68f86SHong Zhang } 85638faf0bSHong Zhang 86cd05a4c0SHong Zhang #undef __FUNCT__ 87cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy" 88c540e29cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommDestroy(PetscSubcomm psubcomm) 89cd05a4c0SHong Zhang { 90cd05a4c0SHong Zhang PetscErrorCode ierr; 91cd05a4c0SHong Zhang 92cd05a4c0SHong Zhang PetscFunctionBegin; 9351e7a19cSBarry Smith ierr = MPI_Comm_free(&psubcomm->dupparent);CHKERRQ(ierr); 9451e7a19cSBarry Smith ierr = MPI_Comm_free(&psubcomm->comm);CHKERRQ(ierr); 95cd05a4c0SHong Zhang ierr = PetscFree(psubcomm);CHKERRQ(ierr); 96cd05a4c0SHong Zhang PetscFunctionReturn(0); 97cd05a4c0SHong Zhang } 98cd05a4c0SHong Zhang 99cd05a4c0SHong Zhang #undef __FUNCT__ 100cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate" 101ab8c242fSMatthew Knepley /*@C 102cd05a4c0SHong Zhang PetscSubcommCreate - Create a PetscSubcomm context. 103cd05a4c0SHong Zhang 104cd05a4c0SHong Zhang Collective on MPI_Comm 105cd05a4c0SHong Zhang 106cd05a4c0SHong Zhang Input Parameter: 107cd05a4c0SHong Zhang + comm - MPI communicator 108638faf0bSHong Zhang . nsubcomm - the number of subcommunicators to be created 109638faf0bSHong Zhang - subcommtype - subcommunicator type 110cd05a4c0SHong Zhang 111cd05a4c0SHong Zhang Output Parameter: 112cd05a4c0SHong Zhang . psubcomm - location to store the PetscSubcomm context 113cd05a4c0SHong Zhang 114638faf0bSHong Zhang Level: advanced 115cd05a4c0SHong Zhang 116638faf0bSHong Zhang .keywords: communicator, create 117638faf0bSHong Zhang 118638faf0bSHong Zhang .seealso: PetscSubcommDestroy() 119638faf0bSHong Zhang @*/ 120*d8a68f86SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm) 121638faf0bSHong Zhang { 122638faf0bSHong Zhang PetscErrorCode ierr; 123638faf0bSHong Zhang 124638faf0bSHong Zhang PetscFunctionBegin; 125638faf0bSHong Zhang 126*d8a68f86SHong Zhang ierr = PetscNew(struct _n_PetscSubcomm,psubcomm);CHKERRQ(ierr); 127*d8a68f86SHong Zhang (*psubcomm)->parent = comm; 128638faf0bSHong Zhang PetscFunctionReturn(0); 129638faf0bSHong Zhang } 130638faf0bSHong Zhang 131638faf0bSHong Zhang #undef __FUNCT__ 13253c77d0aSJed Brown #define __FUNCT__ "PetscSubcommCreate_contiguous" 133*d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm) 134638faf0bSHong Zhang { 135638faf0bSHong Zhang PetscErrorCode ierr; 136d6037b41SHong Zhang PetscMPIInt rank,size,*subsize,duprank=-1,subrank=-1; 137*d8a68f86SHong Zhang PetscInt np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n; 138*d8a68f86SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 139638faf0bSHong Zhang 140638faf0bSHong Zhang PetscFunctionBegin; 14155e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 14255e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 14355e3b8d2SHong Zhang 144638faf0bSHong Zhang /* get size of each subcommunicator */ 145638faf0bSHong Zhang ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr); 146638faf0bSHong Zhang np_subcomm = size/nsubcomm; 147638faf0bSHong Zhang nleftover = size - nsubcomm*np_subcomm; 148638faf0bSHong Zhang for (i=0; i<nsubcomm; i++){ 149638faf0bSHong Zhang subsize[i] = np_subcomm; 150638faf0bSHong Zhang if (i<nleftover) subsize[i]++; 151638faf0bSHong Zhang } 152638faf0bSHong Zhang 153638faf0bSHong Zhang /* get color and subrank of this proc */ 154638faf0bSHong Zhang rankstart = 0; 155638faf0bSHong Zhang for (i=0; i<nsubcomm; i++){ 156638faf0bSHong Zhang if ( rank >= rankstart && rank < rankstart+subsize[i]) { 157638faf0bSHong Zhang color = i; 158638faf0bSHong Zhang subrank = rank - rankstart; 159638faf0bSHong Zhang duprank = rank; 160638faf0bSHong Zhang break; 161638faf0bSHong Zhang } else { 162638faf0bSHong Zhang rankstart += subsize[i]; 163638faf0bSHong Zhang } 164638faf0bSHong Zhang } 165638faf0bSHong Zhang ierr = PetscFree(subsize);CHKERRQ(ierr); 166638faf0bSHong Zhang 167638faf0bSHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 168638faf0bSHong Zhang 169638faf0bSHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 170638faf0bSHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 171638faf0bSHong Zhang 172*d8a68f86SHong Zhang psubcomm->dupparent = dupcomm; 173*d8a68f86SHong Zhang psubcomm->comm = subcomm; 174*d8a68f86SHong Zhang psubcomm->color = color; 175638faf0bSHong Zhang PetscFunctionReturn(0); 176638faf0bSHong Zhang } 177638faf0bSHong Zhang 178638faf0bSHong Zhang #undef __FUNCT__ 179638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced" 180638faf0bSHong Zhang /* 181638faf0bSHong Zhang Note: 182638faf0bSHong Zhang In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators 18345fc02eaSBarry Smith by iteratively taking a process into a subcommunicator. 184cd05a4c0SHong Zhang Example: size=4, nsubcomm=(*psubcomm)->n=3 185cd05a4c0SHong Zhang comm=(*psubcomm)->parent: 186cd05a4c0SHong Zhang rank: [0] [1] [2] [3] 187cd05a4c0SHong Zhang color: 0 1 2 0 188cd05a4c0SHong Zhang 189cd05a4c0SHong Zhang subcomm=(*psubcomm)->comm: 190cd05a4c0SHong Zhang subrank: [0] [0] [0] [1] 191cd05a4c0SHong Zhang 192cd05a4c0SHong Zhang dupcomm=(*psubcomm)->dupparent: 193cd05a4c0SHong Zhang duprank: [0] [2] [3] [1] 194cd05a4c0SHong Zhang 195cd05a4c0SHong Zhang Here, subcomm[color = 0] has subsize=2, owns process [0] and [3] 196cd05a4c0SHong Zhang subcomm[color = 1] has subsize=1, owns process [1] 197cd05a4c0SHong Zhang subcomm[color = 2] has subsize=1, owns process [2] 198cd05a4c0SHong Zhang dupcomm has same number of processes as comm, and its duprank maps 199cd05a4c0SHong Zhang processes in subcomm contiguously into a 1d array: 200cd05a4c0SHong Zhang duprank: [0] [1] [2] [3] 201cd05a4c0SHong Zhang rank: [0] [3] [1] [2] 202cd05a4c0SHong Zhang subcomm[0] subcomm[1] subcomm[2] 203638faf0bSHong Zhang */ 204cd05a4c0SHong Zhang 205*d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm) 206cd05a4c0SHong Zhang { 207cd05a4c0SHong Zhang PetscErrorCode ierr; 208cd05a4c0SHong Zhang PetscMPIInt rank,size,*subsize,duprank,subrank; 209*d8a68f86SHong Zhang PetscInt np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n; 210*d8a68f86SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 211cd05a4c0SHong Zhang 212cd05a4c0SHong Zhang PetscFunctionBegin; 21355e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 21455e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 21555e3b8d2SHong Zhang 216cd05a4c0SHong Zhang /* get size of each subcommunicator */ 217cd05a4c0SHong Zhang ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr); 218cd05a4c0SHong Zhang np_subcomm = size/nsubcomm; 219cd05a4c0SHong Zhang nleftover = size - nsubcomm*np_subcomm; 220cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++){ 221cd05a4c0SHong Zhang subsize[i] = np_subcomm; 222cd05a4c0SHong Zhang if (i<nleftover) subsize[i]++; 223cd05a4c0SHong Zhang } 224cd05a4c0SHong Zhang 225cd05a4c0SHong Zhang /* find color for this proc */ 226cd05a4c0SHong Zhang color = rank%nsubcomm; 227cd05a4c0SHong Zhang subrank = rank/nsubcomm; 228cd05a4c0SHong Zhang 229cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 230cd05a4c0SHong Zhang 231cd05a4c0SHong Zhang j = 0; duprank = 0; 232cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++){ 233cd05a4c0SHong Zhang if (j == color){ 234cd05a4c0SHong Zhang duprank += subrank; 235cd05a4c0SHong Zhang break; 236cd05a4c0SHong Zhang } 237cd05a4c0SHong Zhang duprank += subsize[i]; j++; 238cd05a4c0SHong Zhang } 23945fc02eaSBarry Smith ierr = PetscFree(subsize);CHKERRQ(ierr); 240cd05a4c0SHong Zhang 241cd05a4c0SHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 242cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 243cd05a4c0SHong Zhang 244*d8a68f86SHong Zhang psubcomm->dupparent = dupcomm; 245*d8a68f86SHong Zhang psubcomm->comm = subcomm; 246*d8a68f86SHong Zhang psubcomm->color = color; 247cd05a4c0SHong Zhang PetscFunctionReturn(0); 248cd05a4c0SHong Zhang } 249638faf0bSHong Zhang 250638faf0bSHong Zhang 251