xref: /petsc/src/sys/objects/subcomm.c (revision 422737f6cfaf9dc70ea2c3abb819c59e3a7121c7)
17d0a6c19SBarry Smith 
2cd05a4c0SHong Zhang /*
3cd05a4c0SHong Zhang      Provides utility routines for split MPI communicator.
4cd05a4c0SHong Zhang */
5c6db04a5SJed Brown #include <petscsys.h>    /*I   "petscsys.h"    I*/
6*422737f6SJed Brown #include <petsc-private/threadcommimpl.h> /* Petsc_ThreadComm_keyval */
7cd05a4c0SHong Zhang 
8d8a68f86SHong Zhang extern PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm);
9d8a68f86SHong Zhang extern PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm);
10d8a68f86SHong Zhang 
11d8a68f86SHong Zhang #undef __FUNCT__
12d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetNumber"
13d8a68f86SHong Zhang /*@C
14d8a68f86SHong Zhang   PetscSubcommSetNumber - Set total number of subcommunicators.
15d8a68f86SHong Zhang 
16d8a68f86SHong Zhang    Collective on MPI_Comm
17d8a68f86SHong Zhang 
18d8a68f86SHong Zhang    Input Parameter:
19d8a68f86SHong Zhang +  psubcomm - PetscSubcomm context
20d8a68f86SHong Zhang -  nsubcomm - the total number of subcommunicators in psubcomm
21d8a68f86SHong Zhang 
22d8a68f86SHong Zhang    Level: advanced
23d8a68f86SHong Zhang 
24d8a68f86SHong Zhang .keywords: communicator
25d8a68f86SHong Zhang 
26d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetType(),PetscSubcommSetTypeGeneral()
27d8a68f86SHong Zhang @*/
287087cfbeSBarry Smith PetscErrorCode  PetscSubcommSetNumber(PetscSubcomm psubcomm,PetscInt nsubcomm)
29d8a68f86SHong Zhang {
30d8a68f86SHong Zhang   PetscErrorCode ierr;
31d8a68f86SHong Zhang   MPI_Comm       comm=psubcomm->parent;
32d8a68f86SHong Zhang   PetscMPIInt    rank,size;
33d8a68f86SHong Zhang 
34d8a68f86SHong Zhang   PetscFunctionBegin;
35d8a68f86SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate() first");
36d8a68f86SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
37d8a68f86SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
38d8a68f86SHong 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);
39d8a68f86SHong Zhang 
40d8a68f86SHong Zhang   psubcomm->n = nsubcomm;
41d8a68f86SHong Zhang   PetscFunctionReturn(0);
42d8a68f86SHong Zhang }
43d8a68f86SHong Zhang 
44d8a68f86SHong Zhang #undef __FUNCT__
45d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetType"
46d8a68f86SHong Zhang /*@C
47d8a68f86SHong Zhang   PetscSubcommSetType - Set type of subcommunicators.
48d8a68f86SHong Zhang 
49d8a68f86SHong Zhang    Collective on MPI_Comm
50d8a68f86SHong Zhang 
51d8a68f86SHong Zhang    Input Parameter:
52d8a68f86SHong Zhang +  psubcomm - PetscSubcomm context
531ba920a7SHong Zhang -  subcommtype - subcommunicator type, PETSC_SUBCOMM_CONTIGUOUS,PETSC_SUBCOMM_INTERLACED
54d8a68f86SHong Zhang 
55d8a68f86SHong Zhang    Level: advanced
56d8a68f86SHong Zhang 
57d8a68f86SHong Zhang .keywords: communicator
58d8a68f86SHong Zhang 
59d8a68f86SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral()
60d8a68f86SHong Zhang @*/
617c764164SBarry Smith PetscErrorCode  PetscSubcommSetType(PetscSubcomm psubcomm,PetscSubcommType subcommtype)
62d8a68f86SHong Zhang {
63d8a68f86SHong Zhang   PetscErrorCode ierr;
64d8a68f86SHong Zhang 
65d8a68f86SHong Zhang   PetscFunctionBegin;
66d8a68f86SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
67d8a68f86SHong Zhang   if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n);
68d8a68f86SHong Zhang 
69d8a68f86SHong Zhang   if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS){
70d8a68f86SHong Zhang     ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr);
71d8a68f86SHong Zhang   } else if (subcommtype == PETSC_SUBCOMM_INTERLACED){
72d8a68f86SHong Zhang     ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr);
737c764164SBarry Smith   } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",subcommtype);
74d8a68f86SHong Zhang   PetscFunctionReturn(0);
75d8a68f86SHong Zhang }
76d8a68f86SHong Zhang 
77d8a68f86SHong Zhang #undef __FUNCT__
78d8a68f86SHong Zhang #define __FUNCT__ "PetscSubcommSetTypeGeneral"
791ba920a7SHong Zhang /*@C
801ba920a7SHong Zhang   PetscSubcommSetTypeGeneral - Set type of subcommunicators from user's specifications
811ba920a7SHong Zhang 
821ba920a7SHong Zhang    Collective on MPI_Comm
831ba920a7SHong Zhang 
841ba920a7SHong Zhang    Input Parameter:
851ba920a7SHong Zhang +  psubcomm - PetscSubcomm context
861ba920a7SHong Zhang .  color   - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator.
871ba920a7SHong Zhang .  subrank - rank in the subcommunicator
881ba920a7SHong Zhang -  duprank - rank in the dupparent (see PetscSubcomm)
891ba920a7SHong Zhang 
901ba920a7SHong Zhang    Level: advanced
911ba920a7SHong Zhang 
921ba920a7SHong Zhang .keywords: communicator, create
931ba920a7SHong Zhang 
941ba920a7SHong Zhang .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType()
951ba920a7SHong Zhang @*/
967087cfbeSBarry Smith PetscErrorCode  PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank,PetscMPIInt duprank)
97d8a68f86SHong Zhang {
981ba920a7SHong Zhang   PetscErrorCode ierr;
991ba920a7SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
1001ba920a7SHong Zhang   PetscMPIInt    size;
1011ba920a7SHong Zhang 
102d8a68f86SHong Zhang   PetscFunctionBegin;
1031ba920a7SHong Zhang   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
1041ba920a7SHong Zhang   if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n);
1051ba920a7SHong Zhang 
1061ba920a7SHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
1071ba920a7SHong Zhang 
1081ba920a7SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm
1091ba920a7SHong Zhang      if duprank is not a valid number, then dupcomm is not created - not all applications require dupcomm! */
1101ba920a7SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1117c764164SBarry Smith   if (duprank == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"duprank==PETSC_DECIDE is not supported yet");
1127c764164SBarry Smith   else if (duprank >= 0 && duprank < size){
1131ba920a7SHong Zhang     ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
1141ba920a7SHong Zhang   }
115aa9c1079SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr);
116aa9c1079SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr);
117b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
118b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
1191ba920a7SHong Zhang   psubcomm->color     = color;
120d8a68f86SHong Zhang   PetscFunctionReturn(0);
121d8a68f86SHong Zhang }
122638faf0bSHong Zhang 
123cd05a4c0SHong Zhang #undef __FUNCT__
124cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy"
1256bf464f9SBarry Smith PetscErrorCode  PetscSubcommDestroy(PetscSubcomm *psubcomm)
126cd05a4c0SHong Zhang {
127cd05a4c0SHong Zhang   PetscErrorCode ierr;
128cd05a4c0SHong Zhang 
129cd05a4c0SHong Zhang   PetscFunctionBegin;
1306bf464f9SBarry Smith   if (!*psubcomm) PetscFunctionReturn(0);
131aa9c1079SBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr);
132aa9c1079SBarry Smith   ierr = PetscCommDestroy(&(*psubcomm)->comm);CHKERRQ(ierr);
1336bf464f9SBarry Smith   ierr = PetscFree((*psubcomm));CHKERRQ(ierr);
134cd05a4c0SHong Zhang   PetscFunctionReturn(0);
135cd05a4c0SHong Zhang }
136cd05a4c0SHong Zhang 
137cd05a4c0SHong Zhang #undef __FUNCT__
138cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate"
139ab8c242fSMatthew Knepley /*@C
140cd05a4c0SHong Zhang   PetscSubcommCreate - Create a PetscSubcomm context.
141cd05a4c0SHong Zhang 
142cd05a4c0SHong Zhang    Collective on MPI_Comm
143cd05a4c0SHong Zhang 
144cd05a4c0SHong Zhang    Input Parameter:
1459873d53eSJed Brown .  comm - MPI communicator
146cd05a4c0SHong Zhang 
147cd05a4c0SHong Zhang    Output Parameter:
148cd05a4c0SHong Zhang .  psubcomm - location to store the PetscSubcomm context
149cd05a4c0SHong Zhang 
150638faf0bSHong Zhang    Level: advanced
151cd05a4c0SHong Zhang 
152638faf0bSHong Zhang .keywords: communicator, create
153638faf0bSHong Zhang 
154638faf0bSHong Zhang .seealso: PetscSubcommDestroy()
155638faf0bSHong Zhang @*/
1567087cfbeSBarry Smith PetscErrorCode  PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm)
157638faf0bSHong Zhang {
158638faf0bSHong Zhang   PetscErrorCode ierr;
159638faf0bSHong Zhang 
160638faf0bSHong Zhang   PetscFunctionBegin;
161d8a68f86SHong Zhang   ierr = PetscNew(struct _n_PetscSubcomm,psubcomm);CHKERRQ(ierr);
162d8a68f86SHong Zhang   (*psubcomm)->parent = comm;
163638faf0bSHong Zhang   PetscFunctionReturn(0);
164638faf0bSHong Zhang }
165638faf0bSHong Zhang 
166638faf0bSHong Zhang #undef __FUNCT__
16753c77d0aSJed Brown #define __FUNCT__ "PetscSubcommCreate_contiguous"
168d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm)
169638faf0bSHong Zhang {
170638faf0bSHong Zhang   PetscErrorCode ierr;
171d6037b41SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank=-1,subrank=-1;
172d8a68f86SHong Zhang   PetscInt       np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n;
173d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
174638faf0bSHong Zhang 
175638faf0bSHong Zhang   PetscFunctionBegin;
17655e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
17755e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
17855e3b8d2SHong Zhang 
179638faf0bSHong Zhang   /* get size of each subcommunicator */
180638faf0bSHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
181638faf0bSHong Zhang   np_subcomm = size/nsubcomm;
182638faf0bSHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
183638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++){
184638faf0bSHong Zhang     subsize[i] = np_subcomm;
185638faf0bSHong Zhang     if (i<nleftover) subsize[i]++;
186638faf0bSHong Zhang   }
187638faf0bSHong Zhang 
188638faf0bSHong Zhang   /* get color and subrank of this proc */
189638faf0bSHong Zhang   rankstart = 0;
190638faf0bSHong Zhang   for (i=0; i<nsubcomm; i++){
191638faf0bSHong Zhang     if ( rank >= rankstart && rank < rankstart+subsize[i]) {
192638faf0bSHong Zhang       color   = i;
193638faf0bSHong Zhang       subrank = rank - rankstart;
194638faf0bSHong Zhang       duprank = rank;
195638faf0bSHong Zhang       break;
196638faf0bSHong Zhang     } else {
197638faf0bSHong Zhang       rankstart += subsize[i];
198638faf0bSHong Zhang     }
199638faf0bSHong Zhang   }
200638faf0bSHong Zhang   ierr = PetscFree(subsize);CHKERRQ(ierr);
201638faf0bSHong Zhang 
202638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
203638faf0bSHong Zhang 
204638faf0bSHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
205638faf0bSHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
206638faf0bSHong Zhang 
207aa9c1079SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr);
208aa9c1079SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr);
209b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
210b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
211d8a68f86SHong Zhang   psubcomm->color     = color;
212*422737f6SJed Brown 
213*422737f6SJed Brown #if defined(PETSC_THREADCOMM_ACTIVE)
214*422737f6SJed Brown   {
215*422737f6SJed Brown     PetscThreadComm tcomm;
216*422737f6SJed Brown     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
217*422737f6SJed Brown     ierr = MPI_Attr_put(psubcomm->dupparent,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
218*422737f6SJed Brown     tcomm->refct++;
219*422737f6SJed Brown     ierr = MPI_Attr_put(psubcomm->comm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
220*422737f6SJed Brown     tcomm->refct++;
221*422737f6SJed Brown   }
222*422737f6SJed Brown #endif
223638faf0bSHong Zhang   PetscFunctionReturn(0);
224638faf0bSHong Zhang }
225638faf0bSHong Zhang 
226638faf0bSHong Zhang #undef __FUNCT__
227638faf0bSHong Zhang #define __FUNCT__ "PetscSubcommCreate_interlaced"
228638faf0bSHong Zhang /*
229638faf0bSHong Zhang    Note:
230638faf0bSHong Zhang    In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators
23145fc02eaSBarry Smith    by iteratively taking a process into a subcommunicator.
232cd05a4c0SHong Zhang    Example: size=4, nsubcomm=(*psubcomm)->n=3
233cd05a4c0SHong Zhang      comm=(*psubcomm)->parent:
234cd05a4c0SHong Zhang       rank:     [0]  [1]  [2]  [3]
235cd05a4c0SHong Zhang       color:     0    1    2    0
236cd05a4c0SHong Zhang 
237cd05a4c0SHong Zhang      subcomm=(*psubcomm)->comm:
238cd05a4c0SHong Zhang       subrank:  [0]  [0]  [0]  [1]
239cd05a4c0SHong Zhang 
240cd05a4c0SHong Zhang      dupcomm=(*psubcomm)->dupparent:
241cd05a4c0SHong Zhang       duprank:  [0]  [2]  [3]  [1]
242cd05a4c0SHong Zhang 
243cd05a4c0SHong Zhang      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
244cd05a4c0SHong Zhang            subcomm[color = 1] has subsize=1, owns process [1]
245cd05a4c0SHong Zhang            subcomm[color = 2] has subsize=1, owns process [2]
246cd05a4c0SHong Zhang            dupcomm has same number of processes as comm, and its duprank maps
247cd05a4c0SHong Zhang            processes in subcomm contiguously into a 1d array:
248cd05a4c0SHong Zhang             duprank: [0] [1]      [2]         [3]
249cd05a4c0SHong Zhang             rank:    [0] [3]      [1]         [2]
250cd05a4c0SHong Zhang                     subcomm[0] subcomm[1]  subcomm[2]
251638faf0bSHong Zhang */
252cd05a4c0SHong Zhang 
253d8a68f86SHong Zhang PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm)
254cd05a4c0SHong Zhang {
255cd05a4c0SHong Zhang   PetscErrorCode ierr;
256cd05a4c0SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank,subrank;
257d8a68f86SHong Zhang   PetscInt       np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n;
258d8a68f86SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
259cd05a4c0SHong Zhang 
260cd05a4c0SHong Zhang   PetscFunctionBegin;
26155e3b8d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
26255e3b8d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
26355e3b8d2SHong Zhang 
264cd05a4c0SHong Zhang   /* get size of each subcommunicator */
265cd05a4c0SHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
266cd05a4c0SHong Zhang   np_subcomm = size/nsubcomm;
267cd05a4c0SHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
268cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++){
269cd05a4c0SHong Zhang     subsize[i] = np_subcomm;
270cd05a4c0SHong Zhang     if (i<nleftover) subsize[i]++;
271cd05a4c0SHong Zhang   }
272cd05a4c0SHong Zhang 
273cd05a4c0SHong Zhang   /* find color for this proc */
274cd05a4c0SHong Zhang   color   = rank%nsubcomm;
275cd05a4c0SHong Zhang   subrank = rank/nsubcomm;
276cd05a4c0SHong Zhang 
277cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
278cd05a4c0SHong Zhang 
279cd05a4c0SHong Zhang   j = 0; duprank = 0;
280cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++){
281cd05a4c0SHong Zhang     if (j == color){
282cd05a4c0SHong Zhang       duprank += subrank;
283cd05a4c0SHong Zhang       break;
284cd05a4c0SHong Zhang     }
285cd05a4c0SHong Zhang     duprank += subsize[i]; j++;
286cd05a4c0SHong Zhang   }
28745fc02eaSBarry Smith   ierr = PetscFree(subsize);CHKERRQ(ierr);
288cd05a4c0SHong Zhang 
289cd05a4c0SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
290cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
291cd05a4c0SHong Zhang 
292aa9c1079SBarry Smith   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr);
293aa9c1079SBarry Smith   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr);
294b89831e5SBarry Smith   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
295b89831e5SBarry Smith   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
296d8a68f86SHong Zhang   psubcomm->color     = color;
297*422737f6SJed Brown 
298*422737f6SJed Brown #if defined(PETSC_THREADCOMM_ACTIVE)
299*422737f6SJed Brown   {
300*422737f6SJed Brown     PetscThreadComm tcomm;
301*422737f6SJed Brown     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
302*422737f6SJed Brown     ierr = MPI_Attr_put(psubcomm->dupparent,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
303*422737f6SJed Brown     tcomm->refct++;
304*422737f6SJed Brown     ierr = MPI_Attr_put(psubcomm->comm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
305*422737f6SJed Brown     tcomm->refct++;
306*422737f6SJed Brown   }
307*422737f6SJed Brown #endif
308cd05a4c0SHong Zhang   PetscFunctionReturn(0);
309cd05a4c0SHong Zhang }
310638faf0bSHong Zhang 
311638faf0bSHong Zhang 
312