xref: /petsc/src/sys/objects/subcomm.c (revision 65d9b8f1aa0e7fa6d9b7f735e1aa03b2ff0c056d)
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;
18f68be91cSHong Zhang   PetscInt       type;
19f68be91cSHong Zhang   const char     *psubcommTypes[3] = {"general","contiguous","interlaced"};
20f68be91cSHong Zhang   PetscBool      flg;
21f68be91cSHong Zhang 
22f68be91cSHong Zhang   PetscFunctionBegin;
23f68be91cSHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Must call PetscSubcommCreate firt");
24f68be91cSHong Zhang   ierr = PetscOptionsEList("-psubcomm_type","PETSc subcommunicator","PetscSubcommSetType",psubcommTypes,3,psubcommTypes[2],&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) {
31f68be91cSHong Zhang     case 1:
32f68be91cSHong Zhang       ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr);
33f68be91cSHong Zhang       break;
34f68be91cSHong Zhang     case 2:
35f68be91cSHong Zhang       ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr);
36f68be91cSHong Zhang       break;
37f68be91cSHong Zhang     default:
38f68be91cSHong Zhang       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",type);
39f68be91cSHong Zhang     }
40f68be91cSHong Zhang   }
4119171117SHong Zhang 
4219171117SHong Zhang   ierr = PetscOptionsHasName(NULL, "-psubcomm_view", &flg);CHKERRQ(ierr);
4319171117SHong Zhang   if (flg) {
4419171117SHong Zhang     ierr = PetscSubcommView(psubcomm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
4519171117SHong Zhang   }
46f68be91cSHong Zhang   PetscFunctionReturn(0);
47f68be91cSHong Zhang }
48d8a68f86SHong Zhang 
49d8a68f86SHong Zhang #undef __FUNCT__
50053d1c95SHong Zhang #define __FUNCT__ "PetscSubcommView"
51053d1c95SHong Zhang PetscErrorCode PetscSubcommView(PetscSubcomm psubcomm,PetscViewer viewer)
52053d1c95SHong Zhang {
53053d1c95SHong Zhang   PetscErrorCode    ierr;
54053d1c95SHong Zhang   PetscBool         iascii;
55053d1c95SHong Zhang   PetscViewerFormat format;
56053d1c95SHong Zhang 
57053d1c95SHong Zhang   PetscFunctionBegin;
58053d1c95SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
59053d1c95SHong Zhang   if (iascii) {
60053d1c95SHong Zhang     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
61053d1c95SHong Zhang     if (format == PETSC_VIEWER_DEFAULT) {
62053d1c95SHong Zhang       MPI_Comm    comm=psubcomm->parent;
63053d1c95SHong Zhang       PetscMPIInt rank,size,subsize,subrank,duprank;
64053d1c95SHong Zhang 
65053d1c95SHong Zhang       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
66053d1c95SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"PetscSubcomm type %s with total %D MPI processes:\n",PetscSubcommTypes[psubcomm->type],size);CHKERRQ(ierr);
67053d1c95SHong Zhang       ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
68053d1c95SHong Zhang       ierr = MPI_Comm_size(psubcomm->comm,&subsize);CHKERRQ(ierr);
69053d1c95SHong Zhang       ierr = MPI_Comm_rank(psubcomm->comm,&subrank);CHKERRQ(ierr);
70053d1c95SHong Zhang       ierr = MPI_Comm_rank(psubcomm->dupparent,&duprank);CHKERRQ(ierr);
71053d1c95SHong Zhang       ierr = PetscSynchronizedPrintf(comm,"  [%D], color %D, sub-size %D, sub-rank %D, duprank %D\n",rank,psubcomm->color,subsize,subrank,duprank);
72053d1c95SHong Zhang       ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr);
73053d1c95SHong Zhang     }
74053d1c95SHong Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported yet");
75053d1c95SHong Zhang   PetscFunctionReturn(0);
76053d1c95SHong Zhang }
77053d1c95SHong Zhang 
78053d1c95SHong Zhang #undef __FUNCT__
79d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetNumber"
80d8a68f86SHong Zhang /*@C
81d8a68f86SHong Zhang   PetscSubcommSetNumber - Set total number of subcommunicators.
82d8a68f86SHong Zhang 
83d8a68f86SHong Zhang    Collective on MPI_Comm
84d8a68f86SHong Zhang 
85d8a68f86SHong Zhang    Input Parameter:
86d8a68f86SHong Zhang +  psubcomm - PetscSubcomm context
87d8a68f86SHong Zhang -  nsubcomm - the total number of subcommunicators in psubcomm
88d8a68f86SHong Zhang 
89d8a68f86SHong Zhang    Level: advanced
90d8a68f86SHong Zhang 
91d8a68f86SHong Zhang .keywords: communicator
92d8a68f86SHong Zhang 
93d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetType(),PetscSubcommSetTypeGeneral()
94d8a68f86SHong Zhang @*/
957087cfbeSBarry Smith PetscErrorCode  PetscSubcommSetNumber(PetscSubcomm psubcomm,PetscInt nsubcomm)
96d8a68f86SHong Zhang {
97d8a68f86SHong Zhang   PetscErrorCode ierr;
98d8a68f86SHong Zhang   MPI_Comm       comm=psubcomm->parent;
99d8a68f86SHong Zhang   PetscMPIInt    rank,size;
100d8a68f86SHong Zhang 
101d8a68f86SHong Zhang   PetscFunctionBegin;
102d8a68f86SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate() first");
103d8a68f86SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
104d8a68f86SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
105d8a68f86SHong 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);
106d8a68f86SHong Zhang 
107d8a68f86SHong Zhang   psubcomm->n = nsubcomm;
108d8a68f86SHong Zhang   PetscFunctionReturn(0);
109d8a68f86SHong Zhang }
110d8a68f86SHong Zhang 
111d8a68f86SHong Zhang #undef __FUNCT__
112d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetType"
113d8a68f86SHong Zhang /*@C
114d8a68f86SHong Zhang   PetscSubcommSetType - Set type of subcommunicators.
115d8a68f86SHong Zhang 
116d8a68f86SHong Zhang    Collective on MPI_Comm
117d8a68f86SHong Zhang 
118d8a68f86SHong Zhang    Input Parameter:
119d8a68f86SHong Zhang +  psubcomm - PetscSubcomm context
1201ba920a7SHong Zhang -  subcommtype - subcommunicator type, PETSC_SUBCOMM_CONTIGUOUS,PETSC_SUBCOMM_INTERLACED
121d8a68f86SHong Zhang 
122d8a68f86SHong Zhang    Level: advanced
123d8a68f86SHong Zhang 
124d8a68f86SHong Zhang .keywords: communicator
125d8a68f86SHong Zhang 
126d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral()
127d8a68f86SHong Zhang @*/
1287c764164SBarry Smith PetscErrorCode  PetscSubcommSetType(PetscSubcomm psubcomm,PetscSubcommType subcommtype)
129d8a68f86SHong Zhang {
130d8a68f86SHong Zhang   PetscErrorCode ierr;
131d8a68f86SHong Zhang 
132d8a68f86SHong Zhang   PetscFunctionBegin;
133d8a68f86SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
134d8a68f86SHong Zhang   if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n);
135d8a68f86SHong Zhang 
136d8a68f86SHong Zhang   if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS) {
137d8a68f86SHong Zhang     ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr);
138d8a68f86SHong Zhang   } else if (subcommtype == PETSC_SUBCOMM_INTERLACED) {
139d8a68f86SHong Zhang     ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr);
1407c764164SBarry Smith   } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",subcommtype);
141d8a68f86SHong Zhang   PetscFunctionReturn(0);
142d8a68f86SHong Zhang }
143d8a68f86SHong Zhang 
144d8a68f86SHong Zhang #undef __FUNCT__
145d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetTypeGeneral"
1461ba920a7SHong Zhang /*@C
147*65d9b8f1SHong Zhang   PetscSubcommSetTypeGeneral - Set a PetscSubcomm from user's specifications
1481ba920a7SHong Zhang 
1491ba920a7SHong Zhang    Collective on MPI_Comm
1501ba920a7SHong Zhang 
1511ba920a7SHong Zhang    Input Parameter:
1521ba920a7SHong Zhang +  psubcomm - PetscSubcomm context
1531ba920a7SHong Zhang .  color   - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator.
154*65d9b8f1SHong Zhang -  subrank - rank in the subcommunicator
1551ba920a7SHong Zhang 
1561ba920a7SHong Zhang    Level: advanced
1571ba920a7SHong Zhang 
1581ba920a7SHong Zhang .keywords: communicator, create
1591ba920a7SHong Zhang 
1601ba920a7SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType()
1611ba920a7SHong Zhang @*/
162*65d9b8f1SHong Zhang PetscErrorCode  PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank)
163d8a68f86SHong Zhang {
1641ba920a7SHong Zhang   PetscErrorCode ierr;
1651ba920a7SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
166*65d9b8f1SHong Zhang   PetscMPIInt    size,icolor,duprank,*recvbuf,sendbuf[3],subsize,rank,subsizes[psubcomm->n];
167*65d9b8f1SHong Zhang   PetscInt       i;
1681ba920a7SHong Zhang 
169d8a68f86SHong Zhang   PetscFunctionBegin;
1701ba920a7SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
1711ba920a7SHong Zhang   if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n);
1721ba920a7SHong Zhang 
1731ba920a7SHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
1741ba920a7SHong Zhang 
175*65d9b8f1SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
1761ba920a7SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
177*65d9b8f1SHong Zhang   ierr = PetscMalloc(3*size*sizeof(PetscMPIInt),&recvbuf);CHKERRQ(ierr);
178*65d9b8f1SHong Zhang 
179*65d9b8f1SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
180*65d9b8f1SHong Zhang   ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr);
181*65d9b8f1SHong Zhang 
182*65d9b8f1SHong Zhang   sendbuf[0] = color;
183*65d9b8f1SHong Zhang   sendbuf[1] = subrank;
184*65d9b8f1SHong Zhang   sendbuf[2] = subsize;
185*65d9b8f1SHong Zhang   ierr = MPI_Allgather(sendbuf,3,MPI_INT,recvbuf,3,MPI_INT,comm);CHKERRQ(ierr);
186*65d9b8f1SHong Zhang 
187*65d9b8f1SHong Zhang   /* printf("[%d] (color,subrank,subsize): ",rank); */
188*65d9b8f1SHong Zhang   for (i=0; i<3*size; i+=3) {
189*65d9b8f1SHong Zhang     /* printf(" (%d %d %d);",recvbuf[i],recvbuf[i+1],recvbuf[i+2]); */
190*65d9b8f1SHong Zhang     subsizes[recvbuf[i]] = recvbuf[i+2];
1911ba920a7SHong Zhang   }
192*65d9b8f1SHong Zhang   ierr = PetscFree(recvbuf);CHKERRQ(ierr);
193*65d9b8f1SHong Zhang 
194*65d9b8f1SHong Zhang   duprank = 0;
195*65d9b8f1SHong Zhang   for (icolor=0; icolor<psubcomm->n; icolor++) {
196*65d9b8f1SHong Zhang     if (icolor != color) { /* not color of this process */
197*65d9b8f1SHong Zhang       duprank += subsizes[icolor];
198*65d9b8f1SHong Zhang     } else {
199*65d9b8f1SHong Zhang       duprank += subrank;
200*65d9b8f1SHong Zhang       break;
201*65d9b8f1SHong Zhang     }
202*65d9b8f1SHong Zhang   }
203*65d9b8f1SHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
204*65d9b8f1SHong 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;
211d8a68f86SHong Zhang   PetscFunctionReturn(0);
212d8a68f86SHong Zhang }
213638faf0bSHong Zhang 
214cd05a4c0SHong Zhang #undef __FUNCT__
215cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy"
2166bf464f9SBarry Smith PetscErrorCode  PetscSubcommDestroy(PetscSubcomm *psubcomm)
217cd05a4c0SHong Zhang {
218cd05a4c0SHong Zhang   PetscErrorCode ierr;
219cd05a4c0SHong Zhang 
220cd05a4c0SHong Zhang   PetscFunctionBegin;
2216bf464f9SBarry Smith   if (!*psubcomm) PetscFunctionReturn(0);
222aa9c1079SBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr);
223aa9c1079SBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->comm);CHKERRQ(ierr);
224e37c6257SHong Zhang   ierr = PetscFree((*psubcomm)->subsize);CHKERRQ(ierr);
2256bf464f9SBarry Smith   ierr = PetscFree((*psubcomm));CHKERRQ(ierr);
226cd05a4c0SHong Zhang   PetscFunctionReturn(0);
227cd05a4c0SHong Zhang }
228cd05a4c0SHong Zhang 
229cd05a4c0SHong Zhang #undef __FUNCT__
230cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate"
231ab8c242fSMatthew Knepley /*@C
232cd05a4c0SHong Zhang   PetscSubcommCreate - Create a PetscSubcomm context.
233cd05a4c0SHong Zhang 
234cd05a4c0SHong Zhang    Collective on MPI_Comm
235cd05a4c0SHong Zhang 
236cd05a4c0SHong Zhang    Input Parameter:
2379873d53eSJed Brown .  comm - MPI communicator
238cd05a4c0SHong Zhang 
239cd05a4c0SHong Zhang    Output Parameter:
240cd05a4c0SHong Zhang .  psubcomm - location to store the PetscSubcomm context
241cd05a4c0SHong Zhang 
242638faf0bSHong Zhang    Level: advanced
243cd05a4c0SHong Zhang 
244638faf0bSHong Zhang .keywords: communicator, create
245638faf0bSHong Zhang 
246638faf0bSHong Zhang .seealso: PetscSubcommDestroy()
247638faf0bSHong Zhang @*/
2487087cfbeSBarry Smith PetscErrorCode  PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm)
249638faf0bSHong Zhang {
250638faf0bSHong Zhang   PetscErrorCode ierr;
251d3b23db5SHong Zhang   PetscMPIInt    rank,size;
252638faf0bSHong Zhang 
253638faf0bSHong Zhang   PetscFunctionBegin;
254d8a68f86SHong Zhang   ierr = PetscNew(struct _n_PetscSubcomm,psubcomm);CHKERRQ(ierr);
255a297a907SKarl Rupp 
256d3b23db5SHong Zhang   /* set defaults */
257d3b23db5SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
258d3b23db5SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
259f68be91cSHong Zhang 
260d8a68f86SHong Zhang   (*psubcomm)->parent    = comm;
261d3b23db5SHong Zhang   (*psubcomm)->dupparent = comm;
262d3b23db5SHong Zhang   (*psubcomm)->comm      = PETSC_COMM_SELF;
263d3b23db5SHong Zhang   (*psubcomm)->n         = size;
264d3b23db5SHong Zhang   (*psubcomm)->color     = rank;
265e37c6257SHong Zhang   (*psubcomm)->subsize   = NULL;
266d3b23db5SHong Zhang   (*psubcomm)->type      = PETSC_SUBCOMM_INTERLACED;
267638faf0bSHong Zhang   PetscFunctionReturn(0);
268638faf0bSHong Zhang }
269638faf0bSHong Zhang 
270638faf0bSHong Zhang #undef __FUNCT__
27153c77d0aSJed Brown #define __FUNCT__ "PetscSubcommCreate_contiguous"
272d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm)
273638faf0bSHong Zhang {
274638faf0bSHong Zhang   PetscErrorCode ierr;
275d6037b41SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank=-1,subrank=-1;
276d8a68f86SHong Zhang   PetscInt       np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n;
277d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
278638faf0bSHong Zhang 
279638faf0bSHong Zhang   PetscFunctionBegin;
28055e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
28155e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
28255e3b8d2SHong Zhang 
283638faf0bSHong Zhang   /* get size of each subcommunicator */
284638faf0bSHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
285a297a907SKarl Rupp 
286638faf0bSHong Zhang   np_subcomm = size/nsubcomm;
287638faf0bSHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
288638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++) {
289638faf0bSHong Zhang     subsize[i] = np_subcomm;
290638faf0bSHong Zhang     if (i<nleftover) subsize[i]++;
291638faf0bSHong Zhang   }
292638faf0bSHong Zhang 
293638faf0bSHong Zhang   /* get color and subrank of this proc */
294638faf0bSHong Zhang   rankstart = 0;
295638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++) {
296638faf0bSHong Zhang     if (rank >= rankstart && rank < rankstart+subsize[i]) {
297638faf0bSHong Zhang       color   = i;
298638faf0bSHong Zhang       subrank = rank - rankstart;
299638faf0bSHong Zhang       duprank = rank;
300638faf0bSHong Zhang       break;
301a297a907SKarl Rupp     } else rankstart += subsize[i];
302638faf0bSHong Zhang   }
303638faf0bSHong Zhang 
304638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
305638faf0bSHong Zhang 
306638faf0bSHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
307638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
30852af3f7aSShri Abhyankar   {
30952af3f7aSShri Abhyankar     PetscThreadComm tcomm;
31052af3f7aSShri Abhyankar     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
31152af3f7aSShri Abhyankar     ierr = MPI_Attr_put(dupcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
31252af3f7aSShri Abhyankar     tcomm->refct++;
31352af3f7aSShri Abhyankar     ierr = MPI_Attr_put(subcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
31452af3f7aSShri Abhyankar     tcomm->refct++;
31552af3f7aSShri Abhyankar   }
3160298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
3170298fd71SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr);
318b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
319b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
320a297a907SKarl Rupp 
321d8a68f86SHong Zhang   psubcomm->color   = color;
322e37c6257SHong Zhang   psubcomm->subsize = subsize;
323f38d543fSHong Zhang   psubcomm->type    = PETSC_SUBCOMM_CONTIGUOUS;
324638faf0bSHong Zhang   PetscFunctionReturn(0);
325638faf0bSHong Zhang }
326638faf0bSHong Zhang 
327638faf0bSHong Zhang #undef __FUNCT__
328638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced"
329638faf0bSHong Zhang /*
330638faf0bSHong Zhang    Note:
331638faf0bSHong Zhang    In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators
33245fc02eaSBarry Smith    by iteratively taking a process into a subcommunicator.
333cd05a4c0SHong Zhang    Example: size=4, nsubcomm=(*psubcomm)->n=3
334cd05a4c0SHong Zhang      comm=(*psubcomm)->parent:
335cd05a4c0SHong Zhang       rank:     [0]  [1]  [2]  [3]
336cd05a4c0SHong Zhang       color:     0    1    2    0
337cd05a4c0SHong Zhang 
338cd05a4c0SHong Zhang      subcomm=(*psubcomm)->comm:
339cd05a4c0SHong Zhang       subrank:  [0]  [0]  [0]  [1]
340cd05a4c0SHong Zhang 
341cd05a4c0SHong Zhang      dupcomm=(*psubcomm)->dupparent:
342cd05a4c0SHong Zhang       duprank:  [0]  [2]  [3]  [1]
343cd05a4c0SHong Zhang 
344cd05a4c0SHong Zhang      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
345cd05a4c0SHong Zhang            subcomm[color = 1] has subsize=1, owns process [1]
346cd05a4c0SHong Zhang            subcomm[color = 2] has subsize=1, owns process [2]
347cd05a4c0SHong Zhang            dupcomm has same number of processes as comm, and its duprank maps
348cd05a4c0SHong Zhang            processes in subcomm contiguously into a 1d array:
349cd05a4c0SHong Zhang             duprank: [0] [1]      [2]         [3]
350cd05a4c0SHong Zhang             rank:    [0] [3]      [1]         [2]
351cd05a4c0SHong Zhang                     subcomm[0] subcomm[1]  subcomm[2]
352638faf0bSHong Zhang */
353cd05a4c0SHong Zhang 
354d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm)
355cd05a4c0SHong Zhang {
356cd05a4c0SHong Zhang   PetscErrorCode ierr;
357cd05a4c0SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank,subrank;
358d8a68f86SHong Zhang   PetscInt       np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n;
359d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
360cd05a4c0SHong Zhang 
361cd05a4c0SHong Zhang   PetscFunctionBegin;
36255e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
36355e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
36455e3b8d2SHong Zhang 
365cd05a4c0SHong Zhang   /* get size of each subcommunicator */
366cd05a4c0SHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
367a297a907SKarl Rupp 
368cd05a4c0SHong Zhang   np_subcomm = size/nsubcomm;
369cd05a4c0SHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
370cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++) {
371cd05a4c0SHong Zhang     subsize[i] = np_subcomm;
372cd05a4c0SHong Zhang     if (i<nleftover) subsize[i]++;
373cd05a4c0SHong Zhang   }
374cd05a4c0SHong Zhang 
375cd05a4c0SHong Zhang   /* find color for this proc */
376cd05a4c0SHong Zhang   color   = rank%nsubcomm;
377cd05a4c0SHong Zhang   subrank = rank/nsubcomm;
378cd05a4c0SHong Zhang 
379cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
380cd05a4c0SHong Zhang 
381cd05a4c0SHong Zhang   j = 0; duprank = 0;
382cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++) {
383cd05a4c0SHong Zhang     if (j == color) {
384cd05a4c0SHong Zhang       duprank += subrank;
385cd05a4c0SHong Zhang       break;
386cd05a4c0SHong Zhang     }
387cd05a4c0SHong Zhang     duprank += subsize[i]; j++;
388cd05a4c0SHong Zhang   }
389cd05a4c0SHong Zhang 
390cd05a4c0SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
391cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
39252af3f7aSShri Abhyankar   {
39352af3f7aSShri Abhyankar     PetscThreadComm tcomm;
39452af3f7aSShri Abhyankar     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
39552af3f7aSShri Abhyankar     ierr = MPI_Attr_put(dupcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
39652af3f7aSShri Abhyankar     tcomm->refct++;
39752af3f7aSShri Abhyankar     ierr = MPI_Attr_put(subcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
39852af3f7aSShri Abhyankar     tcomm->refct++;
39952af3f7aSShri Abhyankar   }
4000298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
4010298fd71SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr);
402b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
403b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
404a297a907SKarl Rupp 
405d8a68f86SHong Zhang   psubcomm->color   = color;
406e37c6257SHong Zhang   psubcomm->subsize = subsize;
407f38d543fSHong Zhang   psubcomm->type    = PETSC_SUBCOMM_INTERLACED;
408cd05a4c0SHong Zhang   PetscFunctionReturn(0);
409cd05a4c0SHong Zhang }
410638faf0bSHong Zhang 
411638faf0bSHong Zhang 
412