xref: /petsc/src/sys/objects/subcomm.c (revision 19171117fd94ac73558cd466be253dd691477438)
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   }
41*19171117SHong Zhang 
42*19171117SHong Zhang   ierr = PetscOptionsHasName(NULL, "-psubcomm_view", &flg);CHKERRQ(ierr);
43*19171117SHong Zhang   if (flg) {
44*19171117SHong Zhang     ierr = PetscSubcommView(psubcomm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
45*19171117SHong 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
1471ba920a7SHong Zhang   PetscSubcommSetTypeGeneral - Set type of subcommunicators 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.
1541ba920a7SHong Zhang .  subrank - rank in the subcommunicator
1551ba920a7SHong Zhang -  duprank - rank in the dupparent (see PetscSubcomm)
1561ba920a7SHong Zhang 
1571ba920a7SHong Zhang    Level: advanced
1581ba920a7SHong Zhang 
1591ba920a7SHong Zhang .keywords: communicator, create
1601ba920a7SHong Zhang 
1611ba920a7SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType()
1621ba920a7SHong Zhang @*/
1637087cfbeSBarry Smith PetscErrorCode  PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank,PetscMPIInt duprank)
164d8a68f86SHong Zhang {
1651ba920a7SHong Zhang   PetscErrorCode ierr;
1661ba920a7SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
1671ba920a7SHong Zhang   PetscMPIInt    size;
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 
1751ba920a7SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm
1761ba920a7SHong Zhang      if duprank is not a valid number, then dupcomm is not created - not all applications require dupcomm! */
1771ba920a7SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1787c764164SBarry Smith   if (duprank == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"duprank==PETSC_DECIDE is not supported yet");
1797c764164SBarry Smith   else if (duprank >= 0 && duprank < size) {
1801ba920a7SHong Zhang     ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
1811ba920a7SHong Zhang   }
1820298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
1830298fd71SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr);
184b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
185b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
186a297a907SKarl Rupp 
1871ba920a7SHong Zhang   psubcomm->color = color;
188d8a68f86SHong Zhang   PetscFunctionReturn(0);
189d8a68f86SHong Zhang }
190638faf0bSHong Zhang 
191cd05a4c0SHong Zhang #undef __FUNCT__
192cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy"
1936bf464f9SBarry Smith PetscErrorCode  PetscSubcommDestroy(PetscSubcomm *psubcomm)
194cd05a4c0SHong Zhang {
195cd05a4c0SHong Zhang   PetscErrorCode ierr;
196cd05a4c0SHong Zhang 
197cd05a4c0SHong Zhang   PetscFunctionBegin;
1986bf464f9SBarry Smith   if (!*psubcomm) PetscFunctionReturn(0);
199aa9c1079SBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr);
200aa9c1079SBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->comm);CHKERRQ(ierr);
201e37c6257SHong Zhang   ierr = PetscFree((*psubcomm)->subsize);CHKERRQ(ierr);
2026bf464f9SBarry Smith   ierr = PetscFree((*psubcomm));CHKERRQ(ierr);
203cd05a4c0SHong Zhang   PetscFunctionReturn(0);
204cd05a4c0SHong Zhang }
205cd05a4c0SHong Zhang 
206cd05a4c0SHong Zhang #undef __FUNCT__
207cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate"
208ab8c242fSMatthew Knepley /*@C
209cd05a4c0SHong Zhang   PetscSubcommCreate - Create a PetscSubcomm context.
210cd05a4c0SHong Zhang 
211cd05a4c0SHong Zhang    Collective on MPI_Comm
212cd05a4c0SHong Zhang 
213cd05a4c0SHong Zhang    Input Parameter:
2149873d53eSJed Brown .  comm - MPI communicator
215cd05a4c0SHong Zhang 
216cd05a4c0SHong Zhang    Output Parameter:
217cd05a4c0SHong Zhang .  psubcomm - location to store the PetscSubcomm context
218cd05a4c0SHong Zhang 
219638faf0bSHong Zhang    Level: advanced
220cd05a4c0SHong Zhang 
221638faf0bSHong Zhang .keywords: communicator, create
222638faf0bSHong Zhang 
223638faf0bSHong Zhang .seealso: PetscSubcommDestroy()
224638faf0bSHong Zhang @*/
2257087cfbeSBarry Smith PetscErrorCode  PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm)
226638faf0bSHong Zhang {
227638faf0bSHong Zhang   PetscErrorCode ierr;
228d3b23db5SHong Zhang   PetscMPIInt    rank,size;
229638faf0bSHong Zhang 
230638faf0bSHong Zhang   PetscFunctionBegin;
231d8a68f86SHong Zhang   ierr = PetscNew(struct _n_PetscSubcomm,psubcomm);CHKERRQ(ierr);
232a297a907SKarl Rupp 
233d3b23db5SHong Zhang   /* set defaults */
234d3b23db5SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
235d3b23db5SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
236f68be91cSHong Zhang 
237d8a68f86SHong Zhang   (*psubcomm)->parent    = comm;
238d3b23db5SHong Zhang   (*psubcomm)->dupparent = comm;
239d3b23db5SHong Zhang   (*psubcomm)->comm      = PETSC_COMM_SELF;
240d3b23db5SHong Zhang   (*psubcomm)->n         = size;
241d3b23db5SHong Zhang   (*psubcomm)->color     = rank;
242e37c6257SHong Zhang   (*psubcomm)->subsize   = NULL;
243d3b23db5SHong Zhang   (*psubcomm)->type      = PETSC_SUBCOMM_INTERLACED;
244638faf0bSHong Zhang   PetscFunctionReturn(0);
245638faf0bSHong Zhang }
246638faf0bSHong Zhang 
247638faf0bSHong Zhang #undef __FUNCT__
24853c77d0aSJed Brown #define __FUNCT__ "PetscSubcommCreate_contiguous"
249d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm)
250638faf0bSHong Zhang {
251638faf0bSHong Zhang   PetscErrorCode ierr;
252d6037b41SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank=-1,subrank=-1;
253d8a68f86SHong Zhang   PetscInt       np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n;
254d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
255638faf0bSHong Zhang 
256638faf0bSHong Zhang   PetscFunctionBegin;
25755e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
25855e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
25955e3b8d2SHong Zhang 
260638faf0bSHong Zhang   /* get size of each subcommunicator */
261638faf0bSHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
262a297a907SKarl Rupp 
263638faf0bSHong Zhang   np_subcomm = size/nsubcomm;
264638faf0bSHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
265638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++) {
266638faf0bSHong Zhang     subsize[i] = np_subcomm;
267638faf0bSHong Zhang     if (i<nleftover) subsize[i]++;
268638faf0bSHong Zhang   }
269638faf0bSHong Zhang 
270638faf0bSHong Zhang   /* get color and subrank of this proc */
271638faf0bSHong Zhang   rankstart = 0;
272638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++) {
273638faf0bSHong Zhang     if (rank >= rankstart && rank < rankstart+subsize[i]) {
274638faf0bSHong Zhang       color   = i;
275638faf0bSHong Zhang       subrank = rank - rankstart;
276638faf0bSHong Zhang       duprank = rank;
277638faf0bSHong Zhang       break;
278a297a907SKarl Rupp     } else rankstart += subsize[i];
279638faf0bSHong Zhang   }
280638faf0bSHong Zhang 
281638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
282638faf0bSHong Zhang 
283638faf0bSHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
284638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
28552af3f7aSShri Abhyankar   {
28652af3f7aSShri Abhyankar     PetscThreadComm tcomm;
28752af3f7aSShri Abhyankar     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
28852af3f7aSShri Abhyankar     ierr = MPI_Attr_put(dupcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
28952af3f7aSShri Abhyankar     tcomm->refct++;
29052af3f7aSShri Abhyankar     ierr = MPI_Attr_put(subcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
29152af3f7aSShri Abhyankar     tcomm->refct++;
29252af3f7aSShri Abhyankar   }
2930298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
2940298fd71SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr);
295b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
296b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
297a297a907SKarl Rupp 
298d8a68f86SHong Zhang   psubcomm->color   = color;
299e37c6257SHong Zhang   psubcomm->subsize = subsize;
300f38d543fSHong Zhang   psubcomm->type    = PETSC_SUBCOMM_CONTIGUOUS;
301638faf0bSHong Zhang   PetscFunctionReturn(0);
302638faf0bSHong Zhang }
303638faf0bSHong Zhang 
304638faf0bSHong Zhang #undef __FUNCT__
305638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced"
306638faf0bSHong Zhang /*
307638faf0bSHong Zhang    Note:
308638faf0bSHong Zhang    In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators
30945fc02eaSBarry Smith    by iteratively taking a process into a subcommunicator.
310cd05a4c0SHong Zhang    Example: size=4, nsubcomm=(*psubcomm)->n=3
311cd05a4c0SHong Zhang      comm=(*psubcomm)->parent:
312cd05a4c0SHong Zhang       rank:     [0]  [1]  [2]  [3]
313cd05a4c0SHong Zhang       color:     0    1    2    0
314cd05a4c0SHong Zhang 
315cd05a4c0SHong Zhang      subcomm=(*psubcomm)->comm:
316cd05a4c0SHong Zhang       subrank:  [0]  [0]  [0]  [1]
317cd05a4c0SHong Zhang 
318cd05a4c0SHong Zhang      dupcomm=(*psubcomm)->dupparent:
319cd05a4c0SHong Zhang       duprank:  [0]  [2]  [3]  [1]
320cd05a4c0SHong Zhang 
321cd05a4c0SHong Zhang      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
322cd05a4c0SHong Zhang            subcomm[color = 1] has subsize=1, owns process [1]
323cd05a4c0SHong Zhang            subcomm[color = 2] has subsize=1, owns process [2]
324cd05a4c0SHong Zhang            dupcomm has same number of processes as comm, and its duprank maps
325cd05a4c0SHong Zhang            processes in subcomm contiguously into a 1d array:
326cd05a4c0SHong Zhang             duprank: [0] [1]      [2]         [3]
327cd05a4c0SHong Zhang             rank:    [0] [3]      [1]         [2]
328cd05a4c0SHong Zhang                     subcomm[0] subcomm[1]  subcomm[2]
329638faf0bSHong Zhang */
330cd05a4c0SHong Zhang 
331d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm)
332cd05a4c0SHong Zhang {
333cd05a4c0SHong Zhang   PetscErrorCode ierr;
334cd05a4c0SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank,subrank;
335d8a68f86SHong Zhang   PetscInt       np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n;
336d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
337cd05a4c0SHong Zhang 
338cd05a4c0SHong Zhang   PetscFunctionBegin;
33955e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
34055e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
34155e3b8d2SHong Zhang 
342cd05a4c0SHong Zhang   /* get size of each subcommunicator */
343cd05a4c0SHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
344a297a907SKarl Rupp 
345cd05a4c0SHong Zhang   np_subcomm = size/nsubcomm;
346cd05a4c0SHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
347cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++) {
348cd05a4c0SHong Zhang     subsize[i] = np_subcomm;
349cd05a4c0SHong Zhang     if (i<nleftover) subsize[i]++;
350cd05a4c0SHong Zhang   }
351cd05a4c0SHong Zhang 
352cd05a4c0SHong Zhang   /* find color for this proc */
353cd05a4c0SHong Zhang   color   = rank%nsubcomm;
354cd05a4c0SHong Zhang   subrank = rank/nsubcomm;
355cd05a4c0SHong Zhang 
356cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
357cd05a4c0SHong Zhang 
358cd05a4c0SHong Zhang   j = 0; duprank = 0;
359cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++) {
360cd05a4c0SHong Zhang     if (j == color) {
361cd05a4c0SHong Zhang       duprank += subrank;
362cd05a4c0SHong Zhang       break;
363cd05a4c0SHong Zhang     }
364cd05a4c0SHong Zhang     duprank += subsize[i]; j++;
365cd05a4c0SHong Zhang   }
366cd05a4c0SHong Zhang 
367cd05a4c0SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
368cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
36952af3f7aSShri Abhyankar   {
37052af3f7aSShri Abhyankar     PetscThreadComm tcomm;
37152af3f7aSShri Abhyankar     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
37252af3f7aSShri Abhyankar     ierr = MPI_Attr_put(dupcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
37352af3f7aSShri Abhyankar     tcomm->refct++;
37452af3f7aSShri Abhyankar     ierr = MPI_Attr_put(subcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
37552af3f7aSShri Abhyankar     tcomm->refct++;
37652af3f7aSShri Abhyankar   }
3770298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
3780298fd71SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr);
379b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
380b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
381a297a907SKarl Rupp 
382d8a68f86SHong Zhang   psubcomm->color   = color;
383e37c6257SHong Zhang   psubcomm->subsize = subsize;
384f38d543fSHong Zhang   psubcomm->type    = PETSC_SUBCOMM_INTERLACED;
385cd05a4c0SHong Zhang   PetscFunctionReturn(0);
386cd05a4c0SHong Zhang }
387638faf0bSHong Zhang 
388638faf0bSHong Zhang 
389