17d0a6c19SBarry Smith 2cd05a4c0SHong Zhang /* 3cd05a4c0SHong Zhang Provides utility routines for split MPI communicator. 4cd05a4c0SHong Zhang */ 5c6db04a5SJed Brown #include <petscsys.h> /*I "petscsys.h" I*/ 6053d1c95SHong Zhang #include <petscviewer.h> 7cd05a4c0SHong Zhang 86a6fc655SJed Brown const char *const PetscSubcommTypes[] = {"GENERAL","CONTIGUOUS","INTERLACED","PetscSubcommType","PETSC_SUBCOMM_",0}; 96a6fc655SJed Brown 10*95c0884eSLisandro Dalcin static PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm); 11*95c0884eSLisandro Dalcin static PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm); 12e5acf8a4SHong Zhang 13e5acf8a4SHong Zhang /*@C 14e5acf8a4SHong Zhang PetscSubcommSetFromOptions - Allows setting options from a PetscSubcomm 15e5acf8a4SHong Zhang 16e5acf8a4SHong Zhang Collective on PetscSubcomm 17e5acf8a4SHong Zhang 18e5acf8a4SHong Zhang Input Parameter: 19e5acf8a4SHong Zhang . psubcomm - PetscSubcomm context 20e5acf8a4SHong Zhang 21e5acf8a4SHong Zhang Level: beginner 22e5acf8a4SHong Zhang 23e5acf8a4SHong Zhang @*/ 24f68be91cSHong Zhang PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm psubcomm) 25f68be91cSHong Zhang { 26f68be91cSHong Zhang PetscErrorCode ierr; 2745487dadSJed Brown PetscSubcommType type; 28f68be91cSHong Zhang PetscBool flg; 29f68be91cSHong Zhang 30f68be91cSHong Zhang PetscFunctionBegin; 31f68be91cSHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Must call PetscSubcommCreate firt"); 3260f5f21cSHong Zhang 33e5acf8a4SHong Zhang ierr = PetscOptionsBegin(psubcomm->parent,psubcomm->subcommprefix,"Options for PetscSubcomm",NULL);CHKERRQ(ierr); 3460f5f21cSHong Zhang ierr = PetscOptionsEnum("-psubcomm_type",NULL,NULL,PetscSubcommTypes,(PetscEnum)psubcomm->type,(PetscEnum*)&type,&flg);CHKERRQ(ierr); 35f68be91cSHong Zhang if (flg && psubcomm->type != type) { 36f68be91cSHong Zhang /* free old structures */ 37f68be91cSHong Zhang ierr = PetscCommDestroy(&(psubcomm)->dupparent);CHKERRQ(ierr); 38306c2d5bSBarry Smith ierr = PetscCommDestroy(&(psubcomm)->child);CHKERRQ(ierr); 39f68be91cSHong Zhang ierr = PetscFree((psubcomm)->subsize);CHKERRQ(ierr); 40f68be91cSHong Zhang switch (type) { 4145487dadSJed Brown case PETSC_SUBCOMM_GENERAL: 42b3a4ddeeSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Runtime option PETSC_SUBCOMM_GENERAL is not supported, use PetscSubcommSetTypeGeneral()"); 4345487dadSJed Brown case PETSC_SUBCOMM_CONTIGUOUS: 44f68be91cSHong Zhang ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr); 45f68be91cSHong Zhang break; 4645487dadSJed Brown case PETSC_SUBCOMM_INTERLACED: 47f68be91cSHong Zhang ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr); 48f68be91cSHong Zhang break; 49f68be91cSHong Zhang default: 5045487dadSJed Brown SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSubcommType %s is not supported yet",PetscSubcommTypes[type]); 51f68be91cSHong Zhang } 52f68be91cSHong Zhang } 5319171117SHong Zhang 5460f5f21cSHong Zhang ierr = PetscOptionsName("-psubcomm_view","Triggers display of PetscSubcomm context","PetscSubcommView",&flg);CHKERRQ(ierr); 5519171117SHong Zhang if (flg) { 5619171117SHong Zhang ierr = PetscSubcommView(psubcomm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 5719171117SHong Zhang } 5860f5f21cSHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 59f68be91cSHong Zhang PetscFunctionReturn(0); 60f68be91cSHong Zhang } 61d8a68f86SHong Zhang 62e5acf8a4SHong Zhang /*@C 63e5acf8a4SHong Zhang PetscSubcommSetOptionsPrefix - Sets the prefix used for searching for all 64e5acf8a4SHong Zhang PetscSubcomm items in the options database. 65e5acf8a4SHong Zhang 66e5acf8a4SHong Zhang Logically collective on PetscSubcomm. 67e5acf8a4SHong Zhang 68e5acf8a4SHong Zhang Level: Intermediate 69e5acf8a4SHong Zhang 70e5acf8a4SHong Zhang Input Parameters: 71e5acf8a4SHong Zhang + psubcomm - PetscSubcomm context 72e5acf8a4SHong Zhang - prefix - the prefix to prepend all PetscSubcomm item names with. 73e5acf8a4SHong Zhang 74e5acf8a4SHong Zhang @*/ 75e5acf8a4SHong Zhang PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm psubcomm,const char pre[]) 76e5acf8a4SHong Zhang { 77e5acf8a4SHong Zhang PetscErrorCode ierr; 78e5acf8a4SHong Zhang 79e5acf8a4SHong Zhang PetscFunctionBegin; 80e5acf8a4SHong Zhang if (!pre) { 81e5acf8a4SHong Zhang ierr = PetscFree(psubcomm->subcommprefix);CHKERRQ(ierr); 82e5acf8a4SHong Zhang } else { 83e5acf8a4SHong Zhang if (pre[0] == '-') SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Options prefix should not begin with a hypen"); 84e5acf8a4SHong Zhang ierr = PetscFree(psubcomm->subcommprefix);CHKERRQ(ierr); 85e5acf8a4SHong Zhang ierr = PetscStrallocpy(pre,&(psubcomm->subcommprefix));CHKERRQ(ierr); 86e5acf8a4SHong Zhang } 87e5acf8a4SHong Zhang PetscFunctionReturn(0); 88e5acf8a4SHong Zhang } 89e5acf8a4SHong Zhang 90e5acf8a4SHong Zhang /*@C 9123072422SBarry Smith PetscSubcommView - Views a PetscSubcomm of values as either ASCII text or a binary file 92e5acf8a4SHong Zhang 93e5acf8a4SHong Zhang Collective on PetscSubcomm 94e5acf8a4SHong Zhang 95e5acf8a4SHong Zhang Input Parameter: 96e5acf8a4SHong Zhang + psubcomm - PetscSubcomm context 97e5acf8a4SHong Zhang - viewer - location to view the values 98e5acf8a4SHong Zhang 99e5acf8a4SHong Zhang Level: beginner 100e5acf8a4SHong Zhang @*/ 101053d1c95SHong Zhang PetscErrorCode PetscSubcommView(PetscSubcomm psubcomm,PetscViewer viewer) 102053d1c95SHong Zhang { 103053d1c95SHong Zhang PetscErrorCode ierr; 104053d1c95SHong Zhang PetscBool iascii; 105053d1c95SHong Zhang PetscViewerFormat format; 106053d1c95SHong Zhang 107053d1c95SHong Zhang PetscFunctionBegin; 108053d1c95SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 109053d1c95SHong Zhang if (iascii) { 110053d1c95SHong Zhang ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 111053d1c95SHong Zhang if (format == PETSC_VIEWER_DEFAULT) { 112053d1c95SHong Zhang MPI_Comm comm=psubcomm->parent; 113053d1c95SHong Zhang PetscMPIInt rank,size,subsize,subrank,duprank; 114053d1c95SHong Zhang 115053d1c95SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 11645487dadSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"PetscSubcomm type %s with total %d MPI processes:\n",PetscSubcommTypes[psubcomm->type],size);CHKERRQ(ierr); 117053d1c95SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 118306c2d5bSBarry Smith ierr = MPI_Comm_size(psubcomm->child,&subsize);CHKERRQ(ierr); 119306c2d5bSBarry Smith ierr = MPI_Comm_rank(psubcomm->child,&subrank);CHKERRQ(ierr); 120053d1c95SHong Zhang ierr = MPI_Comm_rank(psubcomm->dupparent,&duprank);CHKERRQ(ierr); 121302440fdSBarry Smith ierr = PetscSynchronizedPrintf(comm," [%d], color %d, sub-size %d, sub-rank %d, duprank %d\n",rank,psubcomm->color,subsize,subrank,duprank);CHKERRQ(ierr); 1220ec8b6e3SBarry Smith ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr); 123053d1c95SHong Zhang } 124053d1c95SHong Zhang } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported yet"); 125053d1c95SHong Zhang PetscFunctionReturn(0); 126053d1c95SHong Zhang } 127053d1c95SHong Zhang 128d8a68f86SHong Zhang /*@C 129d8a68f86SHong Zhang PetscSubcommSetNumber - Set total number of subcommunicators. 130d8a68f86SHong Zhang 131d8a68f86SHong Zhang Collective on MPI_Comm 132d8a68f86SHong Zhang 133d8a68f86SHong Zhang Input Parameter: 134d8a68f86SHong Zhang + psubcomm - PetscSubcomm context 135d8a68f86SHong Zhang - nsubcomm - the total number of subcommunicators in psubcomm 136d8a68f86SHong Zhang 137d8a68f86SHong Zhang Level: advanced 138d8a68f86SHong Zhang 139d8a68f86SHong Zhang .keywords: communicator 140d8a68f86SHong Zhang 141d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetType(),PetscSubcommSetTypeGeneral() 142d8a68f86SHong Zhang @*/ 1437087cfbeSBarry Smith PetscErrorCode PetscSubcommSetNumber(PetscSubcomm psubcomm,PetscInt nsubcomm) 144d8a68f86SHong Zhang { 145d8a68f86SHong Zhang PetscErrorCode ierr; 146d8a68f86SHong Zhang MPI_Comm comm=psubcomm->parent; 14745487dadSJed Brown PetscMPIInt msub,size; 148d8a68f86SHong Zhang 149d8a68f86SHong Zhang PetscFunctionBegin; 150d8a68f86SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate() first"); 151d8a68f86SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 15245487dadSJed Brown ierr = PetscMPIIntCast(nsubcomm,&msub);CHKERRQ(ierr); 15345487dadSJed Brown if (msub < 1 || msub > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Num of subcommunicators %d cannot be < 1 or > input comm size %d",msub,size); 154d8a68f86SHong Zhang 15545487dadSJed Brown psubcomm->n = msub; 156d8a68f86SHong Zhang PetscFunctionReturn(0); 157d8a68f86SHong Zhang } 158d8a68f86SHong Zhang 159d8a68f86SHong Zhang /*@C 160d8a68f86SHong Zhang PetscSubcommSetType - Set type of subcommunicators. 161d8a68f86SHong Zhang 162d8a68f86SHong Zhang Collective on MPI_Comm 163d8a68f86SHong Zhang 164d8a68f86SHong Zhang Input Parameter: 165d8a68f86SHong Zhang + psubcomm - PetscSubcomm context 1661ba920a7SHong Zhang - subcommtype - subcommunicator type, PETSC_SUBCOMM_CONTIGUOUS,PETSC_SUBCOMM_INTERLACED 167d8a68f86SHong Zhang 168d8a68f86SHong Zhang Level: advanced 169d8a68f86SHong Zhang 170d8a68f86SHong Zhang .keywords: communicator 171d8a68f86SHong Zhang 172d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral() 173d8a68f86SHong Zhang @*/ 1747c764164SBarry Smith PetscErrorCode PetscSubcommSetType(PetscSubcomm psubcomm,PetscSubcommType subcommtype) 175d8a68f86SHong Zhang { 176d8a68f86SHong Zhang PetscErrorCode ierr; 177d8a68f86SHong Zhang 178d8a68f86SHong Zhang PetscFunctionBegin; 179d8a68f86SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 18045487dadSJed Brown if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()",psubcomm->n); 181d8a68f86SHong Zhang 182d8a68f86SHong Zhang if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS) { 183d8a68f86SHong Zhang ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr); 184d8a68f86SHong Zhang } else if (subcommtype == PETSC_SUBCOMM_INTERLACED) { 185d8a68f86SHong Zhang ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr); 18645487dadSJed Brown } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %s is not supported yet",PetscSubcommTypes[subcommtype]); 187d8a68f86SHong Zhang PetscFunctionReturn(0); 188d8a68f86SHong Zhang } 189d8a68f86SHong Zhang 1901ba920a7SHong Zhang /*@C 19165d9b8f1SHong Zhang PetscSubcommSetTypeGeneral - Set a PetscSubcomm from user's specifications 1921ba920a7SHong Zhang 1931ba920a7SHong Zhang Collective on MPI_Comm 1941ba920a7SHong Zhang 1951ba920a7SHong Zhang Input Parameter: 1961ba920a7SHong Zhang + psubcomm - PetscSubcomm context 1971ba920a7SHong Zhang . color - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator. 19865d9b8f1SHong Zhang - subrank - rank in the subcommunicator 1991ba920a7SHong Zhang 2001ba920a7SHong Zhang Level: advanced 2011ba920a7SHong Zhang 2021ba920a7SHong Zhang .keywords: communicator, create 2031ba920a7SHong Zhang 2041ba920a7SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType() 2051ba920a7SHong Zhang @*/ 20665d9b8f1SHong Zhang PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank) 207d8a68f86SHong Zhang { 2081ba920a7SHong Zhang PetscErrorCode ierr; 2091ba920a7SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 210c9e2ceb8SHong Zhang PetscMPIInt size,icolor,duprank,*recvbuf,sendbuf[3],mysubsize,rank,*subsize; 21145487dadSJed Brown PetscMPIInt i,nsubcomm=psubcomm->n; 2121ba920a7SHong Zhang 213d8a68f86SHong Zhang PetscFunctionBegin; 2141ba920a7SHong Zhang if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()"); 21545487dadSJed Brown if (nsubcomm < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()",nsubcomm); 2161ba920a7SHong Zhang 2171ba920a7SHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 2181ba920a7SHong Zhang 21965d9b8f1SHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 2209abe469cSDmitry Karpeev /* TODO: this can be done in an ostensibly scalale way (i.e., without allocating an array of size 'size') as is done in PetscObjectsCreateGlobalOrdering(). */ 2211ba920a7SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 222785e854fSJed Brown ierr = PetscMalloc1(2*size,&recvbuf);CHKERRQ(ierr); 22365d9b8f1SHong Zhang 22465d9b8f1SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 225c9e2ceb8SHong Zhang ierr = MPI_Comm_size(subcomm,&mysubsize);CHKERRQ(ierr); 22665d9b8f1SHong Zhang 22765d9b8f1SHong Zhang sendbuf[0] = color; 228c9e2ceb8SHong Zhang sendbuf[1] = mysubsize; 229c9e2ceb8SHong Zhang ierr = MPI_Allgather(sendbuf,2,MPI_INT,recvbuf,2,MPI_INT,comm);CHKERRQ(ierr); 23065d9b8f1SHong Zhang 231eca4a452SHong Zhang ierr = PetscCalloc1(nsubcomm,&subsize);CHKERRQ(ierr); 2329972d2ceSHong Zhang for (i=0; i<2*size; i+=2) { 233c9e2ceb8SHong Zhang subsize[recvbuf[i]] = recvbuf[i+1]; 2341ba920a7SHong Zhang } 23565d9b8f1SHong Zhang ierr = PetscFree(recvbuf);CHKERRQ(ierr); 23665d9b8f1SHong Zhang 23765d9b8f1SHong Zhang duprank = 0; 238c9e2ceb8SHong Zhang for (icolor=0; icolor<nsubcomm; icolor++) { 23965d9b8f1SHong Zhang if (icolor != color) { /* not color of this process */ 240c9e2ceb8SHong Zhang duprank += subsize[icolor]; 24165d9b8f1SHong Zhang } else { 24265d9b8f1SHong Zhang duprank += subrank; 24365d9b8f1SHong Zhang break; 24465d9b8f1SHong Zhang } 24565d9b8f1SHong Zhang } 24665d9b8f1SHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 24765d9b8f1SHong Zhang 2480298fd71SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 249306c2d5bSBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr); 250b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 251b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 252a297a907SKarl Rupp 2531ba920a7SHong Zhang psubcomm->color = color; 254c9e2ceb8SHong Zhang psubcomm->subsize = subsize; 255c9e2ceb8SHong Zhang psubcomm->type = PETSC_SUBCOMM_GENERAL; 256d8a68f86SHong Zhang PetscFunctionReturn(0); 257d8a68f86SHong Zhang } 258638faf0bSHong Zhang 25989587e68SDave May /*@C 26089587e68SDave May PetscSubcommDestroy - Destroys a PetscSubcomm object 26189587e68SDave May 26289587e68SDave May Collective on PetscSubcomm 26389587e68SDave May 26489587e68SDave May Input Parameter: 26589587e68SDave May . psubcomm - the PetscSubcomm context 26689587e68SDave May 26789587e68SDave May Level: advanced 26889587e68SDave May 26989587e68SDave May .seealso: PetscSubcommCreate(),PetscSubcommSetType() 27089587e68SDave May @*/ 2716bf464f9SBarry Smith PetscErrorCode PetscSubcommDestroy(PetscSubcomm *psubcomm) 272cd05a4c0SHong Zhang { 273cd05a4c0SHong Zhang PetscErrorCode ierr; 274cd05a4c0SHong Zhang 275cd05a4c0SHong Zhang PetscFunctionBegin; 2766bf464f9SBarry Smith if (!*psubcomm) PetscFunctionReturn(0); 277aa9c1079SBarry Smith ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr); 278306c2d5bSBarry Smith ierr = PetscCommDestroy(&(*psubcomm)->child);CHKERRQ(ierr); 279e37c6257SHong Zhang ierr = PetscFree((*psubcomm)->subsize);CHKERRQ(ierr); 280e5acf8a4SHong Zhang if ((*psubcomm)->subcommprefix) { ierr = PetscFree((*psubcomm)->subcommprefix);CHKERRQ(ierr); } 2816bf464f9SBarry Smith ierr = PetscFree((*psubcomm));CHKERRQ(ierr); 282cd05a4c0SHong Zhang PetscFunctionReturn(0); 283cd05a4c0SHong Zhang } 284cd05a4c0SHong Zhang 285ab8c242fSMatthew Knepley /*@C 286cd05a4c0SHong Zhang PetscSubcommCreate - Create a PetscSubcomm context. 287cd05a4c0SHong Zhang 288cd05a4c0SHong Zhang Collective on MPI_Comm 289cd05a4c0SHong Zhang 290cd05a4c0SHong Zhang Input Parameter: 2919873d53eSJed Brown . comm - MPI communicator 292cd05a4c0SHong Zhang 293cd05a4c0SHong Zhang Output Parameter: 294cd05a4c0SHong Zhang . psubcomm - location to store the PetscSubcomm context 295cd05a4c0SHong Zhang 296638faf0bSHong Zhang Level: advanced 297cd05a4c0SHong Zhang 298638faf0bSHong Zhang .keywords: communicator, create 299638faf0bSHong Zhang 300638faf0bSHong Zhang .seealso: PetscSubcommDestroy() 301638faf0bSHong Zhang @*/ 3027087cfbeSBarry Smith PetscErrorCode PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm) 303638faf0bSHong Zhang { 304638faf0bSHong Zhang PetscErrorCode ierr; 305d3b23db5SHong Zhang PetscMPIInt rank,size; 306638faf0bSHong Zhang 307638faf0bSHong Zhang PetscFunctionBegin; 308b00a9115SJed Brown ierr = PetscNew(psubcomm);CHKERRQ(ierr); 309a297a907SKarl Rupp 310d3b23db5SHong Zhang /* set defaults */ 311d3b23db5SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 312d3b23db5SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 313f68be91cSHong Zhang 314d8a68f86SHong Zhang (*psubcomm)->parent = comm; 315d3b23db5SHong Zhang (*psubcomm)->dupparent = comm; 316306c2d5bSBarry Smith (*psubcomm)->child = PETSC_COMM_SELF; 317d3b23db5SHong Zhang (*psubcomm)->n = size; 318d3b23db5SHong Zhang (*psubcomm)->color = rank; 319e37c6257SHong Zhang (*psubcomm)->subsize = NULL; 320d3b23db5SHong Zhang (*psubcomm)->type = PETSC_SUBCOMM_INTERLACED; 321638faf0bSHong Zhang PetscFunctionReturn(0); 322638faf0bSHong Zhang } 323638faf0bSHong Zhang 324*95c0884eSLisandro Dalcin static PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm) 325638faf0bSHong Zhang { 326638faf0bSHong Zhang PetscErrorCode ierr; 327d6037b41SHong Zhang PetscMPIInt rank,size,*subsize,duprank=-1,subrank=-1; 32845487dadSJed Brown PetscMPIInt np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n; 329d8a68f86SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 330638faf0bSHong Zhang 331638faf0bSHong Zhang PetscFunctionBegin; 33255e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 33355e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 33455e3b8d2SHong Zhang 335638faf0bSHong Zhang /* get size of each subcommunicator */ 336854ce69bSBarry Smith ierr = PetscMalloc1(1+nsubcomm,&subsize);CHKERRQ(ierr); 337a297a907SKarl Rupp 338638faf0bSHong Zhang np_subcomm = size/nsubcomm; 339638faf0bSHong Zhang nleftover = size - nsubcomm*np_subcomm; 340638faf0bSHong Zhang for (i=0; i<nsubcomm; i++) { 341638faf0bSHong Zhang subsize[i] = np_subcomm; 342638faf0bSHong Zhang if (i<nleftover) subsize[i]++; 343638faf0bSHong Zhang } 344638faf0bSHong Zhang 345638faf0bSHong Zhang /* get color and subrank of this proc */ 346638faf0bSHong Zhang rankstart = 0; 347638faf0bSHong Zhang for (i=0; i<nsubcomm; i++) { 348638faf0bSHong Zhang if (rank >= rankstart && rank < rankstart+subsize[i]) { 349638faf0bSHong Zhang color = i; 350638faf0bSHong Zhang subrank = rank - rankstart; 351638faf0bSHong Zhang duprank = rank; 352638faf0bSHong Zhang break; 353a297a907SKarl Rupp } else rankstart += subsize[i]; 354638faf0bSHong Zhang } 355638faf0bSHong Zhang 356638faf0bSHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 357638faf0bSHong Zhang 358638faf0bSHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 359638faf0bSHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 3600298fd71SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 361306c2d5bSBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr); 362b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 363b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 364a297a907SKarl Rupp 365d8a68f86SHong Zhang psubcomm->color = color; 366e37c6257SHong Zhang psubcomm->subsize = subsize; 367f38d543fSHong Zhang psubcomm->type = PETSC_SUBCOMM_CONTIGUOUS; 368638faf0bSHong Zhang PetscFunctionReturn(0); 369638faf0bSHong Zhang } 370638faf0bSHong Zhang 371638faf0bSHong Zhang /* 372638faf0bSHong Zhang Note: 373638faf0bSHong Zhang In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators 37445fc02eaSBarry Smith by iteratively taking a process into a subcommunicator. 375cd05a4c0SHong Zhang Example: size=4, nsubcomm=(*psubcomm)->n=3 376cd05a4c0SHong Zhang comm=(*psubcomm)->parent: 377cd05a4c0SHong Zhang rank: [0] [1] [2] [3] 378cd05a4c0SHong Zhang color: 0 1 2 0 379cd05a4c0SHong Zhang 380cd05a4c0SHong Zhang subcomm=(*psubcomm)->comm: 381cd05a4c0SHong Zhang subrank: [0] [0] [0] [1] 382cd05a4c0SHong Zhang 383cd05a4c0SHong Zhang dupcomm=(*psubcomm)->dupparent: 384cd05a4c0SHong Zhang duprank: [0] [2] [3] [1] 385cd05a4c0SHong Zhang 386cd05a4c0SHong Zhang Here, subcomm[color = 0] has subsize=2, owns process [0] and [3] 387cd05a4c0SHong Zhang subcomm[color = 1] has subsize=1, owns process [1] 388cd05a4c0SHong Zhang subcomm[color = 2] has subsize=1, owns process [2] 389cd05a4c0SHong Zhang dupcomm has same number of processes as comm, and its duprank maps 390cd05a4c0SHong Zhang processes in subcomm contiguously into a 1d array: 391cd05a4c0SHong Zhang duprank: [0] [1] [2] [3] 392cd05a4c0SHong Zhang rank: [0] [3] [1] [2] 393cd05a4c0SHong Zhang subcomm[0] subcomm[1] subcomm[2] 394638faf0bSHong Zhang */ 395cd05a4c0SHong Zhang 396*95c0884eSLisandro Dalcin static PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm) 397cd05a4c0SHong Zhang { 398cd05a4c0SHong Zhang PetscErrorCode ierr; 399cd05a4c0SHong Zhang PetscMPIInt rank,size,*subsize,duprank,subrank; 40045487dadSJed Brown PetscMPIInt np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n; 401d8a68f86SHong Zhang MPI_Comm subcomm=0,dupcomm=0,comm=psubcomm->parent; 402cd05a4c0SHong Zhang 403cd05a4c0SHong Zhang PetscFunctionBegin; 40455e3b8d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 40555e3b8d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 40655e3b8d2SHong Zhang 407cd05a4c0SHong Zhang /* get size of each subcommunicator */ 408854ce69bSBarry Smith ierr = PetscMalloc1(1+nsubcomm,&subsize);CHKERRQ(ierr); 409a297a907SKarl Rupp 410cd05a4c0SHong Zhang np_subcomm = size/nsubcomm; 411cd05a4c0SHong Zhang nleftover = size - nsubcomm*np_subcomm; 412cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++) { 413cd05a4c0SHong Zhang subsize[i] = np_subcomm; 414cd05a4c0SHong Zhang if (i<nleftover) subsize[i]++; 415cd05a4c0SHong Zhang } 416cd05a4c0SHong Zhang 417cd05a4c0SHong Zhang /* find color for this proc */ 418cd05a4c0SHong Zhang color = rank%nsubcomm; 419cd05a4c0SHong Zhang subrank = rank/nsubcomm; 420cd05a4c0SHong Zhang 421cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr); 422cd05a4c0SHong Zhang 423cd05a4c0SHong Zhang j = 0; duprank = 0; 424cd05a4c0SHong Zhang for (i=0; i<nsubcomm; i++) { 425cd05a4c0SHong Zhang if (j == color) { 426cd05a4c0SHong Zhang duprank += subrank; 427cd05a4c0SHong Zhang break; 428cd05a4c0SHong Zhang } 429cd05a4c0SHong Zhang duprank += subsize[i]; j++; 430cd05a4c0SHong Zhang } 431cd05a4c0SHong Zhang 432cd05a4c0SHong Zhang /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */ 433cd05a4c0SHong Zhang ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr); 4340298fd71SBarry Smith ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr); 435306c2d5bSBarry Smith ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr); 436b89831e5SBarry Smith ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr); 437b89831e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 438a297a907SKarl Rupp 439d8a68f86SHong Zhang psubcomm->color = color; 440e37c6257SHong Zhang psubcomm->subsize = subsize; 441f38d543fSHong Zhang psubcomm->type = PETSC_SUBCOMM_INTERLACED; 442cd05a4c0SHong Zhang PetscFunctionReturn(0); 443cd05a4c0SHong Zhang } 444638faf0bSHong Zhang 445638faf0bSHong Zhang 446