xref: /petsc/src/sys/objects/subcomm.c (revision 6a6fc655e6b1b6c3ff84be1691b1c5064ef5f4c2)
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 
8*6a6fc655SJed Brown const char *const PetscSubcommTypes[] = {"GENERAL","CONTIGUOUS","INTERLACED","PetscSubcommType","PETSC_SUBCOMM_",0};
9*6a6fc655SJed 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   }
117aa9c1079SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr);
118aa9c1079SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr);
119b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
120b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
1211ba920a7SHong Zhang   psubcomm->color     = color;
122d8a68f86SHong Zhang   PetscFunctionReturn(0);
123d8a68f86SHong Zhang }
124638faf0bSHong Zhang 
125cd05a4c0SHong Zhang #undef __FUNCT__
126cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy"
1276bf464f9SBarry Smith PetscErrorCode  PetscSubcommDestroy(PetscSubcomm *psubcomm)
128cd05a4c0SHong Zhang {
129cd05a4c0SHong Zhang   PetscErrorCode ierr;
130cd05a4c0SHong Zhang 
131cd05a4c0SHong Zhang   PetscFunctionBegin;
1326bf464f9SBarry Smith   if (!*psubcomm) PetscFunctionReturn(0);
133aa9c1079SBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr);
134aa9c1079SBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->comm);CHKERRQ(ierr);
1356bf464f9SBarry Smith   ierr = PetscFree((*psubcomm));CHKERRQ(ierr);
136cd05a4c0SHong Zhang   PetscFunctionReturn(0);
137cd05a4c0SHong Zhang }
138cd05a4c0SHong Zhang 
139cd05a4c0SHong Zhang #undef __FUNCT__
140cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate"
141ab8c242fSMatthew Knepley /*@C
142cd05a4c0SHong Zhang   PetscSubcommCreate - Create a PetscSubcomm context.
143cd05a4c0SHong Zhang 
144cd05a4c0SHong Zhang    Collective on MPI_Comm
145cd05a4c0SHong Zhang 
146cd05a4c0SHong Zhang    Input Parameter:
1479873d53eSJed Brown .  comm - MPI communicator
148cd05a4c0SHong Zhang 
149cd05a4c0SHong Zhang    Output Parameter:
150cd05a4c0SHong Zhang .  psubcomm - location to store the PetscSubcomm context
151cd05a4c0SHong Zhang 
152638faf0bSHong Zhang    Level: advanced
153cd05a4c0SHong Zhang 
154638faf0bSHong Zhang .keywords: communicator, create
155638faf0bSHong Zhang 
156638faf0bSHong Zhang .seealso: PetscSubcommDestroy()
157638faf0bSHong Zhang @*/
1587087cfbeSBarry Smith PetscErrorCode  PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm)
159638faf0bSHong Zhang {
160638faf0bSHong Zhang   PetscErrorCode ierr;
161638faf0bSHong Zhang 
162638faf0bSHong Zhang   PetscFunctionBegin;
163d8a68f86SHong Zhang   ierr = PetscNew(struct _n_PetscSubcomm,psubcomm);CHKERRQ(ierr);
164d8a68f86SHong Zhang   (*psubcomm)->parent = comm;
165638faf0bSHong Zhang   PetscFunctionReturn(0);
166638faf0bSHong Zhang }
167638faf0bSHong Zhang 
168638faf0bSHong Zhang #undef __FUNCT__
16953c77d0aSJed Brown #define __FUNCT__ "PetscSubcommCreate_contiguous"
170d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm)
171638faf0bSHong Zhang {
172638faf0bSHong Zhang   PetscErrorCode ierr;
173d6037b41SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank=-1,subrank=-1;
174d8a68f86SHong Zhang   PetscInt       np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n;
175d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
176638faf0bSHong Zhang 
177638faf0bSHong Zhang   PetscFunctionBegin;
17855e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
17955e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
18055e3b8d2SHong Zhang 
181638faf0bSHong Zhang   /* get size of each subcommunicator */
182638faf0bSHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
183638faf0bSHong Zhang   np_subcomm = size/nsubcomm;
184638faf0bSHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
185638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++){
186638faf0bSHong Zhang     subsize[i] = np_subcomm;
187638faf0bSHong Zhang     if (i<nleftover) subsize[i]++;
188638faf0bSHong Zhang   }
189638faf0bSHong Zhang 
190638faf0bSHong Zhang   /* get color and subrank of this proc */
191638faf0bSHong Zhang   rankstart = 0;
192638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++){
193638faf0bSHong Zhang     if ( rank >= rankstart && rank < rankstart+subsize[i]) {
194638faf0bSHong Zhang       color   = i;
195638faf0bSHong Zhang       subrank = rank - rankstart;
196638faf0bSHong Zhang       duprank = rank;
197638faf0bSHong Zhang       break;
198638faf0bSHong Zhang     } else {
199638faf0bSHong Zhang       rankstart += subsize[i];
200638faf0bSHong Zhang     }
201638faf0bSHong Zhang   }
202638faf0bSHong Zhang   ierr = PetscFree(subsize);CHKERRQ(ierr);
203638faf0bSHong Zhang 
204638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
205638faf0bSHong Zhang 
206638faf0bSHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
207638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
208638faf0bSHong Zhang 
209aa9c1079SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr);
210aa9c1079SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr);
211b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
212b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
213d8a68f86SHong Zhang   psubcomm->color     = color;
214422737f6SJed Brown 
215422737f6SJed Brown #if defined(PETSC_THREADCOMM_ACTIVE)
216422737f6SJed Brown   {
217422737f6SJed Brown     PetscThreadComm tcomm;
218422737f6SJed Brown     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
219422737f6SJed Brown     ierr = MPI_Attr_put(psubcomm->dupparent,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
220422737f6SJed Brown     tcomm->refct++;
221422737f6SJed Brown     ierr = MPI_Attr_put(psubcomm->comm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
222422737f6SJed Brown     tcomm->refct++;
223422737f6SJed Brown   }
224422737f6SJed Brown #endif
225638faf0bSHong Zhang   PetscFunctionReturn(0);
226638faf0bSHong Zhang }
227638faf0bSHong Zhang 
228638faf0bSHong Zhang #undef __FUNCT__
229638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced"
230638faf0bSHong Zhang /*
231638faf0bSHong Zhang    Note:
232638faf0bSHong Zhang    In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators
23345fc02eaSBarry Smith    by iteratively taking a process into a subcommunicator.
234cd05a4c0SHong Zhang    Example: size=4, nsubcomm=(*psubcomm)->n=3
235cd05a4c0SHong Zhang      comm=(*psubcomm)->parent:
236cd05a4c0SHong Zhang       rank:     [0]  [1]  [2]  [3]
237cd05a4c0SHong Zhang       color:     0    1    2    0
238cd05a4c0SHong Zhang 
239cd05a4c0SHong Zhang      subcomm=(*psubcomm)->comm:
240cd05a4c0SHong Zhang       subrank:  [0]  [0]  [0]  [1]
241cd05a4c0SHong Zhang 
242cd05a4c0SHong Zhang      dupcomm=(*psubcomm)->dupparent:
243cd05a4c0SHong Zhang       duprank:  [0]  [2]  [3]  [1]
244cd05a4c0SHong Zhang 
245cd05a4c0SHong Zhang      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
246cd05a4c0SHong Zhang            subcomm[color = 1] has subsize=1, owns process [1]
247cd05a4c0SHong Zhang            subcomm[color = 2] has subsize=1, owns process [2]
248cd05a4c0SHong Zhang            dupcomm has same number of processes as comm, and its duprank maps
249cd05a4c0SHong Zhang            processes in subcomm contiguously into a 1d array:
250cd05a4c0SHong Zhang             duprank: [0] [1]      [2]         [3]
251cd05a4c0SHong Zhang             rank:    [0] [3]      [1]         [2]
252cd05a4c0SHong Zhang                     subcomm[0] subcomm[1]  subcomm[2]
253638faf0bSHong Zhang */
254cd05a4c0SHong Zhang 
255d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm)
256cd05a4c0SHong Zhang {
257cd05a4c0SHong Zhang   PetscErrorCode ierr;
258cd05a4c0SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank,subrank;
259d8a68f86SHong Zhang   PetscInt       np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n;
260d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
261cd05a4c0SHong Zhang 
262cd05a4c0SHong Zhang   PetscFunctionBegin;
26355e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
26455e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
26555e3b8d2SHong Zhang 
266cd05a4c0SHong Zhang   /* get size of each subcommunicator */
267cd05a4c0SHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
268cd05a4c0SHong Zhang   np_subcomm = size/nsubcomm;
269cd05a4c0SHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
270cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++){
271cd05a4c0SHong Zhang     subsize[i] = np_subcomm;
272cd05a4c0SHong Zhang     if (i<nleftover) subsize[i]++;
273cd05a4c0SHong Zhang   }
274cd05a4c0SHong Zhang 
275cd05a4c0SHong Zhang   /* find color for this proc */
276cd05a4c0SHong Zhang   color   = rank%nsubcomm;
277cd05a4c0SHong Zhang   subrank = rank/nsubcomm;
278cd05a4c0SHong Zhang 
279cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
280cd05a4c0SHong Zhang 
281cd05a4c0SHong Zhang   j = 0; duprank = 0;
282cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++){
283cd05a4c0SHong Zhang     if (j == color){
284cd05a4c0SHong Zhang       duprank += subrank;
285cd05a4c0SHong Zhang       break;
286cd05a4c0SHong Zhang     }
287cd05a4c0SHong Zhang     duprank += subsize[i]; j++;
288cd05a4c0SHong Zhang   }
28945fc02eaSBarry Smith   ierr = PetscFree(subsize);CHKERRQ(ierr);
290cd05a4c0SHong Zhang 
291cd05a4c0SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
292cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
293cd05a4c0SHong Zhang 
294aa9c1079SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr);
295aa9c1079SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr);
296b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
297b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
298d8a68f86SHong Zhang   psubcomm->color     = color;
299422737f6SJed Brown 
300422737f6SJed Brown #if defined(PETSC_THREADCOMM_ACTIVE)
301422737f6SJed Brown   {
302422737f6SJed Brown     PetscThreadComm tcomm;
303422737f6SJed Brown     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
304422737f6SJed Brown     ierr = MPI_Attr_put(psubcomm->dupparent,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
305422737f6SJed Brown     tcomm->refct++;
306422737f6SJed Brown     ierr = MPI_Attr_put(psubcomm->comm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
307422737f6SJed Brown     tcomm->refct++;
308422737f6SJed Brown   }
309422737f6SJed Brown #endif
310cd05a4c0SHong Zhang   PetscFunctionReturn(0);
311cd05a4c0SHong Zhang }
312638faf0bSHong Zhang 
313638faf0bSHong Zhang 
314