xref: /petsc/src/sys/objects/subcomm.c (revision 574034a94232d912b0dd9c16e5221f02823d3038)
1 
2 /*
3      Provides utility routines for split MPI communicator.
4 */
5 #include <petscsys.h>    /*I   "petscsys.h"    I*/
6 #include <petsc-private/threadcommimpl.h> /* Petsc_ThreadComm_keyval */
7 
8 const char *const PetscSubcommTypes[] = {"GENERAL","CONTIGUOUS","INTERLACED","PetscSubcommType","PETSC_SUBCOMM_",0};
9 
10 extern PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm);
11 extern PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm);
12 
13 #undef __FUNCT__
14 #define __FUNCT__ "PetscSubcommSetNumber"
15 /*@C
16   PetscSubcommSetNumber - Set total number of subcommunicators.
17 
18    Collective on MPI_Comm
19 
20    Input Parameter:
21 +  psubcomm - PetscSubcomm context
22 -  nsubcomm - the total number of subcommunicators in psubcomm
23 
24    Level: advanced
25 
26 .keywords: communicator
27 
28 .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetType(),PetscSubcommSetTypeGeneral()
29 @*/
30 PetscErrorCode  PetscSubcommSetNumber(PetscSubcomm psubcomm,PetscInt nsubcomm)
31 {
32   PetscErrorCode ierr;
33   MPI_Comm       comm=psubcomm->parent;
34   PetscMPIInt    rank,size;
35 
36   PetscFunctionBegin;
37   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate() first");
38   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
39   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
40   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);
41 
42   psubcomm->n = nsubcomm;
43   PetscFunctionReturn(0);
44 }
45 
46 #undef __FUNCT__
47 #define __FUNCT__ "PetscSubcommSetType"
48 /*@C
49   PetscSubcommSetType - Set type of subcommunicators.
50 
51    Collective on MPI_Comm
52 
53    Input Parameter:
54 +  psubcomm - PetscSubcomm context
55 -  subcommtype - subcommunicator type, PETSC_SUBCOMM_CONTIGUOUS,PETSC_SUBCOMM_INTERLACED
56 
57    Level: advanced
58 
59 .keywords: communicator
60 
61 .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral()
62 @*/
63 PetscErrorCode  PetscSubcommSetType(PetscSubcomm psubcomm,PetscSubcommType subcommtype)
64 {
65   PetscErrorCode ierr;
66 
67   PetscFunctionBegin;
68   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
69   if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n);
70 
71   if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS){
72     ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr);
73   } else if (subcommtype == PETSC_SUBCOMM_INTERLACED){
74     ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr);
75   } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",subcommtype);
76   PetscFunctionReturn(0);
77 }
78 
79 #undef __FUNCT__
80 #define __FUNCT__ "PetscSubcommSetTypeGeneral"
81 /*@C
82   PetscSubcommSetTypeGeneral - Set type of subcommunicators from user's specifications
83 
84    Collective on MPI_Comm
85 
86    Input Parameter:
87 +  psubcomm - PetscSubcomm context
88 .  color   - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator.
89 .  subrank - rank in the subcommunicator
90 -  duprank - rank in the dupparent (see PetscSubcomm)
91 
92    Level: advanced
93 
94 .keywords: communicator, create
95 
96 .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType()
97 @*/
98 PetscErrorCode  PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank,PetscMPIInt duprank)
99 {
100   PetscErrorCode ierr;
101   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
102   PetscMPIInt    size;
103 
104   PetscFunctionBegin;
105   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
106   if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %D is incorrect. Call PetscSubcommSetNumber()",psubcomm->n);
107 
108   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
109 
110   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm
111      if duprank is not a valid number, then dupcomm is not created - not all applications require dupcomm! */
112   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
113   if (duprank == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"duprank==PETSC_DECIDE is not supported yet");
114   else if (duprank >= 0 && duprank < size){
115     ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
116   }
117   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr);
118   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr);
119   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
120   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
121   psubcomm->color     = color;
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "PetscSubcommDestroy"
127 PetscErrorCode  PetscSubcommDestroy(PetscSubcomm *psubcomm)
128 {
129   PetscErrorCode ierr;
130 
131   PetscFunctionBegin;
132   if (!*psubcomm) PetscFunctionReturn(0);
133   ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr);
134   ierr = PetscCommDestroy(&(*psubcomm)->comm);CHKERRQ(ierr);
135   ierr = PetscFree((*psubcomm));CHKERRQ(ierr);
136   PetscFunctionReturn(0);
137 }
138 
139 #undef __FUNCT__
140 #define __FUNCT__ "PetscSubcommCreate"
141 /*@C
142   PetscSubcommCreate - Create a PetscSubcomm context.
143 
144    Collective on MPI_Comm
145 
146    Input Parameter:
147 .  comm - MPI communicator
148 
149    Output Parameter:
150 .  psubcomm - location to store the PetscSubcomm context
151 
152    Level: advanced
153 
154 .keywords: communicator, create
155 
156 .seealso: PetscSubcommDestroy()
157 @*/
158 PetscErrorCode  PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm)
159 {
160   PetscErrorCode ierr;
161 
162   PetscFunctionBegin;
163   ierr = PetscNew(struct _n_PetscSubcomm,psubcomm);CHKERRQ(ierr);
164   (*psubcomm)->parent = comm;
165   PetscFunctionReturn(0);
166 }
167 
168 #undef __FUNCT__
169 #define __FUNCT__ "PetscSubcommCreate_contiguous"
170 PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm)
171 {
172   PetscErrorCode ierr;
173   PetscMPIInt    rank,size,*subsize,duprank=-1,subrank=-1;
174   PetscInt       np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n;
175   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
176 
177   PetscFunctionBegin;
178   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
179   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
180 
181   /* get size of each subcommunicator */
182   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
183   np_subcomm = size/nsubcomm;
184   nleftover  = size - nsubcomm*np_subcomm;
185   for (i=0; i<nsubcomm; i++){
186     subsize[i] = np_subcomm;
187     if (i<nleftover) subsize[i]++;
188   }
189 
190   /* get color and subrank of this proc */
191   rankstart = 0;
192   for (i=0; i<nsubcomm; i++){
193     if ( rank >= rankstart && rank < rankstart+subsize[i]) {
194       color   = i;
195       subrank = rank - rankstart;
196       duprank = rank;
197       break;
198     } else {
199       rankstart += subsize[i];
200     }
201   }
202   ierr = PetscFree(subsize);CHKERRQ(ierr);
203 
204   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
205 
206   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
207   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
208 
209   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr);
210   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr);
211   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
212   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
213   psubcomm->color     = color;
214 
215 #if defined(PETSC_THREADCOMM_ACTIVE)
216   {
217     PetscThreadComm tcomm;
218     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
219     ierr = MPI_Attr_put(psubcomm->dupparent,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
220     tcomm->refct++;
221     ierr = MPI_Attr_put(psubcomm->comm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
222     tcomm->refct++;
223   }
224 #endif
225   PetscFunctionReturn(0);
226 }
227 
228 #undef __FUNCT__
229 #define __FUNCT__ "PetscSubcommCreate_interlaced"
230 /*
231    Note:
232    In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators
233    by iteratively taking a process into a subcommunicator.
234    Example: size=4, nsubcomm=(*psubcomm)->n=3
235      comm=(*psubcomm)->parent:
236       rank:     [0]  [1]  [2]  [3]
237       color:     0    1    2    0
238 
239      subcomm=(*psubcomm)->comm:
240       subrank:  [0]  [0]  [0]  [1]
241 
242      dupcomm=(*psubcomm)->dupparent:
243       duprank:  [0]  [2]  [3]  [1]
244 
245      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
246            subcomm[color = 1] has subsize=1, owns process [1]
247            subcomm[color = 2] has subsize=1, owns process [2]
248            dupcomm has same number of processes as comm, and its duprank maps
249            processes in subcomm contiguously into a 1d array:
250             duprank: [0] [1]      [2]         [3]
251             rank:    [0] [3]      [1]         [2]
252                     subcomm[0] subcomm[1]  subcomm[2]
253 */
254 
255 PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm)
256 {
257   PetscErrorCode ierr;
258   PetscMPIInt    rank,size,*subsize,duprank,subrank;
259   PetscInt       np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n;
260   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
261 
262   PetscFunctionBegin;
263   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
264   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
265 
266   /* get size of each subcommunicator */
267   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
268   np_subcomm = size/nsubcomm;
269   nleftover  = size - nsubcomm*np_subcomm;
270   for (i=0; i<nsubcomm; i++){
271     subsize[i] = np_subcomm;
272     if (i<nleftover) subsize[i]++;
273   }
274 
275   /* find color for this proc */
276   color   = rank%nsubcomm;
277   subrank = rank/nsubcomm;
278 
279   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
280 
281   j = 0; duprank = 0;
282   for (i=0; i<nsubcomm; i++){
283     if (j == color){
284       duprank += subrank;
285       break;
286     }
287     duprank += subsize[i]; j++;
288   }
289   ierr = PetscFree(subsize);CHKERRQ(ierr);
290 
291   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
292   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
293 
294   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr);
295   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr);
296   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
297   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
298   psubcomm->color     = color;
299 
300 #if defined(PETSC_THREADCOMM_ACTIVE)
301   {
302     PetscThreadComm tcomm;
303     ierr = PetscCommGetThreadComm(comm,&tcomm);CHKERRQ(ierr);
304     ierr = MPI_Attr_put(psubcomm->dupparent,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
305     tcomm->refct++;
306     ierr = MPI_Attr_put(psubcomm->comm,Petsc_ThreadComm_keyval,tcomm);CHKERRQ(ierr);
307     tcomm->refct++;
308   }
309 #endif
310   PetscFunctionReturn(0);
311 }
312 
313 
314