xref: /petsc/src/sys/objects/subcomm.c (revision d3b23db5a7ee1f155504da2761e4a9fbf8c1ff2c)
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 */
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);
12d8a68f86SHong Zhang 
13d8a68f86SHong Zhang #undef __FUNCT__
14d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetNumber"
15d8a68f86SHong Zhang /*@C
16d8a68f86SHong Zhang   PetscSubcommSetNumber - Set total number of subcommunicators.
17d8a68f86SHong Zhang 
18d8a68f86SHong Zhang    Collective on MPI_Comm
19d8a68f86SHong Zhang 
20d8a68f86SHong Zhang    Input Parameter:
21d8a68f86SHong Zhang +  psubcomm - PetscSubcomm context
22d8a68f86SHong Zhang -  nsubcomm - the total number of subcommunicators in psubcomm
23d8a68f86SHong Zhang 
24d8a68f86SHong Zhang    Level: advanced
25d8a68f86SHong Zhang 
26d8a68f86SHong Zhang .keywords: communicator
27d8a68f86SHong Zhang 
28d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetType(),PetscSubcommSetTypeGeneral()
29d8a68f86SHong Zhang @*/
307087cfbeSBarry Smith PetscErrorCode  PetscSubcommSetNumber(PetscSubcomm psubcomm,PetscInt nsubcomm)
31d8a68f86SHong Zhang {
32d8a68f86SHong Zhang   PetscErrorCode ierr;
33d8a68f86SHong Zhang   MPI_Comm       comm=psubcomm->parent;
34d8a68f86SHong Zhang   PetscMPIInt    rank,size;
35d8a68f86SHong Zhang 
36d8a68f86SHong Zhang   PetscFunctionBegin;
37d8a68f86SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate() first");
38d8a68f86SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
39d8a68f86SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
40d8a68f86SHong 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);
41d8a68f86SHong Zhang 
42d8a68f86SHong Zhang   psubcomm->n = nsubcomm;
43d8a68f86SHong Zhang   PetscFunctionReturn(0);
44d8a68f86SHong Zhang }
45d8a68f86SHong Zhang 
46d8a68f86SHong Zhang #undef __FUNCT__
47d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetType"
48d8a68f86SHong Zhang /*@C
49d8a68f86SHong Zhang   PetscSubcommSetType - Set type of subcommunicators.
50d8a68f86SHong Zhang 
51d8a68f86SHong Zhang    Collective on MPI_Comm
52d8a68f86SHong Zhang 
53d8a68f86SHong Zhang    Input Parameter:
54d8a68f86SHong Zhang +  psubcomm - PetscSubcomm context
551ba920a7SHong Zhang -  subcommtype - subcommunicator type, PETSC_SUBCOMM_CONTIGUOUS,PETSC_SUBCOMM_INTERLACED
56d8a68f86SHong Zhang 
57d8a68f86SHong Zhang    Level: advanced
58d8a68f86SHong Zhang 
59d8a68f86SHong Zhang .keywords: communicator
60d8a68f86SHong Zhang 
61d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral()
62d8a68f86SHong Zhang @*/
637c764164SBarry Smith PetscErrorCode  PetscSubcommSetType(PetscSubcomm psubcomm,PetscSubcommType subcommtype)
64d8a68f86SHong Zhang {
65d8a68f86SHong Zhang   PetscErrorCode ierr;
66d8a68f86SHong Zhang 
67d8a68f86SHong Zhang   PetscFunctionBegin;
68d8a68f86SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
69d8a68f86SHong Zhang   if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n);
70d8a68f86SHong Zhang 
71d8a68f86SHong Zhang   if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS) {
72d8a68f86SHong Zhang     ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr);
73d8a68f86SHong Zhang   } else if (subcommtype == PETSC_SUBCOMM_INTERLACED) {
74d8a68f86SHong Zhang     ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr);
757c764164SBarry Smith   } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",subcommtype);
76d8a68f86SHong Zhang   PetscFunctionReturn(0);
77d8a68f86SHong Zhang }
78d8a68f86SHong Zhang 
79d8a68f86SHong Zhang #undef __FUNCT__
80d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetTypeGeneral"
811ba920a7SHong Zhang /*@C
821ba920a7SHong Zhang   PetscSubcommSetTypeGeneral - Set type of subcommunicators from user's specifications
831ba920a7SHong Zhang 
841ba920a7SHong Zhang    Collective on MPI_Comm
851ba920a7SHong Zhang 
861ba920a7SHong Zhang    Input Parameter:
871ba920a7SHong Zhang +  psubcomm - PetscSubcomm context
881ba920a7SHong Zhang .  color   - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator.
891ba920a7SHong Zhang .  subrank - rank in the subcommunicator
901ba920a7SHong Zhang -  duprank - rank in the dupparent (see PetscSubcomm)
911ba920a7SHong Zhang 
921ba920a7SHong Zhang    Level: advanced
931ba920a7SHong Zhang 
941ba920a7SHong Zhang .keywords: communicator, create
951ba920a7SHong Zhang 
961ba920a7SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType()
971ba920a7SHong Zhang @*/
987087cfbeSBarry Smith PetscErrorCode  PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank,PetscMPIInt duprank)
99d8a68f86SHong Zhang {
1001ba920a7SHong Zhang   PetscErrorCode ierr;
1011ba920a7SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
1021ba920a7SHong Zhang   PetscMPIInt    size;
1031ba920a7SHong Zhang 
104d8a68f86SHong Zhang   PetscFunctionBegin;
1051ba920a7SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
1061ba920a7SHong Zhang   if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n);
1071ba920a7SHong Zhang 
1081ba920a7SHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
1091ba920a7SHong Zhang 
1101ba920a7SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm
1111ba920a7SHong Zhang      if duprank is not a valid number, then dupcomm is not created - not all applications require dupcomm! */
1121ba920a7SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1137c764164SBarry Smith   if (duprank == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"duprank==PETSC_DECIDE is not supported yet");
1147c764164SBarry Smith   else if (duprank >= 0 && duprank < size) {
1151ba920a7SHong Zhang     ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
1161ba920a7SHong Zhang   }
1170298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
1180298fd71SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr);
119b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
120b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
121a297a907SKarl Rupp 
1221ba920a7SHong Zhang   psubcomm->color = color;
123d8a68f86SHong Zhang   PetscFunctionReturn(0);
124d8a68f86SHong Zhang }
125638faf0bSHong Zhang 
126cd05a4c0SHong Zhang #undef __FUNCT__
127cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy"
1286bf464f9SBarry Smith PetscErrorCode  PetscSubcommDestroy(PetscSubcomm *psubcomm)
129cd05a4c0SHong Zhang {
130cd05a4c0SHong Zhang   PetscErrorCode ierr;
131cd05a4c0SHong Zhang 
132cd05a4c0SHong Zhang   PetscFunctionBegin;
1336bf464f9SBarry Smith   if (!*psubcomm) PetscFunctionReturn(0);
134aa9c1079SBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr);
135aa9c1079SBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->comm);CHKERRQ(ierr);
1366bf464f9SBarry Smith   ierr = PetscFree((*psubcomm));CHKERRQ(ierr);
137cd05a4c0SHong Zhang   PetscFunctionReturn(0);
138cd05a4c0SHong Zhang }
139cd05a4c0SHong Zhang 
140cd05a4c0SHong Zhang #undef __FUNCT__
141cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate"
142ab8c242fSMatthew Knepley /*@C
143cd05a4c0SHong Zhang   PetscSubcommCreate - Create a PetscSubcomm context.
144cd05a4c0SHong Zhang 
145cd05a4c0SHong Zhang    Collective on MPI_Comm
146cd05a4c0SHong Zhang 
147cd05a4c0SHong Zhang    Input Parameter:
1489873d53eSJed Brown .  comm - MPI communicator
149cd05a4c0SHong Zhang 
150cd05a4c0SHong Zhang    Output Parameter:
151cd05a4c0SHong Zhang .  psubcomm - location to store the PetscSubcomm context
152cd05a4c0SHong Zhang 
153638faf0bSHong Zhang    Level: advanced
154cd05a4c0SHong Zhang 
155638faf0bSHong Zhang .keywords: communicator, create
156638faf0bSHong Zhang 
157638faf0bSHong Zhang .seealso: PetscSubcommDestroy()
158638faf0bSHong Zhang @*/
1597087cfbeSBarry Smith PetscErrorCode  PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm)
160638faf0bSHong Zhang {
161638faf0bSHong Zhang   PetscErrorCode ierr;
162*d3b23db5SHong Zhang   PetscMPIInt    rank,size;
163638faf0bSHong Zhang 
164638faf0bSHong Zhang   PetscFunctionBegin;
165d8a68f86SHong Zhang   ierr = PetscNew(struct _n_PetscSubcomm,psubcomm);CHKERRQ(ierr);
166a297a907SKarl Rupp 
167*d3b23db5SHong Zhang   /* set defaults */
168*d3b23db5SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
169*d3b23db5SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
170d8a68f86SHong Zhang   (*psubcomm)->parent    = comm;
171*d3b23db5SHong Zhang   (*psubcomm)->dupparent = comm;
172*d3b23db5SHong Zhang   (*psubcomm)->comm      = PETSC_COMM_SELF;
173*d3b23db5SHong Zhang   (*psubcomm)->n         = size;
174*d3b23db5SHong Zhang   (*psubcomm)->color     = rank;
175*d3b23db5SHong Zhang   (*psubcomm)->type      = PETSC_SUBCOMM_INTERLACED;
176638faf0bSHong Zhang   PetscFunctionReturn(0);
177638faf0bSHong Zhang }
178638faf0bSHong Zhang 
179638faf0bSHong Zhang #undef __FUNCT__
18053c77d0aSJed Brown #define __FUNCT__ "PetscSubcommCreate_contiguous"
181d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm)
182638faf0bSHong Zhang {
183638faf0bSHong Zhang   PetscErrorCode ierr;
184d6037b41SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank=-1,subrank=-1;
185d8a68f86SHong Zhang   PetscInt       np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n;
186d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
187638faf0bSHong Zhang 
188638faf0bSHong Zhang   PetscFunctionBegin;
18955e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
19055e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
19155e3b8d2SHong Zhang 
192638faf0bSHong Zhang   /* get size of each subcommunicator */
193638faf0bSHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
194a297a907SKarl Rupp 
195638faf0bSHong Zhang   np_subcomm = size/nsubcomm;
196638faf0bSHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
197638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++) {
198638faf0bSHong Zhang     subsize[i] = np_subcomm;
199638faf0bSHong Zhang     if (i<nleftover) subsize[i]++;
200638faf0bSHong Zhang   }
201638faf0bSHong Zhang 
202638faf0bSHong Zhang   /* get color and subrank of this proc */
203638faf0bSHong Zhang   rankstart = 0;
204638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++) {
205638faf0bSHong Zhang     if (rank >= rankstart && rank < rankstart+subsize[i]) {
206638faf0bSHong Zhang       color   = i;
207638faf0bSHong Zhang       subrank = rank - rankstart;
208638faf0bSHong Zhang       duprank = rank;
209638faf0bSHong Zhang       break;
210a297a907SKarl Rupp     } else rankstart += subsize[i];
211638faf0bSHong Zhang   }
212638faf0bSHong Zhang   ierr = PetscFree(subsize);CHKERRQ(ierr);
213638faf0bSHong Zhang 
214638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
215638faf0bSHong Zhang 
216638faf0bSHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
217638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
21852af3f7aSShri Abhyankar   {
21952af3f7aSShri Abhyankar     PetscThreadComm tcomm;
22052af3f7aSShri Abhyankar     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
22152af3f7aSShri Abhyankar     ierr = MPI_Attr_put(dupcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
22252af3f7aSShri Abhyankar     tcomm->refct++;
22352af3f7aSShri Abhyankar     ierr = MPI_Attr_put(subcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
22452af3f7aSShri Abhyankar     tcomm->refct++;
22552af3f7aSShri Abhyankar   }
2260298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
2270298fd71SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr);
228b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
229b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
230a297a907SKarl Rupp 
231d8a68f86SHong Zhang   psubcomm->color = color;
232f38d543fSHong Zhang   psubcomm->type  = PETSC_SUBCOMM_CONTIGUOUS;
233638faf0bSHong Zhang   PetscFunctionReturn(0);
234638faf0bSHong Zhang }
235638faf0bSHong Zhang 
236638faf0bSHong Zhang #undef __FUNCT__
237638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced"
238638faf0bSHong Zhang /*
239638faf0bSHong Zhang    Note:
240638faf0bSHong Zhang    In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators
24145fc02eaSBarry Smith    by iteratively taking a process into a subcommunicator.
242cd05a4c0SHong Zhang    Example: size=4, nsubcomm=(*psubcomm)->n=3
243cd05a4c0SHong Zhang      comm=(*psubcomm)->parent:
244cd05a4c0SHong Zhang       rank:     [0]  [1]  [2]  [3]
245cd05a4c0SHong Zhang       color:     0    1    2    0
246cd05a4c0SHong Zhang 
247cd05a4c0SHong Zhang      subcomm=(*psubcomm)->comm:
248cd05a4c0SHong Zhang       subrank:  [0]  [0]  [0]  [1]
249cd05a4c0SHong Zhang 
250cd05a4c0SHong Zhang      dupcomm=(*psubcomm)->dupparent:
251cd05a4c0SHong Zhang       duprank:  [0]  [2]  [3]  [1]
252cd05a4c0SHong Zhang 
253cd05a4c0SHong Zhang      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
254cd05a4c0SHong Zhang            subcomm[color = 1] has subsize=1, owns process [1]
255cd05a4c0SHong Zhang            subcomm[color = 2] has subsize=1, owns process [2]
256cd05a4c0SHong Zhang            dupcomm has same number of processes as comm, and its duprank maps
257cd05a4c0SHong Zhang            processes in subcomm contiguously into a 1d array:
258cd05a4c0SHong Zhang             duprank: [0] [1]      [2]         [3]
259cd05a4c0SHong Zhang             rank:    [0] [3]      [1]         [2]
260cd05a4c0SHong Zhang                     subcomm[0] subcomm[1]  subcomm[2]
261638faf0bSHong Zhang */
262cd05a4c0SHong Zhang 
263d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm)
264cd05a4c0SHong Zhang {
265cd05a4c0SHong Zhang   PetscErrorCode ierr;
266cd05a4c0SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank,subrank;
267d8a68f86SHong Zhang   PetscInt       np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n;
268d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
269cd05a4c0SHong Zhang 
270cd05a4c0SHong Zhang   PetscFunctionBegin;
27155e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
27255e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
27355e3b8d2SHong Zhang 
274cd05a4c0SHong Zhang   /* get size of each subcommunicator */
275cd05a4c0SHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
276a297a907SKarl Rupp 
277cd05a4c0SHong Zhang   np_subcomm = size/nsubcomm;
278cd05a4c0SHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
279cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++) {
280cd05a4c0SHong Zhang     subsize[i] = np_subcomm;
281cd05a4c0SHong Zhang     if (i<nleftover) subsize[i]++;
282cd05a4c0SHong Zhang   }
283cd05a4c0SHong Zhang 
284cd05a4c0SHong Zhang   /* find color for this proc */
285cd05a4c0SHong Zhang   color   = rank%nsubcomm;
286cd05a4c0SHong Zhang   subrank = rank/nsubcomm;
287cd05a4c0SHong Zhang 
288cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
289cd05a4c0SHong Zhang 
290cd05a4c0SHong Zhang   j = 0; duprank = 0;
291cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++) {
292cd05a4c0SHong Zhang     if (j == color) {
293cd05a4c0SHong Zhang       duprank += subrank;
294cd05a4c0SHong Zhang       break;
295cd05a4c0SHong Zhang     }
296cd05a4c0SHong Zhang     duprank += subsize[i]; j++;
297cd05a4c0SHong Zhang   }
29845fc02eaSBarry Smith   ierr = PetscFree(subsize);CHKERRQ(ierr);
299cd05a4c0SHong Zhang 
300cd05a4c0SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
301cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
30252af3f7aSShri Abhyankar   {
30352af3f7aSShri Abhyankar     PetscThreadComm tcomm;
30452af3f7aSShri Abhyankar     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
30552af3f7aSShri Abhyankar     ierr = MPI_Attr_put(dupcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
30652af3f7aSShri Abhyankar     tcomm->refct++;
30752af3f7aSShri Abhyankar     ierr = MPI_Attr_put(subcomm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
30852af3f7aSShri Abhyankar     tcomm->refct++;
30952af3f7aSShri Abhyankar   }
3100298fd71SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
3110298fd71SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,NULL);CHKERRQ(ierr);
312b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
313b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
314a297a907SKarl Rupp 
315d8a68f86SHong Zhang   psubcomm->color = color;
316f38d543fSHong Zhang   psubcomm->type  = PETSC_SUBCOMM_INTERLACED;
317cd05a4c0SHong Zhang   PetscFunctionReturn(0);
318cd05a4c0SHong Zhang }
319638faf0bSHong Zhang 
320638faf0bSHong Zhang 
321