17d0a6c19SBarry Smith 2cd05a4c0SHong Zhang /* 3cd05a4c0SHong Zhang Provides utility routines for split MPI communicator. 4cd05a4c0SHong Zhang */ 5c6db04a5SJed Brown #include <petscsys.h> /*I "petscsys.h" I*/ 6cd05a4c0SHong Zhang 7d8a68f86SHong Zhang extern PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm); 8d8a68f86SHong Zhang extern PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm); 9d8a68f86SHong Zhang 10d8a68f86SHong Zhang #undef __FUNCT__ 11d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetNumber" 12d8a68f86SHong Zhang /*@C 13d8a68f86SHong Zhang PetscSubcommSetNumber - Set total number of subcommunicators. 14d8a68f86SHong Zhang 15d8a68f86SHong Zhang Collective on MPI_Comm 16d8a68f86SHong Zhang 17d8a68f86SHong Zhang Input Parameter: 18d8a68f86SHong Zhang + psubcomm - PetscSubcomm context 19d8a68f86SHong Zhang - nsubcomm - the total number of subcommunicators in psubcomm 20d8a68f86SHong Zhang 21d8a68f86SHong Zhang Level: advanced 22d8a68f86SHong Zhang 23d8a68f86SHong Zhang .keywords: communicator 24d8a68f86SHong Zhang 25d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetType(),PetscSubcommSetTypeGeneral() 26d8a68f86SHong Zhang @*/ 277087cfbeSBarry Smith PetscErrorCode PetscSubcommSetNumber(PetscSubcomm psubcomm,PetscInt nsubcomm) 28d8a68f86SHong Zhang { 29d8a68f86SHong Zhang PetscErrorCode ierr; 30d8a68f86SHong Zhang MPI_Comm comm=psubcomm->parent; 31d8a68f86SHong Zhang PetscMPIInt rank,size; 32d8a68f86SHong Zhang 33d8a68f86SHong Zhang PetscFunctionBegin; 34d8a68f86SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate() first"); 35d8a68f86SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 36d8a68f86SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 37d8a68f86SHong 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); 38d8a68f86SHong Zhang 39d8a68f86SHong Zhang psubcomm->n = nsubcomm; 40d8a68f86SHong Zhang PetscFunctionReturn(0); 41d8a68f86SHong Zhang } 42d8a68f86SHong Zhang 43d8a68f86SHong Zhang #undef __FUNCT__ 44d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetType" 45d8a68f86SHong Zhang /*@C 46d8a68f86SHong Zhang PetscSubcommSetType - Set type of subcommunicators. 47d8a68f86SHong Zhang 48d8a68f86SHong Zhang Collective on MPI_Comm 49d8a68f86SHong Zhang 50d8a68f86SHong Zhang Input Parameter: 51d8a68f86SHong Zhang + psubcomm - PetscSubcomm context 521ba920a7SHong Zhang - subcommtype - subcommunicator type, PETSC_SUBCOMM_CONTIGUOUS,PETSC_SUBCOMM_INTERLACED 53d8a68f86SHong Zhang 54d8a68f86SHong Zhang Level: advanced 55d8a68f86SHong Zhang 56d8a68f86SHong Zhang .keywords: communicator 57d8a68f86SHong Zhang 58d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral() 59d8a68f86SHong Zhang @*/ 607087cfbeSBarry Smith PetscErrorCode PetscSubcommSetType(PetscSubcomm psubcomm,const PetscSubcommType subcommtype) 61d8a68f86SHong Zhang { 62d8a68f86SHong Zhang PetscErrorCode ierr; 63d8a68f86SHong Zhang 64d8a68f86SHong Zhang PetscFunctionBegin; 65d8a68f86SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 66d8a68f86SHong Zhang if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n); 67d8a68f86SHong Zhang 68d8a68f86SHong Zhang if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS){ 69d8a68f86SHong Zhang ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr); 70d8a68f86SHong Zhang } else if (subcommtype == PETSC_SUBCOMM_INTERLACED){ 71d8a68f86SHong Zhang ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr); 72d8a68f86SHong Zhang } else { 73d8a68f86SHong Zhang SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",subcommtype); 74d8a68f86SHong Zhang } 75d8a68f86SHong Zhang PetscFunctionReturn(0); 76d8a68f86SHong Zhang } 77d8a68f86SHong Zhang 78d8a68f86SHong Zhang #undef __FUNCT__ 79d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetTypeGeneral" 801ba920a7SHong Zhang /*@C 811ba920a7SHong Zhang PetscSubcommSetTypeGeneral - Set type of subcommunicators from user's specifications 821ba920a7SHong Zhang 831ba920a7SHong Zhang Collective on MPI_Comm 841ba920a7SHong Zhang 851ba920a7SHong Zhang Input Parameter: 861ba920a7SHong Zhang + psubcomm - PetscSubcomm context 871ba920a7SHong Zhang . color - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator. 881ba920a7SHong Zhang . subrank - rank in the subcommunicator 891ba920a7SHong Zhang - duprank - rank in the dupparent (see PetscSubcomm) 901ba920a7SHong Zhang 911ba920a7SHong Zhang Level: advanced 921ba920a7SHong Zhang 931ba920a7SHong Zhang .keywords: communicator, create 941ba920a7SHong Zhang 951ba920a7SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType() 961ba920a7SHong Zhang @*/ 977087cfbeSBarry Smith PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank,PetscMPIInt duprank) 98d8a68f86SHong Zhang { 991ba920a7SHong Zhang PetscErrorCode ierr; 1001ba920a7SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 1011ba920a7SHong Zhang PetscMPIInt size; 1021ba920a7SHong Zhang 103d8a68f86SHong Zhang PetscFunctionBegin; 1041ba920a7SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 1051ba920a7SHong Zhang if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n); 1061ba920a7SHong Zhang 1071ba920a7SHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 1081ba920a7SHong Zhang 1091ba920a7SHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm 1101ba920a7SHong Zhang if duprank is not a valid number, then dupcomm is not created - not all applications require dupcomm! */ 1111ba920a7SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1121ba920a7SHong Zhang if (duprank == PETSC_DECIDE){ 1131ba920a7SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"duprank==PETSC_DECIDE is not supported yet"); 1141ba920a7SHong Zhang } else if (duprank >= 0 && duprank < size){ 1151ba920a7SHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 1161ba920a7SHong Zhang } 117aa9c1079SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr); 118aa9c1079SBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr); 119*b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 120*b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 1211ba920a7SHong Zhang psubcomm->color = color; 122d8a68f86SHong Zhang PetscFunctionReturn(0); 123d8a68f86SHong Zhang } 124638faf0bSHong Zhang 125cd05a4c0SHong Zhang #undef __FUNCT__ 126cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy" 1276bf464f9SBarry Smith PetscErrorCode PetscSubcommDestroy(PetscSubcomm *psubcomm) 128cd05a4c0SHong Zhang { 129cd05a4c0SHong Zhang PetscErrorCode ierr; 130cd05a4c0SHong Zhang 131cd05a4c0SHong Zhang PetscFunctionBegin; 1326bf464f9SBarry Smith if (!*psubcomm) PetscFunctionReturn(0); 133aa9c1079SBarry Smith ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr); 134aa9c1079SBarry Smith ierr = PetscCommDestroy(&(*psubcomm)->comm);CHKERRQ(ierr); 1356bf464f9SBarry Smith ierr = PetscFree((*psubcomm));CHKERRQ(ierr); 136cd05a4c0SHong Zhang PetscFunctionReturn(0); 137cd05a4c0SHong Zhang } 138cd05a4c0SHong Zhang 139cd05a4c0SHong Zhang #undef __FUNCT__ 140cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate" 141ab8c242fSMatthew Knepley /*@C 142cd05a4c0SHong Zhang PetscSubcommCreate - Create a PetscSubcomm context. 143cd05a4c0SHong Zhang 144cd05a4c0SHong Zhang Collective on MPI_Comm 145cd05a4c0SHong Zhang 146cd05a4c0SHong Zhang Input Parameter: 147cd05a4c0SHong Zhang + comm - MPI communicator 148638faf0bSHong Zhang . nsubcomm - the number of subcommunicators to be created 149638faf0bSHong Zhang - subcommtype - subcommunicator type 150cd05a4c0SHong Zhang 151cd05a4c0SHong Zhang Output Parameter: 152cd05a4c0SHong Zhang . psubcomm - location to store the PetscSubcomm context 153cd05a4c0SHong Zhang 154638faf0bSHong Zhang Level: advanced 155cd05a4c0SHong Zhang 156638faf0bSHong Zhang .keywords: communicator, create 157638faf0bSHong Zhang 158638faf0bSHong Zhang .seealso: PetscSubcommDestroy() 159638faf0bSHong Zhang @*/ 1607087cfbeSBarry Smith PetscErrorCode PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm) 161638faf0bSHong Zhang { 162638faf0bSHong Zhang PetscErrorCode ierr; 163638faf0bSHong Zhang 164638faf0bSHong Zhang PetscFunctionBegin; 165638faf0bSHong Zhang 166d8a68f86SHong Zhang ierr = PetscNew(struct _n_PetscSubcomm,psubcomm);CHKERRQ(ierr); 167d8a68f86SHong Zhang (*psubcomm)->parent = comm; 168638faf0bSHong Zhang PetscFunctionReturn(0); 169638faf0bSHong Zhang } 170638faf0bSHong Zhang 171638faf0bSHong Zhang #undef __FUNCT__ 17253c77d0aSJed Brown #define __FUNCT__ "PetscSubcommCreate_contiguous" 173d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm) 174638faf0bSHong Zhang { 175638faf0bSHong Zhang PetscErrorCode ierr; 176d6037b41SHong Zhang PetscMPIInt rank,size,*subsize,duprank=-1,subrank=-1; 177d8a68f86SHong Zhang PetscInt np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n; 178d8a68f86SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 179638faf0bSHong Zhang 180638faf0bSHong Zhang PetscFunctionBegin; 18155e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 18255e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 18355e3b8d2SHong Zhang 184638faf0bSHong Zhang /* get size of each subcommunicator */ 185638faf0bSHong Zhang ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr); 186638faf0bSHong Zhang np_subcomm = size/nsubcomm; 187638faf0bSHong Zhang nleftover = size - nsubcomm*np_subcomm; 188638faf0bSHong Zhang for (i=0; i<nsubcomm; i++){ 189638faf0bSHong Zhang subsize[i] = np_subcomm; 190638faf0bSHong Zhang if (i<nleftover) subsize[i]++; 191638faf0bSHong Zhang } 192638faf0bSHong Zhang 193638faf0bSHong Zhang /* get color and subrank of this proc */ 194638faf0bSHong Zhang rankstart = 0; 195638faf0bSHong Zhang for (i=0; i<nsubcomm; i++){ 196638faf0bSHong Zhang if ( rank >= rankstart && rank < rankstart+subsize[i]) { 197638faf0bSHong Zhang color = i; 198638faf0bSHong Zhang subrank = rank - rankstart; 199638faf0bSHong Zhang duprank = rank; 200638faf0bSHong Zhang break; 201638faf0bSHong Zhang } else { 202638faf0bSHong Zhang rankstart += subsize[i]; 203638faf0bSHong Zhang } 204638faf0bSHong Zhang } 205638faf0bSHong Zhang ierr = PetscFree(subsize);CHKERRQ(ierr); 206638faf0bSHong Zhang 207638faf0bSHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 208638faf0bSHong Zhang 209638faf0bSHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 210638faf0bSHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 211638faf0bSHong Zhang 212aa9c1079SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr); 213aa9c1079SBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr); 214*b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 215*b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 216d8a68f86SHong Zhang psubcomm->color = color; 217638faf0bSHong Zhang PetscFunctionReturn(0); 218638faf0bSHong Zhang } 219638faf0bSHong Zhang 220638faf0bSHong Zhang #undef __FUNCT__ 221638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced" 222638faf0bSHong Zhang /* 223638faf0bSHong Zhang Note: 224638faf0bSHong Zhang In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators 22545fc02eaSBarry Smith by iteratively taking a process into a subcommunicator. 226cd05a4c0SHong Zhang Example: size=4, nsubcomm=(*psubcomm)->n=3 227cd05a4c0SHong Zhang comm=(*psubcomm)->parent: 228cd05a4c0SHong Zhang rank: [0] [1] [2] [3] 229cd05a4c0SHong Zhang color: 0 1 2 0 230cd05a4c0SHong Zhang 231cd05a4c0SHong Zhang subcomm=(*psubcomm)->comm: 232cd05a4c0SHong Zhang subrank: [0] [0] [0] [1] 233cd05a4c0SHong Zhang 234cd05a4c0SHong Zhang dupcomm=(*psubcomm)->dupparent: 235cd05a4c0SHong Zhang duprank: [0] [2] [3] [1] 236cd05a4c0SHong Zhang 237cd05a4c0SHong Zhang Here, subcomm[color = 0] has subsize=2, owns process [0] and [3] 238cd05a4c0SHong Zhang subcomm[color = 1] has subsize=1, owns process [1] 239cd05a4c0SHong Zhang subcomm[color = 2] has subsize=1, owns process [2] 240cd05a4c0SHong Zhang dupcomm has same number of processes as comm, and its duprank maps 241cd05a4c0SHong Zhang processes in subcomm contiguously into a 1d array: 242cd05a4c0SHong Zhang duprank: [0] [1] [2] [3] 243cd05a4c0SHong Zhang rank: [0] [3] [1] [2] 244cd05a4c0SHong Zhang subcomm[0] subcomm[1] subcomm[2] 245638faf0bSHong Zhang */ 246cd05a4c0SHong Zhang 247d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm) 248cd05a4c0SHong Zhang { 249cd05a4c0SHong Zhang PetscErrorCode ierr; 250cd05a4c0SHong Zhang PetscMPIInt rank,size,*subsize,duprank,subrank; 251d8a68f86SHong Zhang PetscInt np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n; 252d8a68f86SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 253cd05a4c0SHong Zhang 254cd05a4c0SHong Zhang PetscFunctionBegin; 25555e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 25655e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 25755e3b8d2SHong Zhang 258cd05a4c0SHong Zhang /* get size of each subcommunicator */ 259cd05a4c0SHong Zhang ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr); 260cd05a4c0SHong Zhang np_subcomm = size/nsubcomm; 261cd05a4c0SHong Zhang nleftover = size - nsubcomm*np_subcomm; 262cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++){ 263cd05a4c0SHong Zhang subsize[i] = np_subcomm; 264cd05a4c0SHong Zhang if (i<nleftover) subsize[i]++; 265cd05a4c0SHong Zhang } 266cd05a4c0SHong Zhang 267cd05a4c0SHong Zhang /* find color for this proc */ 268cd05a4c0SHong Zhang color = rank%nsubcomm; 269cd05a4c0SHong Zhang subrank = rank/nsubcomm; 270cd05a4c0SHong Zhang 271cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 272cd05a4c0SHong Zhang 273cd05a4c0SHong Zhang j = 0; duprank = 0; 274cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++){ 275cd05a4c0SHong Zhang if (j == color){ 276cd05a4c0SHong Zhang duprank += subrank; 277cd05a4c0SHong Zhang break; 278cd05a4c0SHong Zhang } 279cd05a4c0SHong Zhang duprank += subsize[i]; j++; 280cd05a4c0SHong Zhang } 28145fc02eaSBarry Smith ierr = PetscFree(subsize);CHKERRQ(ierr); 282cd05a4c0SHong Zhang 283cd05a4c0SHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 284cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 285cd05a4c0SHong Zhang 286aa9c1079SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr); 287aa9c1079SBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr); 288*b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 289*b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 290d8a68f86SHong Zhang psubcomm->color = color; 291cd05a4c0SHong Zhang PetscFunctionReturn(0); 292cd05a4c0SHong Zhang } 293638faf0bSHong Zhang 294638faf0bSHong Zhang 295