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); 13*f68be91cSHong Zhang #undef __FUNCT__ 14*f68be91cSHong Zhang #define __FUNCT__ "PetscSubcommSetFromOptions" 15*f68be91cSHong Zhang PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm psubcomm) 16*f68be91cSHong Zhang { 17*f68be91cSHong Zhang PetscErrorCode ierr; 18*f68be91cSHong Zhang PetscInt type; 19*f68be91cSHong Zhang const char *psubcommTypes[3] = {"general","contiguous","interlaced"}; 20*f68be91cSHong Zhang PetscBool flg; 21*f68be91cSHong Zhang 22*f68be91cSHong Zhang PetscFunctionBegin; 23*f68be91cSHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Must call PetscSubcommCreate firt"); 24*f68be91cSHong Zhang ierr = PetscOptionsEList("-psubcomm_type","PETSc subcommunicator","PetscSubcommSetType",psubcommTypes,3,psubcommTypes[2],&type,&flg);CHKERRQ(ierr); 25*f68be91cSHong Zhang if (flg && psubcomm->type != type) { 26*f68be91cSHong Zhang /* free old structures */ 27*f68be91cSHong Zhang ierr = PetscCommDestroy(&(psubcomm)->dupparent);CHKERRQ(ierr); 28*f68be91cSHong Zhang ierr = PetscCommDestroy(&(psubcomm)->comm);CHKERRQ(ierr); 29*f68be91cSHong Zhang ierr = PetscFree((psubcomm)->subsize);CHKERRQ(ierr); 30*f68be91cSHong Zhang switch (type) { 31*f68be91cSHong Zhang case 1: 32*f68be91cSHong Zhang ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr); 33*f68be91cSHong Zhang break; 34*f68be91cSHong Zhang case 2: 35*f68be91cSHong Zhang ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr); 36*f68be91cSHong Zhang break; 37*f68be91cSHong Zhang default: 38*f68be91cSHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",type); 39*f68be91cSHong Zhang } 40*f68be91cSHong Zhang } 41*f68be91cSHong Zhang PetscFunctionReturn(0); 42*f68be91cSHong Zhang } 43d8a68f86SHong Zhang 44d8a68f86SHong Zhang #undef __FUNCT__ 45053d1c95SHong Zhang #define __FUNCT__ "PetscSubcommView" 46053d1c95SHong Zhang PetscErrorCode PetscSubcommView(PetscSubcomm psubcomm,PetscViewer viewer) 47053d1c95SHong Zhang { 48053d1c95SHong Zhang PetscErrorCode ierr; 49053d1c95SHong Zhang PetscBool iascii; 50053d1c95SHong Zhang PetscViewerFormat format; 51053d1c95SHong Zhang 52053d1c95SHong Zhang PetscFunctionBegin; 53053d1c95SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 54053d1c95SHong Zhang if (iascii) { 55053d1c95SHong Zhang ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 56053d1c95SHong Zhang if (format == PETSC_VIEWER_DEFAULT) { 57053d1c95SHong Zhang MPI_Comm comm=psubcomm->parent; 58053d1c95SHong Zhang PetscMPIInt rank,size,subsize,subrank,duprank; 59053d1c95SHong Zhang 60053d1c95SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 61053d1c95SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"PetscSubcomm type %s with total %D MPI processes:\n",PetscSubcommTypes[psubcomm->type],size);CHKERRQ(ierr); 62053d1c95SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 63053d1c95SHong Zhang ierr = MPI_Comm_size(psubcomm->comm,&subsize);CHKERRQ(ierr); 64053d1c95SHong Zhang ierr = MPI_Comm_rank(psubcomm->comm,&subrank);CHKERRQ(ierr); 65053d1c95SHong Zhang ierr = MPI_Comm_rank(psubcomm->dupparent,&duprank);CHKERRQ(ierr); 66053d1c95SHong Zhang ierr = PetscSynchronizedPrintf(comm," [%D], color %D, sub-size %D, sub-rank %D, duprank %D\n",rank,psubcomm->color,subsize,subrank,duprank); 67053d1c95SHong Zhang ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 68053d1c95SHong Zhang } 69053d1c95SHong Zhang } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported yet"); 70053d1c95SHong Zhang PetscFunctionReturn(0); 71053d1c95SHong Zhang } 72053d1c95SHong Zhang 73053d1c95SHong Zhang #undef __FUNCT__ 74d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetNumber" 75d8a68f86SHong Zhang /*@C 76d8a68f86SHong Zhang PetscSubcommSetNumber - Set total number of subcommunicators. 77d8a68f86SHong Zhang 78d8a68f86SHong Zhang Collective on MPI_Comm 79d8a68f86SHong Zhang 80d8a68f86SHong Zhang Input Parameter: 81d8a68f86SHong Zhang + psubcomm - PetscSubcomm context 82d8a68f86SHong Zhang - nsubcomm - the total number of subcommunicators in psubcomm 83d8a68f86SHong Zhang 84d8a68f86SHong Zhang Level: advanced 85d8a68f86SHong Zhang 86d8a68f86SHong Zhang .keywords: communicator 87d8a68f86SHong Zhang 88d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetType(),PetscSubcommSetTypeGeneral() 89d8a68f86SHong Zhang @*/ 907087cfbeSBarry Smith PetscErrorCode PetscSubcommSetNumber(PetscSubcomm psubcomm,PetscInt nsubcomm) 91d8a68f86SHong Zhang { 92d8a68f86SHong Zhang PetscErrorCode ierr; 93d8a68f86SHong Zhang MPI_Comm comm=psubcomm->parent; 94d8a68f86SHong Zhang PetscMPIInt rank,size; 95d8a68f86SHong Zhang 96d8a68f86SHong Zhang PetscFunctionBegin; 97d8a68f86SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate() first"); 98d8a68f86SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 99d8a68f86SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 100d8a68f86SHong 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); 101d8a68f86SHong Zhang 102d8a68f86SHong Zhang psubcomm->n = nsubcomm; 103d8a68f86SHong Zhang PetscFunctionReturn(0); 104d8a68f86SHong Zhang } 105d8a68f86SHong Zhang 106d8a68f86SHong Zhang #undef __FUNCT__ 107d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetType" 108d8a68f86SHong Zhang /*@C 109d8a68f86SHong Zhang PetscSubcommSetType - Set type of subcommunicators. 110d8a68f86SHong Zhang 111d8a68f86SHong Zhang Collective on MPI_Comm 112d8a68f86SHong Zhang 113d8a68f86SHong Zhang Input Parameter: 114d8a68f86SHong Zhang + psubcomm - PetscSubcomm context 1151ba920a7SHong Zhang - subcommtype - subcommunicator type, PETSC_SUBCOMM_CONTIGUOUS,PETSC_SUBCOMM_INTERLACED 116d8a68f86SHong Zhang 117d8a68f86SHong Zhang Level: advanced 118d8a68f86SHong Zhang 119d8a68f86SHong Zhang .keywords: communicator 120d8a68f86SHong Zhang 121d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral() 122d8a68f86SHong Zhang @*/ 1237c764164SBarry Smith PetscErrorCode PetscSubcommSetType(PetscSubcomm psubcomm,PetscSubcommType subcommtype) 124d8a68f86SHong Zhang { 125d8a68f86SHong Zhang PetscErrorCode ierr; 126d8a68f86SHong Zhang 127d8a68f86SHong Zhang PetscFunctionBegin; 128d8a68f86SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 129d8a68f86SHong Zhang if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n); 130d8a68f86SHong Zhang 131d8a68f86SHong Zhang if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS) { 132d8a68f86SHong Zhang ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr); 133d8a68f86SHong Zhang } else if (subcommtype == PETSC_SUBCOMM_INTERLACED) { 134d8a68f86SHong Zhang ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr); 1357c764164SBarry Smith } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",subcommtype); 136d8a68f86SHong Zhang PetscFunctionReturn(0); 137d8a68f86SHong Zhang } 138d8a68f86SHong Zhang 139d8a68f86SHong Zhang #undef __FUNCT__ 140d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetTypeGeneral" 1411ba920a7SHong Zhang /*@C 1421ba920a7SHong Zhang PetscSubcommSetTypeGeneral - Set type of subcommunicators from user's specifications 1431ba920a7SHong Zhang 1441ba920a7SHong Zhang Collective on MPI_Comm 1451ba920a7SHong Zhang 1461ba920a7SHong Zhang Input Parameter: 1471ba920a7SHong Zhang + psubcomm - PetscSubcomm context 1481ba920a7SHong Zhang . color - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator. 1491ba920a7SHong Zhang . subrank - rank in the subcommunicator 1501ba920a7SHong Zhang - duprank - rank in the dupparent (see PetscSubcomm) 1511ba920a7SHong Zhang 1521ba920a7SHong Zhang Level: advanced 1531ba920a7SHong Zhang 1541ba920a7SHong Zhang .keywords: communicator, create 1551ba920a7SHong Zhang 1561ba920a7SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType() 1571ba920a7SHong Zhang @*/ 1587087cfbeSBarry Smith PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank,PetscMPIInt duprank) 159d8a68f86SHong Zhang { 1601ba920a7SHong Zhang PetscErrorCode ierr; 1611ba920a7SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 1621ba920a7SHong Zhang PetscMPIInt size; 1631ba920a7SHong Zhang 164d8a68f86SHong Zhang PetscFunctionBegin; 1651ba920a7SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 1661ba920a7SHong Zhang if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n); 1671ba920a7SHong Zhang 1681ba920a7SHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 1691ba920a7SHong Zhang 1701ba920a7SHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm 1711ba920a7SHong Zhang if duprank is not a valid number, then dupcomm is not created - not all applications require dupcomm! */ 1721ba920a7SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1737c764164SBarry Smith if (duprank == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"duprank==PETSC_DECIDE is not supported yet"); 1747c764164SBarry Smith else if (duprank >= 0 && duprank < size) { 1751ba920a7SHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 1761ba920a7SHong Zhang } 1770298fd71SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 1780298fd71SBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr); 179b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 180b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 181a297a907SKarl Rupp 1821ba920a7SHong Zhang psubcomm->color = color; 183d8a68f86SHong Zhang PetscFunctionReturn(0); 184d8a68f86SHong Zhang } 185638faf0bSHong Zhang 186cd05a4c0SHong Zhang #undef __FUNCT__ 187cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy" 1886bf464f9SBarry Smith PetscErrorCode PetscSubcommDestroy(PetscSubcomm *psubcomm) 189cd05a4c0SHong Zhang { 190cd05a4c0SHong Zhang PetscErrorCode ierr; 191cd05a4c0SHong Zhang 192cd05a4c0SHong Zhang PetscFunctionBegin; 1936bf464f9SBarry Smith if (!*psubcomm) PetscFunctionReturn(0); 194aa9c1079SBarry Smith ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr); 195aa9c1079SBarry Smith ierr = PetscCommDestroy(&(*psubcomm)->comm);CHKERRQ(ierr); 196e37c6257SHong Zhang ierr = PetscFree((*psubcomm)->subsize);CHKERRQ(ierr); 1976bf464f9SBarry Smith ierr = PetscFree((*psubcomm));CHKERRQ(ierr); 198cd05a4c0SHong Zhang PetscFunctionReturn(0); 199cd05a4c0SHong Zhang } 200cd05a4c0SHong Zhang 201cd05a4c0SHong Zhang #undef __FUNCT__ 202cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate" 203ab8c242fSMatthew Knepley /*@C 204cd05a4c0SHong Zhang PetscSubcommCreate - Create a PetscSubcomm context. 205cd05a4c0SHong Zhang 206cd05a4c0SHong Zhang Collective on MPI_Comm 207cd05a4c0SHong Zhang 208cd05a4c0SHong Zhang Input Parameter: 2099873d53eSJed Brown . comm - MPI communicator 210cd05a4c0SHong Zhang 211cd05a4c0SHong Zhang Output Parameter: 212cd05a4c0SHong Zhang . psubcomm - location to store the PetscSubcomm context 213cd05a4c0SHong Zhang 214638faf0bSHong Zhang Level: advanced 215cd05a4c0SHong Zhang 216638faf0bSHong Zhang .keywords: communicator, create 217638faf0bSHong Zhang 218638faf0bSHong Zhang .seealso: PetscSubcommDestroy() 219638faf0bSHong Zhang @*/ 2207087cfbeSBarry Smith PetscErrorCode PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm) 221638faf0bSHong Zhang { 222638faf0bSHong Zhang PetscErrorCode ierr; 223d3b23db5SHong Zhang PetscMPIInt rank,size; 224638faf0bSHong Zhang 225638faf0bSHong Zhang PetscFunctionBegin; 226d8a68f86SHong Zhang ierr = PetscNew(struct _n_PetscSubcomm,psubcomm);CHKERRQ(ierr); 227a297a907SKarl Rupp 228d3b23db5SHong Zhang /* set defaults */ 229d3b23db5SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 230d3b23db5SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 231*f68be91cSHong Zhang 232d8a68f86SHong Zhang (*psubcomm)->parent = comm; 233d3b23db5SHong Zhang (*psubcomm)->dupparent = comm; 234d3b23db5SHong Zhang (*psubcomm)->comm = PETSC_COMM_SELF; 235d3b23db5SHong Zhang (*psubcomm)->n = size; 236d3b23db5SHong Zhang (*psubcomm)->color = rank; 237e37c6257SHong Zhang (*psubcomm)->subsize = NULL; 238d3b23db5SHong Zhang (*psubcomm)->type = PETSC_SUBCOMM_INTERLACED; 239638faf0bSHong Zhang PetscFunctionReturn(0); 240638faf0bSHong Zhang } 241638faf0bSHong Zhang 242638faf0bSHong Zhang #undef __FUNCT__ 24353c77d0aSJed Brown #define __FUNCT__ "PetscSubcommCreate_contiguous" 244d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm) 245638faf0bSHong Zhang { 246638faf0bSHong Zhang PetscErrorCode ierr; 247d6037b41SHong Zhang PetscMPIInt rank,size,*subsize,duprank=-1,subrank=-1; 248d8a68f86SHong Zhang PetscInt np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n; 249d8a68f86SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 250638faf0bSHong Zhang 251638faf0bSHong Zhang PetscFunctionBegin; 25255e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 25355e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 25455e3b8d2SHong Zhang 255638faf0bSHong Zhang /* get size of each subcommunicator */ 256638faf0bSHong Zhang ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr); 257a297a907SKarl Rupp 258638faf0bSHong Zhang np_subcomm = size/nsubcomm; 259638faf0bSHong Zhang nleftover = size - nsubcomm*np_subcomm; 260638faf0bSHong Zhang for (i=0; i<nsubcomm; i++) { 261638faf0bSHong Zhang subsize[i] = np_subcomm; 262638faf0bSHong Zhang if (i<nleftover) subsize[i]++; 263638faf0bSHong Zhang } 264638faf0bSHong Zhang 265638faf0bSHong Zhang /* get color and subrank of this proc */ 266638faf0bSHong Zhang rankstart = 0; 267638faf0bSHong Zhang for (i=0; i<nsubcomm; i++) { 268638faf0bSHong Zhang if (rank >= rankstart && rank < rankstart+subsize[i]) { 269638faf0bSHong Zhang color = i; 270638faf0bSHong Zhang subrank = rank - rankstart; 271638faf0bSHong Zhang duprank = rank; 272638faf0bSHong Zhang break; 273a297a907SKarl Rupp } else rankstart += subsize[i]; 274638faf0bSHong Zhang } 275638faf0bSHong Zhang 276638faf0bSHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 277638faf0bSHong Zhang 278638faf0bSHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 279638faf0bSHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 28052af3f7aSShri Abhyankar { 28152af3f7aSShri Abhyankar PetscThreadComm tcomm; 28252af3f7aSShri Abhyankar ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr); 28352af3f7aSShri Abhyankar ierr = MPI_Attr_put(dupcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr); 28452af3f7aSShri Abhyankar tcomm->refct++; 28552af3f7aSShri Abhyankar ierr = MPI_Attr_put(subcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr); 28652af3f7aSShri Abhyankar tcomm->refct++; 28752af3f7aSShri Abhyankar } 2880298fd71SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 2890298fd71SBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr); 290b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 291b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 292a297a907SKarl Rupp 293d8a68f86SHong Zhang psubcomm->color = color; 294e37c6257SHong Zhang psubcomm->subsize = subsize; 295f38d543fSHong Zhang psubcomm->type = PETSC_SUBCOMM_CONTIGUOUS; 296638faf0bSHong Zhang PetscFunctionReturn(0); 297638faf0bSHong Zhang } 298638faf0bSHong Zhang 299638faf0bSHong Zhang #undef __FUNCT__ 300638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced" 301638faf0bSHong Zhang /* 302638faf0bSHong Zhang Note: 303638faf0bSHong Zhang In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators 30445fc02eaSBarry Smith by iteratively taking a process into a subcommunicator. 305cd05a4c0SHong Zhang Example: size=4, nsubcomm=(*psubcomm)->n=3 306cd05a4c0SHong Zhang comm=(*psubcomm)->parent: 307cd05a4c0SHong Zhang rank: [0] [1] [2] [3] 308cd05a4c0SHong Zhang color: 0 1 2 0 309cd05a4c0SHong Zhang 310cd05a4c0SHong Zhang subcomm=(*psubcomm)->comm: 311cd05a4c0SHong Zhang subrank: [0] [0] [0] [1] 312cd05a4c0SHong Zhang 313cd05a4c0SHong Zhang dupcomm=(*psubcomm)->dupparent: 314cd05a4c0SHong Zhang duprank: [0] [2] [3] [1] 315cd05a4c0SHong Zhang 316cd05a4c0SHong Zhang Here, subcomm[color = 0] has subsize=2, owns process [0] and [3] 317cd05a4c0SHong Zhang subcomm[color = 1] has subsize=1, owns process [1] 318cd05a4c0SHong Zhang subcomm[color = 2] has subsize=1, owns process [2] 319cd05a4c0SHong Zhang dupcomm has same number of processes as comm, and its duprank maps 320cd05a4c0SHong Zhang processes in subcomm contiguously into a 1d array: 321cd05a4c0SHong Zhang duprank: [0] [1] [2] [3] 322cd05a4c0SHong Zhang rank: [0] [3] [1] [2] 323cd05a4c0SHong Zhang subcomm[0] subcomm[1] subcomm[2] 324638faf0bSHong Zhang */ 325cd05a4c0SHong Zhang 326d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm) 327cd05a4c0SHong Zhang { 328cd05a4c0SHong Zhang PetscErrorCode ierr; 329cd05a4c0SHong Zhang PetscMPIInt rank,size,*subsize,duprank,subrank; 330d8a68f86SHong Zhang PetscInt np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n; 331d8a68f86SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 332cd05a4c0SHong Zhang 333cd05a4c0SHong Zhang PetscFunctionBegin; 33455e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 33555e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 33655e3b8d2SHong Zhang 337cd05a4c0SHong Zhang /* get size of each subcommunicator */ 338cd05a4c0SHong Zhang ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr); 339a297a907SKarl Rupp 340cd05a4c0SHong Zhang np_subcomm = size/nsubcomm; 341cd05a4c0SHong Zhang nleftover = size - nsubcomm*np_subcomm; 342cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++) { 343cd05a4c0SHong Zhang subsize[i] = np_subcomm; 344cd05a4c0SHong Zhang if (i<nleftover) subsize[i]++; 345cd05a4c0SHong Zhang } 346cd05a4c0SHong Zhang 347cd05a4c0SHong Zhang /* find color for this proc */ 348cd05a4c0SHong Zhang color = rank%nsubcomm; 349cd05a4c0SHong Zhang subrank = rank/nsubcomm; 350cd05a4c0SHong Zhang 351cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 352cd05a4c0SHong Zhang 353cd05a4c0SHong Zhang j = 0; duprank = 0; 354cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++) { 355cd05a4c0SHong Zhang if (j == color) { 356cd05a4c0SHong Zhang duprank += subrank; 357cd05a4c0SHong Zhang break; 358cd05a4c0SHong Zhang } 359cd05a4c0SHong Zhang duprank += subsize[i]; j++; 360cd05a4c0SHong Zhang } 361cd05a4c0SHong Zhang 362cd05a4c0SHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 363cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 36452af3f7aSShri Abhyankar { 36552af3f7aSShri Abhyankar PetscThreadComm tcomm; 36652af3f7aSShri Abhyankar ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr); 36752af3f7aSShri Abhyankar ierr = MPI_Attr_put(dupcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr); 36852af3f7aSShri Abhyankar tcomm->refct++; 36952af3f7aSShri Abhyankar ierr = MPI_Attr_put(subcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr); 37052af3f7aSShri Abhyankar tcomm->refct++; 37152af3f7aSShri Abhyankar } 3720298fd71SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 3730298fd71SBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr); 374b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 375b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 376a297a907SKarl Rupp 377d8a68f86SHong Zhang psubcomm->color = color; 378e37c6257SHong Zhang psubcomm->subsize = subsize; 379f38d543fSHong Zhang psubcomm->type = PETSC_SUBCOMM_INTERLACED; 380cd05a4c0SHong Zhang PetscFunctionReturn(0); 381cd05a4c0SHong Zhang } 382638faf0bSHong Zhang 383638faf0bSHong Zhang 384