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 @*/ 60*7c764164SBarry Smith PetscErrorCode PetscSubcommSetType(PetscSubcomm psubcomm,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); 72*7c764164SBarry Smith } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",subcommtype); 73d8a68f86SHong Zhang PetscFunctionReturn(0); 74d8a68f86SHong Zhang } 75d8a68f86SHong Zhang 76d8a68f86SHong Zhang #undef __FUNCT__ 77d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetTypeGeneral" 781ba920a7SHong Zhang /*@C 791ba920a7SHong Zhang PetscSubcommSetTypeGeneral - Set type of subcommunicators from user's specifications 801ba920a7SHong Zhang 811ba920a7SHong Zhang Collective on MPI_Comm 821ba920a7SHong Zhang 831ba920a7SHong Zhang Input Parameter: 841ba920a7SHong Zhang + psubcomm - PetscSubcomm context 851ba920a7SHong Zhang . color - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator. 861ba920a7SHong Zhang . subrank - rank in the subcommunicator 871ba920a7SHong Zhang - duprank - rank in the dupparent (see PetscSubcomm) 881ba920a7SHong Zhang 891ba920a7SHong Zhang Level: advanced 901ba920a7SHong Zhang 911ba920a7SHong Zhang .keywords: communicator, create 921ba920a7SHong Zhang 931ba920a7SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType() 941ba920a7SHong Zhang @*/ 957087cfbeSBarry Smith PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank,PetscMPIInt duprank) 96d8a68f86SHong Zhang { 971ba920a7SHong Zhang PetscErrorCode ierr; 981ba920a7SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 991ba920a7SHong Zhang PetscMPIInt size; 1001ba920a7SHong Zhang 101d8a68f86SHong Zhang PetscFunctionBegin; 1021ba920a7SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 1031ba920a7SHong Zhang if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n); 1041ba920a7SHong Zhang 1051ba920a7SHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 1061ba920a7SHong Zhang 1071ba920a7SHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm 1081ba920a7SHong Zhang if duprank is not a valid number, then dupcomm is not created - not all applications require dupcomm! */ 1091ba920a7SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 110*7c764164SBarry Smith if (duprank == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"duprank==PETSC_DECIDE is not supported yet"); 111*7c764164SBarry Smith else if (duprank >= 0 && duprank < size){ 1121ba920a7SHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 1131ba920a7SHong Zhang } 114aa9c1079SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr); 115aa9c1079SBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr); 116b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 117b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 1181ba920a7SHong Zhang psubcomm->color = color; 119d8a68f86SHong Zhang PetscFunctionReturn(0); 120d8a68f86SHong Zhang } 121638faf0bSHong Zhang 122cd05a4c0SHong Zhang #undef __FUNCT__ 123cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy" 1246bf464f9SBarry Smith PetscErrorCode PetscSubcommDestroy(PetscSubcomm *psubcomm) 125cd05a4c0SHong Zhang { 126cd05a4c0SHong Zhang PetscErrorCode ierr; 127cd05a4c0SHong Zhang 128cd05a4c0SHong Zhang PetscFunctionBegin; 1296bf464f9SBarry Smith if (!*psubcomm) PetscFunctionReturn(0); 130aa9c1079SBarry Smith ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr); 131aa9c1079SBarry Smith ierr = PetscCommDestroy(&(*psubcomm)->comm);CHKERRQ(ierr); 1326bf464f9SBarry Smith ierr = PetscFree((*psubcomm));CHKERRQ(ierr); 133cd05a4c0SHong Zhang PetscFunctionReturn(0); 134cd05a4c0SHong Zhang } 135cd05a4c0SHong Zhang 136cd05a4c0SHong Zhang #undef __FUNCT__ 137cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate" 138ab8c242fSMatthew Knepley /*@C 139cd05a4c0SHong Zhang PetscSubcommCreate - Create a PetscSubcomm context. 140cd05a4c0SHong Zhang 141cd05a4c0SHong Zhang Collective on MPI_Comm 142cd05a4c0SHong Zhang 143cd05a4c0SHong Zhang Input Parameter: 1449873d53eSJed Brown . comm - MPI communicator 145cd05a4c0SHong Zhang 146cd05a4c0SHong Zhang Output Parameter: 147cd05a4c0SHong Zhang . psubcomm - location to store the PetscSubcomm context 148cd05a4c0SHong Zhang 149638faf0bSHong Zhang Level: advanced 150cd05a4c0SHong Zhang 151638faf0bSHong Zhang .keywords: communicator, create 152638faf0bSHong Zhang 153638faf0bSHong Zhang .seealso: PetscSubcommDestroy() 154638faf0bSHong Zhang @*/ 1557087cfbeSBarry Smith PetscErrorCode PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm) 156638faf0bSHong Zhang { 157638faf0bSHong Zhang PetscErrorCode ierr; 158638faf0bSHong Zhang 159638faf0bSHong Zhang PetscFunctionBegin; 160d8a68f86SHong Zhang ierr = PetscNew(struct _n_PetscSubcomm,psubcomm);CHKERRQ(ierr); 161d8a68f86SHong Zhang (*psubcomm)->parent = comm; 162638faf0bSHong Zhang PetscFunctionReturn(0); 163638faf0bSHong Zhang } 164638faf0bSHong Zhang 165638faf0bSHong Zhang #undef __FUNCT__ 16653c77d0aSJed Brown #define __FUNCT__ "PetscSubcommCreate_contiguous" 167d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm) 168638faf0bSHong Zhang { 169638faf0bSHong Zhang PetscErrorCode ierr; 170d6037b41SHong Zhang PetscMPIInt rank,size,*subsize,duprank=-1,subrank=-1; 171d8a68f86SHong Zhang PetscInt np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n; 172d8a68f86SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 173638faf0bSHong Zhang 174638faf0bSHong Zhang PetscFunctionBegin; 17555e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 17655e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 17755e3b8d2SHong Zhang 178638faf0bSHong Zhang /* get size of each subcommunicator */ 179638faf0bSHong Zhang ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr); 180638faf0bSHong Zhang np_subcomm = size/nsubcomm; 181638faf0bSHong Zhang nleftover = size - nsubcomm*np_subcomm; 182638faf0bSHong Zhang for (i=0; i<nsubcomm; i++){ 183638faf0bSHong Zhang subsize[i] = np_subcomm; 184638faf0bSHong Zhang if (i<nleftover) subsize[i]++; 185638faf0bSHong Zhang } 186638faf0bSHong Zhang 187638faf0bSHong Zhang /* get color and subrank of this proc */ 188638faf0bSHong Zhang rankstart = 0; 189638faf0bSHong Zhang for (i=0; i<nsubcomm; i++){ 190638faf0bSHong Zhang if ( rank >= rankstart && rank < rankstart+subsize[i]) { 191638faf0bSHong Zhang color = i; 192638faf0bSHong Zhang subrank = rank - rankstart; 193638faf0bSHong Zhang duprank = rank; 194638faf0bSHong Zhang break; 195638faf0bSHong Zhang } else { 196638faf0bSHong Zhang rankstart += subsize[i]; 197638faf0bSHong Zhang } 198638faf0bSHong Zhang } 199638faf0bSHong Zhang ierr = PetscFree(subsize);CHKERRQ(ierr); 200638faf0bSHong Zhang 201638faf0bSHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 202638faf0bSHong Zhang 203638faf0bSHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 204638faf0bSHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 205638faf0bSHong Zhang 206aa9c1079SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr); 207aa9c1079SBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr); 208b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 209b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 210d8a68f86SHong Zhang psubcomm->color = color; 211638faf0bSHong Zhang PetscFunctionReturn(0); 212638faf0bSHong Zhang } 213638faf0bSHong Zhang 214638faf0bSHong Zhang #undef __FUNCT__ 215638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced" 216638faf0bSHong Zhang /* 217638faf0bSHong Zhang Note: 218638faf0bSHong Zhang In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators 21945fc02eaSBarry Smith by iteratively taking a process into a subcommunicator. 220cd05a4c0SHong Zhang Example: size=4, nsubcomm=(*psubcomm)->n=3 221cd05a4c0SHong Zhang comm=(*psubcomm)->parent: 222cd05a4c0SHong Zhang rank: [0] [1] [2] [3] 223cd05a4c0SHong Zhang color: 0 1 2 0 224cd05a4c0SHong Zhang 225cd05a4c0SHong Zhang subcomm=(*psubcomm)->comm: 226cd05a4c0SHong Zhang subrank: [0] [0] [0] [1] 227cd05a4c0SHong Zhang 228cd05a4c0SHong Zhang dupcomm=(*psubcomm)->dupparent: 229cd05a4c0SHong Zhang duprank: [0] [2] [3] [1] 230cd05a4c0SHong Zhang 231cd05a4c0SHong Zhang Here, subcomm[color = 0] has subsize=2, owns process [0] and [3] 232cd05a4c0SHong Zhang subcomm[color = 1] has subsize=1, owns process [1] 233cd05a4c0SHong Zhang subcomm[color = 2] has subsize=1, owns process [2] 234cd05a4c0SHong Zhang dupcomm has same number of processes as comm, and its duprank maps 235cd05a4c0SHong Zhang processes in subcomm contiguously into a 1d array: 236cd05a4c0SHong Zhang duprank: [0] [1] [2] [3] 237cd05a4c0SHong Zhang rank: [0] [3] [1] [2] 238cd05a4c0SHong Zhang subcomm[0] subcomm[1] subcomm[2] 239638faf0bSHong Zhang */ 240cd05a4c0SHong Zhang 241d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm) 242cd05a4c0SHong Zhang { 243cd05a4c0SHong Zhang PetscErrorCode ierr; 244cd05a4c0SHong Zhang PetscMPIInt rank,size,*subsize,duprank,subrank; 245d8a68f86SHong Zhang PetscInt np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n; 246d8a68f86SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 247cd05a4c0SHong Zhang 248cd05a4c0SHong Zhang PetscFunctionBegin; 24955e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 25055e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 25155e3b8d2SHong Zhang 252cd05a4c0SHong Zhang /* get size of each subcommunicator */ 253cd05a4c0SHong Zhang ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr); 254cd05a4c0SHong Zhang np_subcomm = size/nsubcomm; 255cd05a4c0SHong Zhang nleftover = size - nsubcomm*np_subcomm; 256cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++){ 257cd05a4c0SHong Zhang subsize[i] = np_subcomm; 258cd05a4c0SHong Zhang if (i<nleftover) subsize[i]++; 259cd05a4c0SHong Zhang } 260cd05a4c0SHong Zhang 261cd05a4c0SHong Zhang /* find color for this proc */ 262cd05a4c0SHong Zhang color = rank%nsubcomm; 263cd05a4c0SHong Zhang subrank = rank/nsubcomm; 264cd05a4c0SHong Zhang 265cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 266cd05a4c0SHong Zhang 267cd05a4c0SHong Zhang j = 0; duprank = 0; 268cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++){ 269cd05a4c0SHong Zhang if (j == color){ 270cd05a4c0SHong Zhang duprank += subrank; 271cd05a4c0SHong Zhang break; 272cd05a4c0SHong Zhang } 273cd05a4c0SHong Zhang duprank += subsize[i]; j++; 274cd05a4c0SHong Zhang } 27545fc02eaSBarry Smith ierr = PetscFree(subsize);CHKERRQ(ierr); 276cd05a4c0SHong Zhang 277cd05a4c0SHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 278cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 279cd05a4c0SHong Zhang 280aa9c1079SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr); 281aa9c1079SBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr); 282b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 283b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 284d8a68f86SHong Zhang psubcomm->color = color; 285cd05a4c0SHong Zhang PetscFunctionReturn(0); 286cd05a4c0SHong Zhang } 287638faf0bSHong Zhang 288638faf0bSHong Zhang 289