17d0a6c19SBarry Smith 2cd05a4c0SHong Zhang /* 3cd05a4c0SHong Zhang Provides utility routines for split MPI communicator. 4cd05a4c0SHong Zhang */ 5c6db04a5SJed Brown #include <petscsys.h> /*I "petscsys.h" I*/ 6422737f6SJed Brown #include <petsc-private/threadcommimpl.h> /* Petsc_ThreadComm_keyval */ 7053d1c95SHong Zhang #include <petscviewer.h> 8cd05a4c0SHong Zhang 96a6fc655SJed Brown const char *const PetscSubcommTypes[] = {"GENERAL","CONTIGUOUS","INTERLACED","PetscSubcommType","PETSC_SUBCOMM_",0}; 106a6fc655SJed Brown 11d8a68f86SHong Zhang extern PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm); 12d8a68f86SHong Zhang extern PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm); 13d8a68f86SHong Zhang 14d8a68f86SHong Zhang #undef __FUNCT__ 15053d1c95SHong Zhang #define __FUNCT__ "PetscSubcommView" 16053d1c95SHong Zhang PetscErrorCode PetscSubcommView(PetscSubcomm psubcomm,PetscViewer viewer) 17053d1c95SHong Zhang { 18053d1c95SHong Zhang PetscErrorCode ierr; 19053d1c95SHong Zhang PetscBool iascii; 20053d1c95SHong Zhang PetscViewerFormat format; 21053d1c95SHong Zhang 22053d1c95SHong Zhang PetscFunctionBegin; 23053d1c95SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 24053d1c95SHong Zhang if (iascii) { 25053d1c95SHong Zhang ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 26053d1c95SHong Zhang if (format == PETSC_VIEWER_DEFAULT) { 27053d1c95SHong Zhang MPI_Comm comm=psubcomm->parent; 28053d1c95SHong Zhang PetscMPIInt rank,size,subsize,subrank,duprank; 29053d1c95SHong Zhang 30053d1c95SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 31053d1c95SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"PetscSubcomm type %s with total %D MPI processes:\n",PetscSubcommTypes[psubcomm->type],size);CHKERRQ(ierr); 32053d1c95SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 33053d1c95SHong Zhang ierr = MPI_Comm_size(psubcomm->comm,&subsize);CHKERRQ(ierr); 34053d1c95SHong Zhang ierr = MPI_Comm_rank(psubcomm->comm,&subrank);CHKERRQ(ierr); 35053d1c95SHong Zhang ierr = MPI_Comm_rank(psubcomm->dupparent,&duprank);CHKERRQ(ierr); 36053d1c95SHong Zhang ierr = PetscSynchronizedPrintf(comm," [%D], color %D, sub-size %D, sub-rank %D, duprank %D\n",rank,psubcomm->color,subsize,subrank,duprank); 37053d1c95SHong Zhang ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 38053d1c95SHong Zhang } 39053d1c95SHong Zhang } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported yet"); 40053d1c95SHong Zhang PetscFunctionReturn(0); 41053d1c95SHong Zhang } 42053d1c95SHong Zhang 43053d1c95SHong Zhang #undef __FUNCT__ 44d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetNumber" 45d8a68f86SHong Zhang /*@C 46d8a68f86SHong Zhang PetscSubcommSetNumber - Set total number of subcommunicators. 47d8a68f86SHong Zhang 48d8a68f86SHong Zhang Collective on MPI_Comm 49d8a68f86SHong Zhang 50d8a68f86SHong Zhang Input Parameter: 51d8a68f86SHong Zhang + psubcomm - PetscSubcomm context 52d8a68f86SHong Zhang - nsubcomm - the total number of subcommunicators in psubcomm 53d8a68f86SHong Zhang 54d8a68f86SHong Zhang Level: advanced 55d8a68f86SHong Zhang 56d8a68f86SHong Zhang .keywords: communicator 57d8a68f86SHong Zhang 58d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetType(),PetscSubcommSetTypeGeneral() 59d8a68f86SHong Zhang @*/ 607087cfbeSBarry Smith PetscErrorCode PetscSubcommSetNumber(PetscSubcomm psubcomm,PetscInt nsubcomm) 61d8a68f86SHong Zhang { 62d8a68f86SHong Zhang PetscErrorCode ierr; 63d8a68f86SHong Zhang MPI_Comm comm=psubcomm->parent; 64d8a68f86SHong Zhang PetscMPIInt rank,size; 65d8a68f86SHong Zhang 66d8a68f86SHong Zhang PetscFunctionBegin; 67d8a68f86SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate() first"); 68d8a68f86SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 69d8a68f86SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 70d8a68f86SHong 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); 71d8a68f86SHong Zhang 72d8a68f86SHong Zhang psubcomm->n = nsubcomm; 73d8a68f86SHong Zhang PetscFunctionReturn(0); 74d8a68f86SHong Zhang } 75d8a68f86SHong Zhang 76d8a68f86SHong Zhang #undef __FUNCT__ 77d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetType" 78d8a68f86SHong Zhang /*@C 79d8a68f86SHong Zhang PetscSubcommSetType - Set type of subcommunicators. 80d8a68f86SHong Zhang 81d8a68f86SHong Zhang Collective on MPI_Comm 82d8a68f86SHong Zhang 83d8a68f86SHong Zhang Input Parameter: 84d8a68f86SHong Zhang + psubcomm - PetscSubcomm context 851ba920a7SHong Zhang - subcommtype - subcommunicator type, PETSC_SUBCOMM_CONTIGUOUS,PETSC_SUBCOMM_INTERLACED 86d8a68f86SHong Zhang 87d8a68f86SHong Zhang Level: advanced 88d8a68f86SHong Zhang 89d8a68f86SHong Zhang .keywords: communicator 90d8a68f86SHong Zhang 91d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral() 92d8a68f86SHong Zhang @*/ 937c764164SBarry Smith PetscErrorCode PetscSubcommSetType(PetscSubcomm psubcomm,PetscSubcommType subcommtype) 94d8a68f86SHong Zhang { 95d8a68f86SHong Zhang PetscErrorCode ierr; 96d8a68f86SHong Zhang 97d8a68f86SHong Zhang PetscFunctionBegin; 98d8a68f86SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 99d8a68f86SHong Zhang if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n); 100d8a68f86SHong Zhang 101d8a68f86SHong Zhang if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS) { 102d8a68f86SHong Zhang ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr); 103d8a68f86SHong Zhang } else if (subcommtype == PETSC_SUBCOMM_INTERLACED) { 104d8a68f86SHong Zhang ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr); 1057c764164SBarry Smith } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",subcommtype); 106d8a68f86SHong Zhang PetscFunctionReturn(0); 107d8a68f86SHong Zhang } 108d8a68f86SHong Zhang 109d8a68f86SHong Zhang #undef __FUNCT__ 110d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetTypeGeneral" 1111ba920a7SHong Zhang /*@C 1121ba920a7SHong Zhang PetscSubcommSetTypeGeneral - Set type of subcommunicators from user's specifications 1131ba920a7SHong Zhang 1141ba920a7SHong Zhang Collective on MPI_Comm 1151ba920a7SHong Zhang 1161ba920a7SHong Zhang Input Parameter: 1171ba920a7SHong Zhang + psubcomm - PetscSubcomm context 1181ba920a7SHong Zhang . color - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator. 1191ba920a7SHong Zhang . subrank - rank in the subcommunicator 1201ba920a7SHong Zhang - duprank - rank in the dupparent (see PetscSubcomm) 1211ba920a7SHong Zhang 1221ba920a7SHong Zhang Level: advanced 1231ba920a7SHong Zhang 1241ba920a7SHong Zhang .keywords: communicator, create 1251ba920a7SHong Zhang 1261ba920a7SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType() 1271ba920a7SHong Zhang @*/ 1287087cfbeSBarry Smith PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank,PetscMPIInt duprank) 129d8a68f86SHong Zhang { 1301ba920a7SHong Zhang PetscErrorCode ierr; 1311ba920a7SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 1321ba920a7SHong Zhang PetscMPIInt size; 1331ba920a7SHong Zhang 134d8a68f86SHong Zhang PetscFunctionBegin; 1351ba920a7SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 1361ba920a7SHong Zhang if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n); 1371ba920a7SHong Zhang 1381ba920a7SHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 1391ba920a7SHong Zhang 1401ba920a7SHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm 1411ba920a7SHong Zhang if duprank is not a valid number, then dupcomm is not created - not all applications require dupcomm! */ 1421ba920a7SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1437c764164SBarry Smith if (duprank == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"duprank==PETSC_DECIDE is not supported yet"); 1447c764164SBarry Smith else if (duprank >= 0 && duprank < size) { 1451ba920a7SHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 1461ba920a7SHong Zhang } 1470298fd71SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 1480298fd71SBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr); 149b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 150b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 151a297a907SKarl Rupp 1521ba920a7SHong Zhang psubcomm->color = color; 153d8a68f86SHong Zhang PetscFunctionReturn(0); 154d8a68f86SHong Zhang } 155638faf0bSHong Zhang 156cd05a4c0SHong Zhang #undef __FUNCT__ 157cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy" 1586bf464f9SBarry Smith PetscErrorCode PetscSubcommDestroy(PetscSubcomm *psubcomm) 159cd05a4c0SHong Zhang { 160cd05a4c0SHong Zhang PetscErrorCode ierr; 161cd05a4c0SHong Zhang 162cd05a4c0SHong Zhang PetscFunctionBegin; 1636bf464f9SBarry Smith if (!*psubcomm) PetscFunctionReturn(0); 164aa9c1079SBarry Smith ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr); 165aa9c1079SBarry Smith ierr = PetscCommDestroy(&(*psubcomm)->comm);CHKERRQ(ierr); 166*e37c6257SHong Zhang ierr = PetscFree((*psubcomm)->subsize);CHKERRQ(ierr); 1676bf464f9SBarry Smith ierr = PetscFree((*psubcomm));CHKERRQ(ierr); 168cd05a4c0SHong Zhang PetscFunctionReturn(0); 169cd05a4c0SHong Zhang } 170cd05a4c0SHong Zhang 171cd05a4c0SHong Zhang #undef __FUNCT__ 172cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate" 173ab8c242fSMatthew Knepley /*@C 174cd05a4c0SHong Zhang PetscSubcommCreate - Create a PetscSubcomm context. 175cd05a4c0SHong Zhang 176cd05a4c0SHong Zhang Collective on MPI_Comm 177cd05a4c0SHong Zhang 178cd05a4c0SHong Zhang Input Parameter: 1799873d53eSJed Brown . comm - MPI communicator 180cd05a4c0SHong Zhang 181cd05a4c0SHong Zhang Output Parameter: 182cd05a4c0SHong Zhang . psubcomm - location to store the PetscSubcomm context 183cd05a4c0SHong Zhang 184638faf0bSHong Zhang Level: advanced 185cd05a4c0SHong Zhang 186638faf0bSHong Zhang .keywords: communicator, create 187638faf0bSHong Zhang 188638faf0bSHong Zhang .seealso: PetscSubcommDestroy() 189638faf0bSHong Zhang @*/ 1907087cfbeSBarry Smith PetscErrorCode PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm) 191638faf0bSHong Zhang { 192638faf0bSHong Zhang PetscErrorCode ierr; 193d3b23db5SHong Zhang PetscMPIInt rank,size; 194638faf0bSHong Zhang 195638faf0bSHong Zhang PetscFunctionBegin; 196d8a68f86SHong Zhang ierr = PetscNew(struct _n_PetscSubcomm,psubcomm);CHKERRQ(ierr); 197a297a907SKarl Rupp 198d3b23db5SHong Zhang /* set defaults */ 199d3b23db5SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 200d3b23db5SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 201d8a68f86SHong Zhang (*psubcomm)->parent = comm; 202d3b23db5SHong Zhang (*psubcomm)->dupparent = comm; 203d3b23db5SHong Zhang (*psubcomm)->comm = PETSC_COMM_SELF; 204d3b23db5SHong Zhang (*psubcomm)->n = size; 205d3b23db5SHong Zhang (*psubcomm)->color = rank; 206*e37c6257SHong Zhang (*psubcomm)->subsize = NULL; 207d3b23db5SHong Zhang (*psubcomm)->type = PETSC_SUBCOMM_INTERLACED; 208638faf0bSHong Zhang PetscFunctionReturn(0); 209638faf0bSHong Zhang } 210638faf0bSHong Zhang 211638faf0bSHong Zhang #undef __FUNCT__ 21253c77d0aSJed Brown #define __FUNCT__ "PetscSubcommCreate_contiguous" 213d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm) 214638faf0bSHong Zhang { 215638faf0bSHong Zhang PetscErrorCode ierr; 216d6037b41SHong Zhang PetscMPIInt rank,size,*subsize,duprank=-1,subrank=-1; 217d8a68f86SHong Zhang PetscInt np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n; 218d8a68f86SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 219638faf0bSHong Zhang 220638faf0bSHong Zhang PetscFunctionBegin; 22155e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 22255e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 22355e3b8d2SHong Zhang 224638faf0bSHong Zhang /* get size of each subcommunicator */ 225638faf0bSHong Zhang ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr); 226a297a907SKarl Rupp 227638faf0bSHong Zhang np_subcomm = size/nsubcomm; 228638faf0bSHong Zhang nleftover = size - nsubcomm*np_subcomm; 229638faf0bSHong Zhang for (i=0; i<nsubcomm; i++) { 230638faf0bSHong Zhang subsize[i] = np_subcomm; 231638faf0bSHong Zhang if (i<nleftover) subsize[i]++; 232638faf0bSHong Zhang } 233638faf0bSHong Zhang 234638faf0bSHong Zhang /* get color and subrank of this proc */ 235638faf0bSHong Zhang rankstart = 0; 236638faf0bSHong Zhang for (i=0; i<nsubcomm; i++) { 237638faf0bSHong Zhang if (rank >= rankstart && rank < rankstart+subsize[i]) { 238638faf0bSHong Zhang color = i; 239638faf0bSHong Zhang subrank = rank - rankstart; 240638faf0bSHong Zhang duprank = rank; 241638faf0bSHong Zhang break; 242a297a907SKarl Rupp } else rankstart += subsize[i]; 243638faf0bSHong Zhang } 244638faf0bSHong Zhang 245638faf0bSHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 246638faf0bSHong Zhang 247638faf0bSHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 248638faf0bSHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 24952af3f7aSShri Abhyankar { 25052af3f7aSShri Abhyankar PetscThreadComm tcomm; 25152af3f7aSShri Abhyankar ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr); 25252af3f7aSShri Abhyankar ierr = MPI_Attr_put(dupcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr); 25352af3f7aSShri Abhyankar tcomm->refct++; 25452af3f7aSShri Abhyankar ierr = MPI_Attr_put(subcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr); 25552af3f7aSShri Abhyankar tcomm->refct++; 25652af3f7aSShri Abhyankar } 2570298fd71SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 2580298fd71SBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr); 259b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 260b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 261a297a907SKarl Rupp 262d8a68f86SHong Zhang psubcomm->color = color; 263*e37c6257SHong Zhang psubcomm->subsize = subsize; 264f38d543fSHong Zhang psubcomm->type = PETSC_SUBCOMM_CONTIGUOUS; 265638faf0bSHong Zhang PetscFunctionReturn(0); 266638faf0bSHong Zhang } 267638faf0bSHong Zhang 268638faf0bSHong Zhang #undef __FUNCT__ 269638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced" 270638faf0bSHong Zhang /* 271638faf0bSHong Zhang Note: 272638faf0bSHong Zhang In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators 27345fc02eaSBarry Smith by iteratively taking a process into a subcommunicator. 274cd05a4c0SHong Zhang Example: size=4, nsubcomm=(*psubcomm)->n=3 275cd05a4c0SHong Zhang comm=(*psubcomm)->parent: 276cd05a4c0SHong Zhang rank: [0] [1] [2] [3] 277cd05a4c0SHong Zhang color: 0 1 2 0 278cd05a4c0SHong Zhang 279cd05a4c0SHong Zhang subcomm=(*psubcomm)->comm: 280cd05a4c0SHong Zhang subrank: [0] [0] [0] [1] 281cd05a4c0SHong Zhang 282cd05a4c0SHong Zhang dupcomm=(*psubcomm)->dupparent: 283cd05a4c0SHong Zhang duprank: [0] [2] [3] [1] 284cd05a4c0SHong Zhang 285cd05a4c0SHong Zhang Here, subcomm[color = 0] has subsize=2, owns process [0] and [3] 286cd05a4c0SHong Zhang subcomm[color = 1] has subsize=1, owns process [1] 287cd05a4c0SHong Zhang subcomm[color = 2] has subsize=1, owns process [2] 288cd05a4c0SHong Zhang dupcomm has same number of processes as comm, and its duprank maps 289cd05a4c0SHong Zhang processes in subcomm contiguously into a 1d array: 290cd05a4c0SHong Zhang duprank: [0] [1] [2] [3] 291cd05a4c0SHong Zhang rank: [0] [3] [1] [2] 292cd05a4c0SHong Zhang subcomm[0] subcomm[1] subcomm[2] 293638faf0bSHong Zhang */ 294cd05a4c0SHong Zhang 295d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm) 296cd05a4c0SHong Zhang { 297cd05a4c0SHong Zhang PetscErrorCode ierr; 298cd05a4c0SHong Zhang PetscMPIInt rank,size,*subsize,duprank,subrank; 299d8a68f86SHong Zhang PetscInt np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n; 300d8a68f86SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 301cd05a4c0SHong Zhang 302cd05a4c0SHong Zhang PetscFunctionBegin; 30355e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 30455e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 30555e3b8d2SHong Zhang 306cd05a4c0SHong Zhang /* get size of each subcommunicator */ 307cd05a4c0SHong Zhang ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr); 308a297a907SKarl Rupp 309cd05a4c0SHong Zhang np_subcomm = size/nsubcomm; 310cd05a4c0SHong Zhang nleftover = size - nsubcomm*np_subcomm; 311cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++) { 312cd05a4c0SHong Zhang subsize[i] = np_subcomm; 313cd05a4c0SHong Zhang if (i<nleftover) subsize[i]++; 314cd05a4c0SHong Zhang } 315cd05a4c0SHong Zhang 316cd05a4c0SHong Zhang /* find color for this proc */ 317cd05a4c0SHong Zhang color = rank%nsubcomm; 318cd05a4c0SHong Zhang subrank = rank/nsubcomm; 319cd05a4c0SHong Zhang 320cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 321cd05a4c0SHong Zhang 322cd05a4c0SHong Zhang j = 0; duprank = 0; 323cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++) { 324cd05a4c0SHong Zhang if (j == color) { 325cd05a4c0SHong Zhang duprank += subrank; 326cd05a4c0SHong Zhang break; 327cd05a4c0SHong Zhang } 328cd05a4c0SHong Zhang duprank += subsize[i]; j++; 329cd05a4c0SHong Zhang } 330cd05a4c0SHong Zhang 331cd05a4c0SHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 332cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 33352af3f7aSShri Abhyankar { 33452af3f7aSShri Abhyankar PetscThreadComm tcomm; 33552af3f7aSShri Abhyankar ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr); 33652af3f7aSShri Abhyankar ierr = MPI_Attr_put(dupcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr); 33752af3f7aSShri Abhyankar tcomm->refct++; 33852af3f7aSShri Abhyankar ierr = MPI_Attr_put(subcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr); 33952af3f7aSShri Abhyankar tcomm->refct++; 34052af3f7aSShri Abhyankar } 3410298fd71SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 3420298fd71SBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr); 343b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 344b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 345a297a907SKarl Rupp 346d8a68f86SHong Zhang psubcomm->color = color; 347*e37c6257SHong Zhang psubcomm->subsize = subsize; 348f38d543fSHong Zhang psubcomm->type = PETSC_SUBCOMM_INTERLACED; 349cd05a4c0SHong Zhang PetscFunctionReturn(0); 350cd05a4c0SHong Zhang } 351638faf0bSHong Zhang 352638faf0bSHong Zhang 353