xref: /petsc/src/sys/objects/subcomm.c (revision 9972d2cea69f576862c8d7ac556aa3a0a50ab130)
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);
13f68be91cSHong Zhang #undef __FUNCT__
14f68be91cSHong Zhang #define __FUNCT__ "PetscSubcommSetFromOptions"
15f68be91cSHong Zhang PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm psubcomm)
16f68be91cSHong Zhang {
17f68be91cSHong Zhang   PetscErrorCode ierr;
1845487dadSJed Brown   PetscSubcommType type;
19f68be91cSHong Zhang   PetscBool      flg;
20f68be91cSHong Zhang 
21f68be91cSHong Zhang   PetscFunctionBegin;
22f68be91cSHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Must call PetscSubcommCreate firt");
2345487dadSJed Brown   type = psubcomm->type;
2445487dadSJed Brown   ierr = PetscOptionsEnum("-psubcomm_type","PETSc subcommunicator","PetscSubcommSetType",PetscSubcommTypes,(PetscEnum)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);
28f68be91cSHong Zhang     ierr = PetscCommDestroy(&(psubcomm)->comm);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 
4419171117SHong Zhang   ierr = PetscOptionsHasName(NULL, "-psubcomm_view", &flg);CHKERRQ(ierr);
4519171117SHong Zhang   if (flg) {
4619171117SHong Zhang     ierr = PetscSubcommView(psubcomm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
4719171117SHong Zhang   }
48f68be91cSHong Zhang   PetscFunctionReturn(0);
49f68be91cSHong Zhang }
50d8a68f86SHong Zhang 
51d8a68f86SHong Zhang #undef __FUNCT__
52053d1c95SHong Zhang #define __FUNCT__ "PetscSubcommView"
53053d1c95SHong Zhang PetscErrorCode PetscSubcommView(PetscSubcomm psubcomm,PetscViewer viewer)
54053d1c95SHong Zhang {
55053d1c95SHong Zhang   PetscErrorCode    ierr;
56053d1c95SHong Zhang   PetscBool         iascii;
57053d1c95SHong Zhang   PetscViewerFormat format;
58053d1c95SHong Zhang 
59053d1c95SHong Zhang   PetscFunctionBegin;
60053d1c95SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
61053d1c95SHong Zhang   if (iascii) {
62053d1c95SHong Zhang     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
63053d1c95SHong Zhang     if (format == PETSC_VIEWER_DEFAULT) {
64053d1c95SHong Zhang       MPI_Comm    comm=psubcomm->parent;
65053d1c95SHong Zhang       PetscMPIInt rank,size,subsize,subrank,duprank;
66053d1c95SHong Zhang 
67053d1c95SHong Zhang       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
6845487dadSJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"PetscSubcomm type %s with total %d MPI processes:\n",PetscSubcommTypes[psubcomm->type],size);CHKERRQ(ierr);
69053d1c95SHong Zhang       ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
70053d1c95SHong Zhang       ierr = MPI_Comm_size(psubcomm->comm,&subsize);CHKERRQ(ierr);
71053d1c95SHong Zhang       ierr = MPI_Comm_rank(psubcomm->comm,&subrank);CHKERRQ(ierr);
72053d1c95SHong Zhang       ierr = MPI_Comm_rank(psubcomm->dupparent,&duprank);CHKERRQ(ierr);
7345487dadSJed Brown       ierr = PetscSynchronizedPrintf(comm,"  [%d], color %d, sub-size %d, sub-rank %d, duprank %d\n",rank,psubcomm->color,subsize,subrank,duprank);
74053d1c95SHong Zhang       ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr);
75053d1c95SHong Zhang     }
76053d1c95SHong Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported yet");
77053d1c95SHong Zhang   PetscFunctionReturn(0);
78053d1c95SHong Zhang }
79053d1c95SHong Zhang 
80053d1c95SHong Zhang #undef __FUNCT__
81d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetNumber"
82d8a68f86SHong Zhang /*@C
83d8a68f86SHong Zhang   PetscSubcommSetNumber - Set total number of subcommunicators.
84d8a68f86SHong Zhang 
85d8a68f86SHong Zhang    Collective on MPI_Comm
86d8a68f86SHong Zhang 
87d8a68f86SHong Zhang    Input Parameter:
88d8a68f86SHong Zhang +  psubcomm - PetscSubcomm context
89d8a68f86SHong Zhang -  nsubcomm - the total number of subcommunicators in psubcomm
90d8a68f86SHong Zhang 
91d8a68f86SHong Zhang    Level: advanced
92d8a68f86SHong Zhang 
93d8a68f86SHong Zhang .keywords: communicator
94d8a68f86SHong Zhang 
95d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetType(),PetscSubcommSetTypeGeneral()
96d8a68f86SHong Zhang @*/
977087cfbeSBarry Smith PetscErrorCode  PetscSubcommSetNumber(PetscSubcomm psubcomm,PetscInt nsubcomm)
98d8a68f86SHong Zhang {
99d8a68f86SHong Zhang   PetscErrorCode ierr;
100d8a68f86SHong Zhang   MPI_Comm       comm=psubcomm->parent;
10145487dadSJed Brown   PetscMPIInt    msub,size;
102d8a68f86SHong Zhang 
103d8a68f86SHong Zhang   PetscFunctionBegin;
104d8a68f86SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate() first");
105d8a68f86SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
10645487dadSJed Brown   ierr = PetscMPIIntCast(nsubcomm,&msub);CHKERRQ(ierr);
10745487dadSJed 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);
108d8a68f86SHong Zhang 
10945487dadSJed Brown   psubcomm->n = msub;
110d8a68f86SHong Zhang   PetscFunctionReturn(0);
111d8a68f86SHong Zhang }
112d8a68f86SHong Zhang 
113d8a68f86SHong Zhang #undef __FUNCT__
114d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetType"
115d8a68f86SHong Zhang /*@C
116d8a68f86SHong Zhang   PetscSubcommSetType - Set type of subcommunicators.
117d8a68f86SHong Zhang 
118d8a68f86SHong Zhang    Collective on MPI_Comm
119d8a68f86SHong Zhang 
120d8a68f86SHong Zhang    Input Parameter:
121d8a68f86SHong Zhang +  psubcomm - PetscSubcomm context
1221ba920a7SHong Zhang -  subcommtype - subcommunicator type, PETSC_SUBCOMM_CONTIGUOUS,PETSC_SUBCOMM_INTERLACED
123d8a68f86SHong Zhang 
124d8a68f86SHong Zhang    Level: advanced
125d8a68f86SHong Zhang 
126d8a68f86SHong Zhang .keywords: communicator
127d8a68f86SHong Zhang 
128d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral()
129d8a68f86SHong Zhang @*/
1307c764164SBarry Smith PetscErrorCode  PetscSubcommSetType(PetscSubcomm psubcomm,PetscSubcommType subcommtype)
131d8a68f86SHong Zhang {
132d8a68f86SHong Zhang   PetscErrorCode ierr;
133d8a68f86SHong Zhang 
134d8a68f86SHong Zhang   PetscFunctionBegin;
135d8a68f86SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
13645487dadSJed Brown   if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()",psubcomm->n);
137d8a68f86SHong Zhang 
138d8a68f86SHong Zhang   if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS) {
139d8a68f86SHong Zhang     ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr);
140d8a68f86SHong Zhang   } else if (subcommtype == PETSC_SUBCOMM_INTERLACED) {
141d8a68f86SHong Zhang     ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr);
14245487dadSJed Brown   } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %s is not supported yet",PetscSubcommTypes[subcommtype]);
143d8a68f86SHong Zhang   PetscFunctionReturn(0);
144d8a68f86SHong Zhang }
145d8a68f86SHong Zhang 
146d8a68f86SHong Zhang #undef __FUNCT__
147d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetTypeGeneral"
1481ba920a7SHong Zhang /*@C
14965d9b8f1SHong Zhang   PetscSubcommSetTypeGeneral - Set a PetscSubcomm from user's specifications
1501ba920a7SHong Zhang 
1511ba920a7SHong Zhang    Collective on MPI_Comm
1521ba920a7SHong Zhang 
1531ba920a7SHong Zhang    Input Parameter:
1541ba920a7SHong Zhang +  psubcomm - PetscSubcomm context
1551ba920a7SHong Zhang .  color   - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator.
15665d9b8f1SHong Zhang -  subrank - rank in the subcommunicator
1571ba920a7SHong Zhang 
1581ba920a7SHong Zhang    Level: advanced
1591ba920a7SHong Zhang 
1601ba920a7SHong Zhang .keywords: communicator, create
1611ba920a7SHong Zhang 
1621ba920a7SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType()
1631ba920a7SHong Zhang @*/
16465d9b8f1SHong Zhang PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank)
165d8a68f86SHong Zhang {
1661ba920a7SHong Zhang   PetscErrorCode ierr;
1671ba920a7SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
168c9e2ceb8SHong Zhang   PetscMPIInt    size,icolor,duprank,*recvbuf,sendbuf[3],mysubsize,rank,*subsize;
16945487dadSJed Brown   PetscMPIInt    i,nsubcomm=psubcomm->n;
1701ba920a7SHong Zhang 
171d8a68f86SHong Zhang   PetscFunctionBegin;
1721ba920a7SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
17345487dadSJed Brown   if (nsubcomm < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()",nsubcomm);
1741ba920a7SHong Zhang 
1751ba920a7SHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
1761ba920a7SHong Zhang 
17765d9b8f1SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
1781ba920a7SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
179c9e2ceb8SHong Zhang   ierr = PetscMalloc(2*size*sizeof(PetscMPIInt),&recvbuf);CHKERRQ(ierr);
18065d9b8f1SHong Zhang 
18165d9b8f1SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
182c9e2ceb8SHong Zhang   ierr = MPI_Comm_size(subcomm,&mysubsize);CHKERRQ(ierr);
18365d9b8f1SHong Zhang 
18465d9b8f1SHong Zhang   sendbuf[0] = color;
185c9e2ceb8SHong Zhang   sendbuf[1] = mysubsize;
186c9e2ceb8SHong Zhang   ierr = MPI_Allgather(sendbuf,2,MPI_INT,recvbuf,2,MPI_INT,comm);CHKERRQ(ierr);
18765d9b8f1SHong Zhang 
188c9e2ceb8SHong Zhang   ierr = PetscMalloc(nsubcomm*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
189*9972d2ceSHong Zhang   for (i=0; i<2*size; i+=2) {
190c9e2ceb8SHong Zhang     subsize[recvbuf[i]] = recvbuf[i+1];
1911ba920a7SHong Zhang   }
19265d9b8f1SHong Zhang   ierr = PetscFree(recvbuf);CHKERRQ(ierr);
19365d9b8f1SHong Zhang 
19465d9b8f1SHong Zhang   duprank = 0;
195c9e2ceb8SHong Zhang   for (icolor=0; icolor<nsubcomm; icolor++) {
19665d9b8f1SHong Zhang     if (icolor != color) { /* not color of this process */
197c9e2ceb8SHong Zhang       duprank += subsize[icolor];
19865d9b8f1SHong Zhang     } else {
19965d9b8f1SHong Zhang       duprank += subrank;
20065d9b8f1SHong Zhang       break;
20165d9b8f1SHong Zhang     }
20265d9b8f1SHong Zhang   }
20365d9b8f1SHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
20465d9b8f1SHong Zhang 
2050298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
2060298fd71SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr);
207b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
208b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
209a297a907SKarl Rupp 
2101ba920a7SHong Zhang   psubcomm->color   = color;
211c9e2ceb8SHong Zhang   psubcomm->subsize = subsize;
212c9e2ceb8SHong Zhang   psubcomm->type    = PETSC_SUBCOMM_GENERAL;
213d8a68f86SHong Zhang   PetscFunctionReturn(0);
214d8a68f86SHong Zhang }
215638faf0bSHong Zhang 
216cd05a4c0SHong Zhang #undef __FUNCT__
217cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy"
2186bf464f9SBarry Smith PetscErrorCode  PetscSubcommDestroy(PetscSubcomm *psubcomm)
219cd05a4c0SHong Zhang {
220cd05a4c0SHong Zhang   PetscErrorCode ierr;
221cd05a4c0SHong Zhang 
222cd05a4c0SHong Zhang   PetscFunctionBegin;
2236bf464f9SBarry Smith   if (!*psubcomm) PetscFunctionReturn(0);
224aa9c1079SBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr);
225aa9c1079SBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->comm);CHKERRQ(ierr);
226e37c6257SHong Zhang   ierr = PetscFree((*psubcomm)->subsize);CHKERRQ(ierr);
2276bf464f9SBarry Smith   ierr = PetscFree((*psubcomm));CHKERRQ(ierr);
228cd05a4c0SHong Zhang   PetscFunctionReturn(0);
229cd05a4c0SHong Zhang }
230cd05a4c0SHong Zhang 
231cd05a4c0SHong Zhang #undef __FUNCT__
232cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate"
233ab8c242fSMatthew Knepley /*@C
234cd05a4c0SHong Zhang   PetscSubcommCreate - Create a PetscSubcomm context.
235cd05a4c0SHong Zhang 
236cd05a4c0SHong Zhang    Collective on MPI_Comm
237cd05a4c0SHong Zhang 
238cd05a4c0SHong Zhang    Input Parameter:
2399873d53eSJed Brown .  comm - MPI communicator
240cd05a4c0SHong Zhang 
241cd05a4c0SHong Zhang    Output Parameter:
242cd05a4c0SHong Zhang .  psubcomm - location to store the PetscSubcomm context
243cd05a4c0SHong Zhang 
244638faf0bSHong Zhang    Level: advanced
245cd05a4c0SHong Zhang 
246638faf0bSHong Zhang .keywords: communicator, create
247638faf0bSHong Zhang 
248638faf0bSHong Zhang .seealso: PetscSubcommDestroy()
249638faf0bSHong Zhang @*/
2507087cfbeSBarry Smith PetscErrorCode  PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm)
251638faf0bSHong Zhang {
252638faf0bSHong Zhang   PetscErrorCode ierr;
253d3b23db5SHong Zhang   PetscMPIInt    rank,size;
254638faf0bSHong Zhang 
255638faf0bSHong Zhang   PetscFunctionBegin;
256d8a68f86SHong Zhang   ierr = PetscNew(struct _n_PetscSubcomm,psubcomm);CHKERRQ(ierr);
257a297a907SKarl Rupp 
258d3b23db5SHong Zhang   /* set defaults */
259d3b23db5SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
260d3b23db5SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
261f68be91cSHong Zhang 
262d8a68f86SHong Zhang   (*psubcomm)->parent    = comm;
263d3b23db5SHong Zhang   (*psubcomm)->dupparent = comm;
264d3b23db5SHong Zhang   (*psubcomm)->comm      = PETSC_COMM_SELF;
265d3b23db5SHong Zhang   (*psubcomm)->n         = size;
266d3b23db5SHong Zhang   (*psubcomm)->color     = rank;
267e37c6257SHong Zhang   (*psubcomm)->subsize   = NULL;
268d3b23db5SHong Zhang   (*psubcomm)->type      = PETSC_SUBCOMM_INTERLACED;
269638faf0bSHong Zhang   PetscFunctionReturn(0);
270638faf0bSHong Zhang }
271638faf0bSHong Zhang 
272638faf0bSHong Zhang #undef __FUNCT__
27353c77d0aSJed Brown #define __FUNCT__ "PetscSubcommCreate_contiguous"
274d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm)
275638faf0bSHong Zhang {
276638faf0bSHong Zhang   PetscErrorCode ierr;
277d6037b41SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank=-1,subrank=-1;
27845487dadSJed Brown   PetscMPIInt    np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n;
279d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
280638faf0bSHong Zhang 
281638faf0bSHong Zhang   PetscFunctionBegin;
28255e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
28355e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
28455e3b8d2SHong Zhang 
285638faf0bSHong Zhang   /* get size of each subcommunicator */
286638faf0bSHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
287a297a907SKarl Rupp 
288638faf0bSHong Zhang   np_subcomm = size/nsubcomm;
289638faf0bSHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
290638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++) {
291638faf0bSHong Zhang     subsize[i] = np_subcomm;
292638faf0bSHong Zhang     if (i<nleftover) subsize[i]++;
293638faf0bSHong Zhang   }
294638faf0bSHong Zhang 
295638faf0bSHong Zhang   /* get color and subrank of this proc */
296638faf0bSHong Zhang   rankstart = 0;
297638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++) {
298638faf0bSHong Zhang     if (rank >= rankstart && rank < rankstart+subsize[i]) {
299638faf0bSHong Zhang       color   = i;
300638faf0bSHong Zhang       subrank = rank - rankstart;
301638faf0bSHong Zhang       duprank = rank;
302638faf0bSHong Zhang       break;
303a297a907SKarl Rupp     } else rankstart += subsize[i];
304638faf0bSHong Zhang   }
305638faf0bSHong Zhang 
306638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
307638faf0bSHong Zhang 
308638faf0bSHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
309638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
31052af3f7aSShri Abhyankar   {
31152af3f7aSShri Abhyankar     PetscThreadComm tcomm;
31252af3f7aSShri Abhyankar     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
31352af3f7aSShri Abhyankar     ierr = MPI_Attr_put(dupcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
31452af3f7aSShri Abhyankar     tcomm->refct++;
31552af3f7aSShri Abhyankar     ierr = MPI_Attr_put(subcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
31652af3f7aSShri Abhyankar     tcomm->refct++;
31752af3f7aSShri Abhyankar   }
3180298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
3190298fd71SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr);
320b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
321b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
322a297a907SKarl Rupp 
323d8a68f86SHong Zhang   psubcomm->color   = color;
324e37c6257SHong Zhang   psubcomm->subsize = subsize;
325f38d543fSHong Zhang   psubcomm->type    = PETSC_SUBCOMM_CONTIGUOUS;
326638faf0bSHong Zhang   PetscFunctionReturn(0);
327638faf0bSHong Zhang }
328638faf0bSHong Zhang 
329638faf0bSHong Zhang #undef __FUNCT__
330638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced"
331638faf0bSHong Zhang /*
332638faf0bSHong Zhang    Note:
333638faf0bSHong Zhang    In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators
33445fc02eaSBarry Smith    by iteratively taking a process into a subcommunicator.
335cd05a4c0SHong Zhang    Example: size=4, nsubcomm=(*psubcomm)->n=3
336cd05a4c0SHong Zhang      comm=(*psubcomm)->parent:
337cd05a4c0SHong Zhang       rank:     [0]  [1]  [2]  [3]
338cd05a4c0SHong Zhang       color:     0    1    2    0
339cd05a4c0SHong Zhang 
340cd05a4c0SHong Zhang      subcomm=(*psubcomm)->comm:
341cd05a4c0SHong Zhang       subrank:  [0]  [0]  [0]  [1]
342cd05a4c0SHong Zhang 
343cd05a4c0SHong Zhang      dupcomm=(*psubcomm)->dupparent:
344cd05a4c0SHong Zhang       duprank:  [0]  [2]  [3]  [1]
345cd05a4c0SHong Zhang 
346cd05a4c0SHong Zhang      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
347cd05a4c0SHong Zhang            subcomm[color = 1] has subsize=1, owns process [1]
348cd05a4c0SHong Zhang            subcomm[color = 2] has subsize=1, owns process [2]
349cd05a4c0SHong Zhang            dupcomm has same number of processes as comm, and its duprank maps
350cd05a4c0SHong Zhang            processes in subcomm contiguously into a 1d array:
351cd05a4c0SHong Zhang             duprank: [0] [1]      [2]         [3]
352cd05a4c0SHong Zhang             rank:    [0] [3]      [1]         [2]
353cd05a4c0SHong Zhang                     subcomm[0] subcomm[1]  subcomm[2]
354638faf0bSHong Zhang */
355cd05a4c0SHong Zhang 
356d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm)
357cd05a4c0SHong Zhang {
358cd05a4c0SHong Zhang   PetscErrorCode ierr;
359cd05a4c0SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank,subrank;
36045487dadSJed Brown   PetscMPIInt    np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n;
361d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
362cd05a4c0SHong Zhang 
363cd05a4c0SHong Zhang   PetscFunctionBegin;
36455e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
36555e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
36655e3b8d2SHong Zhang 
367cd05a4c0SHong Zhang   /* get size of each subcommunicator */
368cd05a4c0SHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
369a297a907SKarl Rupp 
370cd05a4c0SHong Zhang   np_subcomm = size/nsubcomm;
371cd05a4c0SHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
372cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++) {
373cd05a4c0SHong Zhang     subsize[i] = np_subcomm;
374cd05a4c0SHong Zhang     if (i<nleftover) subsize[i]++;
375cd05a4c0SHong Zhang   }
376cd05a4c0SHong Zhang 
377cd05a4c0SHong Zhang   /* find color for this proc */
378cd05a4c0SHong Zhang   color   = rank%nsubcomm;
379cd05a4c0SHong Zhang   subrank = rank/nsubcomm;
380cd05a4c0SHong Zhang 
381cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
382cd05a4c0SHong Zhang 
383cd05a4c0SHong Zhang   j = 0; duprank = 0;
384cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++) {
385cd05a4c0SHong Zhang     if (j == color) {
386cd05a4c0SHong Zhang       duprank += subrank;
387cd05a4c0SHong Zhang       break;
388cd05a4c0SHong Zhang     }
389cd05a4c0SHong Zhang     duprank += subsize[i]; j++;
390cd05a4c0SHong Zhang   }
391cd05a4c0SHong Zhang 
392cd05a4c0SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
393cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
39452af3f7aSShri Abhyankar   {
39552af3f7aSShri Abhyankar     PetscThreadComm tcomm;
39652af3f7aSShri Abhyankar     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
39752af3f7aSShri Abhyankar     ierr = MPI_Attr_put(dupcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
39852af3f7aSShri Abhyankar     tcomm->refct++;
39952af3f7aSShri Abhyankar     ierr = MPI_Attr_put(subcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
40052af3f7aSShri Abhyankar     tcomm->refct++;
40152af3f7aSShri Abhyankar   }
4020298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
4030298fd71SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr);
404b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
405b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
406a297a907SKarl Rupp 
407d8a68f86SHong Zhang   psubcomm->color   = color;
408e37c6257SHong Zhang   psubcomm->subsize = subsize;
409f38d543fSHong Zhang   psubcomm->type    = PETSC_SUBCOMM_INTERLACED;
410cd05a4c0SHong Zhang   PetscFunctionReturn(0);
411cd05a4c0SHong Zhang }
412638faf0bSHong Zhang 
413638faf0bSHong Zhang 
414