xref: /petsc/src/sys/objects/subcomm.c (revision 60f5f21c9986be7e96bc5f6a72a0920da4c6def8)
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 
10d8a68f86SHong Zhang extern PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm);
11d8a68f86SHong Zhang extern PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm);
12f68be91cSHong Zhang #undef __FUNCT__
13f68be91cSHong Zhang #define __FUNCT__ "PetscSubcommSetFromOptions"
14f68be91cSHong Zhang PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm psubcomm)
15f68be91cSHong Zhang {
16f68be91cSHong Zhang   PetscErrorCode   ierr;
1745487dadSJed Brown   PetscSubcommType type;
18f68be91cSHong Zhang   PetscBool        flg;
19f68be91cSHong Zhang 
20f68be91cSHong Zhang   PetscFunctionBegin;
21f68be91cSHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Must call PetscSubcommCreate firt");
22*60f5f21cSHong Zhang 
23*60f5f21cSHong Zhang   ierr = PetscOptionsBegin(psubcomm->parent,NULL,"Options for PetscSubcomm",NULL);CHKERRQ(ierr);
24*60f5f21cSHong Zhang   ierr = PetscOptionsEnum("-psubcomm_type",NULL,NULL,PetscSubcommTypes,(PetscEnum)psubcomm->type,(PetscEnum*)&type,&flg);CHKERRQ(ierr);
25f68be91cSHong Zhang   if (flg && psubcomm->type != type) {
26f68be91cSHong Zhang     /* free old structures */
27f68be91cSHong Zhang     ierr = PetscCommDestroy(&(psubcomm)->dupparent);CHKERRQ(ierr);
28306c2d5bSBarry Smith     ierr = PetscCommDestroy(&(psubcomm)->child);CHKERRQ(ierr);
29f68be91cSHong Zhang     ierr = PetscFree((psubcomm)->subsize);CHKERRQ(ierr);
30f68be91cSHong Zhang     switch (type) {
3145487dadSJed Brown     case PETSC_SUBCOMM_GENERAL:
32b3a4ddeeSHong Zhang       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Runtime option PETSC_SUBCOMM_GENERAL is not supported, use PetscSubcommSetTypeGeneral()");
3345487dadSJed Brown     case PETSC_SUBCOMM_CONTIGUOUS:
34f68be91cSHong Zhang       ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr);
35f68be91cSHong Zhang       break;
3645487dadSJed Brown     case PETSC_SUBCOMM_INTERLACED:
37f68be91cSHong Zhang       ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr);
38f68be91cSHong Zhang       break;
39f68be91cSHong Zhang     default:
4045487dadSJed Brown       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSubcommType %s is not supported yet",PetscSubcommTypes[type]);
41f68be91cSHong Zhang     }
42f68be91cSHong Zhang   }
4319171117SHong Zhang 
44*60f5f21cSHong Zhang   ierr = PetscOptionsName("-psubcomm_view","Triggers display of PetscSubcomm context","PetscSubcommView",&flg);CHKERRQ(ierr);
4519171117SHong Zhang   if (flg) {
4619171117SHong Zhang     ierr = PetscSubcommView(psubcomm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
4719171117SHong Zhang   }
48*60f5f21cSHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
49f68be91cSHong Zhang   PetscFunctionReturn(0);
50f68be91cSHong Zhang }
51d8a68f86SHong Zhang 
52d8a68f86SHong Zhang #undef __FUNCT__
53053d1c95SHong Zhang #define __FUNCT__ "PetscSubcommView"
54053d1c95SHong Zhang PetscErrorCode PetscSubcommView(PetscSubcomm psubcomm,PetscViewer viewer)
55053d1c95SHong Zhang {
56053d1c95SHong Zhang   PetscErrorCode    ierr;
57053d1c95SHong Zhang   PetscBool         iascii;
58053d1c95SHong Zhang   PetscViewerFormat format;
59053d1c95SHong Zhang 
60053d1c95SHong Zhang   PetscFunctionBegin;
61053d1c95SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
62053d1c95SHong Zhang   if (iascii) {
63053d1c95SHong Zhang     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
64053d1c95SHong Zhang     if (format == PETSC_VIEWER_DEFAULT) {
65053d1c95SHong Zhang       MPI_Comm    comm=psubcomm->parent;
66053d1c95SHong Zhang       PetscMPIInt rank,size,subsize,subrank,duprank;
67053d1c95SHong Zhang 
68053d1c95SHong Zhang       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
6945487dadSJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"PetscSubcomm type %s with total %d MPI processes:\n",PetscSubcommTypes[psubcomm->type],size);CHKERRQ(ierr);
70053d1c95SHong Zhang       ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
71306c2d5bSBarry Smith       ierr = MPI_Comm_size(psubcomm->child,&subsize);CHKERRQ(ierr);
72306c2d5bSBarry Smith       ierr = MPI_Comm_rank(psubcomm->child,&subrank);CHKERRQ(ierr);
73053d1c95SHong Zhang       ierr = MPI_Comm_rank(psubcomm->dupparent,&duprank);CHKERRQ(ierr);
74302440fdSBarry Smith       ierr = PetscSynchronizedPrintf(comm,"  [%d], color %d, sub-size %d, sub-rank %d, duprank %d\n",rank,psubcomm->color,subsize,subrank,duprank);CHKERRQ(ierr);
750ec8b6e3SBarry Smith       ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
76053d1c95SHong Zhang     }
77053d1c95SHong Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported yet");
78053d1c95SHong Zhang   PetscFunctionReturn(0);
79053d1c95SHong Zhang }
80053d1c95SHong Zhang 
81053d1c95SHong Zhang #undef __FUNCT__
82d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetNumber"
83d8a68f86SHong Zhang /*@C
84d8a68f86SHong Zhang   PetscSubcommSetNumber - Set total number of subcommunicators.
85d8a68f86SHong Zhang 
86d8a68f86SHong Zhang    Collective on MPI_Comm
87d8a68f86SHong Zhang 
88d8a68f86SHong Zhang    Input Parameter:
89d8a68f86SHong Zhang +  psubcomm - PetscSubcomm context
90d8a68f86SHong Zhang -  nsubcomm - the total number of subcommunicators in psubcomm
91d8a68f86SHong Zhang 
92d8a68f86SHong Zhang    Level: advanced
93d8a68f86SHong Zhang 
94d8a68f86SHong Zhang .keywords: communicator
95d8a68f86SHong Zhang 
96d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetType(),PetscSubcommSetTypeGeneral()
97d8a68f86SHong Zhang @*/
987087cfbeSBarry Smith PetscErrorCode  PetscSubcommSetNumber(PetscSubcomm psubcomm,PetscInt nsubcomm)
99d8a68f86SHong Zhang {
100d8a68f86SHong Zhang   PetscErrorCode ierr;
101d8a68f86SHong Zhang   MPI_Comm       comm=psubcomm->parent;
10245487dadSJed Brown   PetscMPIInt    msub,size;
103d8a68f86SHong Zhang 
104d8a68f86SHong Zhang   PetscFunctionBegin;
105d8a68f86SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate() first");
106d8a68f86SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
10745487dadSJed Brown   ierr = PetscMPIIntCast(nsubcomm,&msub);CHKERRQ(ierr);
10845487dadSJed 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);
109d8a68f86SHong Zhang 
11045487dadSJed Brown   psubcomm->n = msub;
111d8a68f86SHong Zhang   PetscFunctionReturn(0);
112d8a68f86SHong Zhang }
113d8a68f86SHong Zhang 
114d8a68f86SHong Zhang #undef __FUNCT__
115d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetType"
116d8a68f86SHong Zhang /*@C
117d8a68f86SHong Zhang   PetscSubcommSetType - Set type of subcommunicators.
118d8a68f86SHong Zhang 
119d8a68f86SHong Zhang    Collective on MPI_Comm
120d8a68f86SHong Zhang 
121d8a68f86SHong Zhang    Input Parameter:
122d8a68f86SHong Zhang +  psubcomm - PetscSubcomm context
1231ba920a7SHong Zhang -  subcommtype - subcommunicator type, PETSC_SUBCOMM_CONTIGUOUS,PETSC_SUBCOMM_INTERLACED
124d8a68f86SHong Zhang 
125d8a68f86SHong Zhang    Level: advanced
126d8a68f86SHong Zhang 
127d8a68f86SHong Zhang .keywords: communicator
128d8a68f86SHong Zhang 
129d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral()
130d8a68f86SHong Zhang @*/
1317c764164SBarry Smith PetscErrorCode  PetscSubcommSetType(PetscSubcomm psubcomm,PetscSubcommType subcommtype)
132d8a68f86SHong Zhang {
133d8a68f86SHong Zhang   PetscErrorCode ierr;
134d8a68f86SHong Zhang 
135d8a68f86SHong Zhang   PetscFunctionBegin;
136d8a68f86SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
13745487dadSJed Brown   if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()",psubcomm->n);
138d8a68f86SHong Zhang 
139d8a68f86SHong Zhang   if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS) {
140d8a68f86SHong Zhang     ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr);
141d8a68f86SHong Zhang   } else if (subcommtype == PETSC_SUBCOMM_INTERLACED) {
142d8a68f86SHong Zhang     ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr);
14345487dadSJed Brown   } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %s is not supported yet",PetscSubcommTypes[subcommtype]);
144d8a68f86SHong Zhang   PetscFunctionReturn(0);
145d8a68f86SHong Zhang }
146d8a68f86SHong Zhang 
147d8a68f86SHong Zhang #undef __FUNCT__
148d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetTypeGeneral"
1491ba920a7SHong Zhang /*@C
15065d9b8f1SHong Zhang   PetscSubcommSetTypeGeneral - Set a PetscSubcomm from user's specifications
1511ba920a7SHong Zhang 
1521ba920a7SHong Zhang    Collective on MPI_Comm
1531ba920a7SHong Zhang 
1541ba920a7SHong Zhang    Input Parameter:
1551ba920a7SHong Zhang +  psubcomm - PetscSubcomm context
1561ba920a7SHong Zhang .  color   - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator.
15765d9b8f1SHong Zhang -  subrank - rank in the subcommunicator
1581ba920a7SHong Zhang 
1591ba920a7SHong Zhang    Level: advanced
1601ba920a7SHong Zhang 
1611ba920a7SHong Zhang .keywords: communicator, create
1621ba920a7SHong Zhang 
1631ba920a7SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType()
1641ba920a7SHong Zhang @*/
16565d9b8f1SHong Zhang PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank)
166d8a68f86SHong Zhang {
1671ba920a7SHong Zhang   PetscErrorCode ierr;
1681ba920a7SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
169c9e2ceb8SHong Zhang   PetscMPIInt    size,icolor,duprank,*recvbuf,sendbuf[3],mysubsize,rank,*subsize;
17045487dadSJed Brown   PetscMPIInt    i,nsubcomm=psubcomm->n;
1711ba920a7SHong Zhang 
172d8a68f86SHong Zhang   PetscFunctionBegin;
1731ba920a7SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
17445487dadSJed Brown   if (nsubcomm < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()",nsubcomm);
1751ba920a7SHong Zhang 
1761ba920a7SHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
1771ba920a7SHong Zhang 
17865d9b8f1SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
1799abe469cSDmitry 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(). */
1801ba920a7SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
181785e854fSJed Brown   ierr = PetscMalloc1(2*size,&recvbuf);CHKERRQ(ierr);
18265d9b8f1SHong Zhang 
18365d9b8f1SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
184c9e2ceb8SHong Zhang   ierr = MPI_Comm_size(subcomm,&mysubsize);CHKERRQ(ierr);
18565d9b8f1SHong Zhang 
18665d9b8f1SHong Zhang   sendbuf[0] = color;
187c9e2ceb8SHong Zhang   sendbuf[1] = mysubsize;
188c9e2ceb8SHong Zhang   ierr = MPI_Allgather(sendbuf,2,MPI_INT,recvbuf,2,MPI_INT,comm);CHKERRQ(ierr);
18965d9b8f1SHong Zhang 
190eca4a452SHong Zhang   ierr = PetscCalloc1(nsubcomm,&subsize);CHKERRQ(ierr);
1919972d2ceSHong Zhang   for (i=0; i<2*size; i+=2) {
192c9e2ceb8SHong Zhang     subsize[recvbuf[i]] = recvbuf[i+1];
1931ba920a7SHong Zhang   }
19465d9b8f1SHong Zhang   ierr = PetscFree(recvbuf);CHKERRQ(ierr);
19565d9b8f1SHong Zhang 
19665d9b8f1SHong Zhang   duprank = 0;
197c9e2ceb8SHong Zhang   for (icolor=0; icolor<nsubcomm; icolor++) {
19865d9b8f1SHong Zhang     if (icolor != color) { /* not color of this process */
199c9e2ceb8SHong Zhang       duprank += subsize[icolor];
20065d9b8f1SHong Zhang     } else {
20165d9b8f1SHong Zhang       duprank += subrank;
20265d9b8f1SHong Zhang       break;
20365d9b8f1SHong Zhang     }
20465d9b8f1SHong Zhang   }
20565d9b8f1SHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
20665d9b8f1SHong Zhang 
2070298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
208306c2d5bSBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr);
209b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
210b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
211a297a907SKarl Rupp 
2121ba920a7SHong Zhang   psubcomm->color   = color;
213c9e2ceb8SHong Zhang   psubcomm->subsize = subsize;
214c9e2ceb8SHong Zhang   psubcomm->type    = PETSC_SUBCOMM_GENERAL;
215d8a68f86SHong Zhang   PetscFunctionReturn(0);
216d8a68f86SHong Zhang }
217638faf0bSHong Zhang 
218cd05a4c0SHong Zhang #undef __FUNCT__
219cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy"
2206bf464f9SBarry Smith PetscErrorCode  PetscSubcommDestroy(PetscSubcomm *psubcomm)
221cd05a4c0SHong Zhang {
222cd05a4c0SHong Zhang   PetscErrorCode ierr;
223cd05a4c0SHong Zhang 
224cd05a4c0SHong Zhang   PetscFunctionBegin;
2256bf464f9SBarry Smith   if (!*psubcomm) PetscFunctionReturn(0);
226aa9c1079SBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr);
227306c2d5bSBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->child);CHKERRQ(ierr);
228e37c6257SHong Zhang   ierr = PetscFree((*psubcomm)->subsize);CHKERRQ(ierr);
2296bf464f9SBarry Smith   ierr = PetscFree((*psubcomm));CHKERRQ(ierr);
230cd05a4c0SHong Zhang   PetscFunctionReturn(0);
231cd05a4c0SHong Zhang }
232cd05a4c0SHong Zhang 
233cd05a4c0SHong Zhang #undef __FUNCT__
234cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate"
235ab8c242fSMatthew Knepley /*@C
236cd05a4c0SHong Zhang   PetscSubcommCreate - Create a PetscSubcomm context.
237cd05a4c0SHong Zhang 
238cd05a4c0SHong Zhang    Collective on MPI_Comm
239cd05a4c0SHong Zhang 
240cd05a4c0SHong Zhang    Input Parameter:
2419873d53eSJed Brown .  comm - MPI communicator
242cd05a4c0SHong Zhang 
243cd05a4c0SHong Zhang    Output Parameter:
244cd05a4c0SHong Zhang .  psubcomm - location to store the PetscSubcomm context
245cd05a4c0SHong Zhang 
246638faf0bSHong Zhang    Level: advanced
247cd05a4c0SHong Zhang 
248638faf0bSHong Zhang .keywords: communicator, create
249638faf0bSHong Zhang 
250638faf0bSHong Zhang .seealso: PetscSubcommDestroy()
251638faf0bSHong Zhang @*/
2527087cfbeSBarry Smith PetscErrorCode  PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm)
253638faf0bSHong Zhang {
254638faf0bSHong Zhang   PetscErrorCode ierr;
255d3b23db5SHong Zhang   PetscMPIInt    rank,size;
256638faf0bSHong Zhang 
257638faf0bSHong Zhang   PetscFunctionBegin;
258b00a9115SJed Brown   ierr = PetscNew(psubcomm);CHKERRQ(ierr);
259a297a907SKarl Rupp 
260d3b23db5SHong Zhang   /* set defaults */
261d3b23db5SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
262d3b23db5SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
263f68be91cSHong Zhang 
264d8a68f86SHong Zhang   (*psubcomm)->parent    = comm;
265d3b23db5SHong Zhang   (*psubcomm)->dupparent = comm;
266306c2d5bSBarry Smith   (*psubcomm)->child     = PETSC_COMM_SELF;
267d3b23db5SHong Zhang   (*psubcomm)->n         = size;
268d3b23db5SHong Zhang   (*psubcomm)->color     = rank;
269e37c6257SHong Zhang   (*psubcomm)->subsize   = NULL;
270d3b23db5SHong Zhang   (*psubcomm)->type      = PETSC_SUBCOMM_INTERLACED;
271638faf0bSHong Zhang   PetscFunctionReturn(0);
272638faf0bSHong Zhang }
273638faf0bSHong Zhang 
274638faf0bSHong Zhang #undef __FUNCT__
27553c77d0aSJed Brown #define __FUNCT__ "PetscSubcommCreate_contiguous"
276d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm)
277638faf0bSHong Zhang {
278638faf0bSHong Zhang   PetscErrorCode ierr;
279d6037b41SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank=-1,subrank=-1;
28045487dadSJed Brown   PetscMPIInt    np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n;
281d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
282638faf0bSHong Zhang 
283638faf0bSHong Zhang   PetscFunctionBegin;
28455e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
28555e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
28655e3b8d2SHong Zhang 
287638faf0bSHong Zhang   /* get size of each subcommunicator */
288854ce69bSBarry Smith   ierr = PetscMalloc1(1+nsubcomm,&subsize);CHKERRQ(ierr);
289a297a907SKarl Rupp 
290638faf0bSHong Zhang   np_subcomm = size/nsubcomm;
291638faf0bSHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
292638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++) {
293638faf0bSHong Zhang     subsize[i] = np_subcomm;
294638faf0bSHong Zhang     if (i<nleftover) subsize[i]++;
295638faf0bSHong Zhang   }
296638faf0bSHong Zhang 
297638faf0bSHong Zhang   /* get color and subrank of this proc */
298638faf0bSHong Zhang   rankstart = 0;
299638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++) {
300638faf0bSHong Zhang     if (rank >= rankstart && rank < rankstart+subsize[i]) {
301638faf0bSHong Zhang       color   = i;
302638faf0bSHong Zhang       subrank = rank - rankstart;
303638faf0bSHong Zhang       duprank = rank;
304638faf0bSHong Zhang       break;
305a297a907SKarl Rupp     } else rankstart += subsize[i];
306638faf0bSHong Zhang   }
307638faf0bSHong Zhang 
308638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
309638faf0bSHong Zhang 
310638faf0bSHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
311638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
3120298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
313306c2d5bSBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr);
314b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
315b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
316a297a907SKarl Rupp 
317d8a68f86SHong Zhang   psubcomm->color   = color;
318e37c6257SHong Zhang   psubcomm->subsize = subsize;
319f38d543fSHong Zhang   psubcomm->type    = PETSC_SUBCOMM_CONTIGUOUS;
320638faf0bSHong Zhang   PetscFunctionReturn(0);
321638faf0bSHong Zhang }
322638faf0bSHong Zhang 
323638faf0bSHong Zhang #undef __FUNCT__
324638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced"
325638faf0bSHong Zhang /*
326638faf0bSHong Zhang    Note:
327638faf0bSHong Zhang    In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators
32845fc02eaSBarry Smith    by iteratively taking a process into a subcommunicator.
329cd05a4c0SHong Zhang    Example: size=4, nsubcomm=(*psubcomm)->n=3
330cd05a4c0SHong Zhang      comm=(*psubcomm)->parent:
331cd05a4c0SHong Zhang       rank:     [0]  [1]  [2]  [3]
332cd05a4c0SHong Zhang       color:     0    1    2    0
333cd05a4c0SHong Zhang 
334cd05a4c0SHong Zhang      subcomm=(*psubcomm)->comm:
335cd05a4c0SHong Zhang       subrank:  [0]  [0]  [0]  [1]
336cd05a4c0SHong Zhang 
337cd05a4c0SHong Zhang      dupcomm=(*psubcomm)->dupparent:
338cd05a4c0SHong Zhang       duprank:  [0]  [2]  [3]  [1]
339cd05a4c0SHong Zhang 
340cd05a4c0SHong Zhang      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
341cd05a4c0SHong Zhang            subcomm[color = 1] has subsize=1, owns process [1]
342cd05a4c0SHong Zhang            subcomm[color = 2] has subsize=1, owns process [2]
343cd05a4c0SHong Zhang            dupcomm has same number of processes as comm, and its duprank maps
344cd05a4c0SHong Zhang            processes in subcomm contiguously into a 1d array:
345cd05a4c0SHong Zhang             duprank: [0] [1]      [2]         [3]
346cd05a4c0SHong Zhang             rank:    [0] [3]      [1]         [2]
347cd05a4c0SHong Zhang                     subcomm[0] subcomm[1]  subcomm[2]
348638faf0bSHong Zhang */
349cd05a4c0SHong Zhang 
350d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm)
351cd05a4c0SHong Zhang {
352cd05a4c0SHong Zhang   PetscErrorCode ierr;
353cd05a4c0SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank,subrank;
35445487dadSJed Brown   PetscMPIInt    np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n;
355d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
356cd05a4c0SHong Zhang 
357cd05a4c0SHong Zhang   PetscFunctionBegin;
35855e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
35955e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
36055e3b8d2SHong Zhang 
361cd05a4c0SHong Zhang   /* get size of each subcommunicator */
362854ce69bSBarry Smith   ierr = PetscMalloc1(1+nsubcomm,&subsize);CHKERRQ(ierr);
363a297a907SKarl Rupp 
364cd05a4c0SHong Zhang   np_subcomm = size/nsubcomm;
365cd05a4c0SHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
366cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++) {
367cd05a4c0SHong Zhang     subsize[i] = np_subcomm;
368cd05a4c0SHong Zhang     if (i<nleftover) subsize[i]++;
369cd05a4c0SHong Zhang   }
370cd05a4c0SHong Zhang 
371cd05a4c0SHong Zhang   /* find color for this proc */
372cd05a4c0SHong Zhang   color   = rank%nsubcomm;
373cd05a4c0SHong Zhang   subrank = rank/nsubcomm;
374cd05a4c0SHong Zhang 
375cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
376cd05a4c0SHong Zhang 
377cd05a4c0SHong Zhang   j = 0; duprank = 0;
378cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++) {
379cd05a4c0SHong Zhang     if (j == color) {
380cd05a4c0SHong Zhang       duprank += subrank;
381cd05a4c0SHong Zhang       break;
382cd05a4c0SHong Zhang     }
383cd05a4c0SHong Zhang     duprank += subsize[i]; j++;
384cd05a4c0SHong Zhang   }
385cd05a4c0SHong Zhang 
386cd05a4c0SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
387cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
3880298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
389306c2d5bSBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr);
390b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
391b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
392a297a907SKarl Rupp 
393d8a68f86SHong Zhang   psubcomm->color   = color;
394e37c6257SHong Zhang   psubcomm->subsize = subsize;
395f38d543fSHong Zhang   psubcomm->type    = PETSC_SUBCOMM_INTERLACED;
396cd05a4c0SHong Zhang   PetscFunctionReturn(0);
397cd05a4c0SHong Zhang }
398638faf0bSHong Zhang 
399638faf0bSHong Zhang 
400