xref: /petsc/src/sys/objects/subcomm.c (revision 45f3bb6e2ab689450ec389543af70463a1ef8274)
1 
2 /*
3      Provides utility routines for split MPI communicator.
4 */
5 #include <petscsys.h>    /*I   "petscsys.h"    I*/
6 #include <petscviewer.h>
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__ "PetscSubcommSetFromOptions"
15 /*@C
16    PetscSubcommSetFromOptions - Allows setting options from a PetscSubcomm
17 
18    Collective on PetscSubcomm
19 
20    Input Parameter:
21 .  psubcomm - PetscSubcomm context
22 
23    Level: beginner
24 
25 @*/
26 PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm psubcomm)
27 {
28   PetscErrorCode   ierr;
29   PetscSubcommType type;
30   PetscBool        flg;
31 
32   PetscFunctionBegin;
33   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Must call PetscSubcommCreate firt");
34 
35   ierr = PetscOptionsBegin(psubcomm->parent,psubcomm->subcommprefix,"Options for PetscSubcomm",NULL);CHKERRQ(ierr);
36   ierr = PetscOptionsEnum("-psubcomm_type",NULL,NULL,PetscSubcommTypes,(PetscEnum)psubcomm->type,(PetscEnum*)&type,&flg);CHKERRQ(ierr);
37   if (flg && psubcomm->type != type) {
38     /* free old structures */
39     ierr = PetscCommDestroy(&(psubcomm)->dupparent);CHKERRQ(ierr);
40     ierr = PetscCommDestroy(&(psubcomm)->child);CHKERRQ(ierr);
41     ierr = PetscFree((psubcomm)->subsize);CHKERRQ(ierr);
42     switch (type) {
43     case PETSC_SUBCOMM_GENERAL:
44       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Runtime option PETSC_SUBCOMM_GENERAL is not supported, use PetscSubcommSetTypeGeneral()");
45     case PETSC_SUBCOMM_CONTIGUOUS:
46       ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr);
47       break;
48     case PETSC_SUBCOMM_INTERLACED:
49       ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr);
50       break;
51     default:
52       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSubcommType %s is not supported yet",PetscSubcommTypes[type]);
53     }
54   }
55 
56   ierr = PetscOptionsName("-psubcomm_view","Triggers display of PetscSubcomm context","PetscSubcommView",&flg);CHKERRQ(ierr);
57   if (flg) {
58     ierr = PetscSubcommView(psubcomm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
59   }
60   ierr = PetscOptionsEnd();CHKERRQ(ierr);
61   PetscFunctionReturn(0);
62 }
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "PetscSubcommSetOptionsPrefix"
66 /*@C
67   PetscSubcommSetOptionsPrefix - Sets the prefix used for searching for all
68   PetscSubcomm items in the options database.
69 
70   Logically collective on PetscSubcomm.
71 
72   Level: Intermediate
73 
74   Input Parameters:
75 +   psubcomm - PetscSubcomm context
76 -   prefix - the prefix to prepend all PetscSubcomm item names with.
77 
78 @*/
79 PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm psubcomm,const char pre[])
80 {
81   PetscErrorCode   ierr;
82 
83   PetscFunctionBegin;
84    if (!pre) {
85     ierr = PetscFree(psubcomm->subcommprefix);CHKERRQ(ierr);
86   } else {
87     if (pre[0] == '-') SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Options prefix should not begin with a hypen");
88     ierr = PetscFree(psubcomm->subcommprefix);CHKERRQ(ierr);
89     ierr = PetscStrallocpy(pre,&(psubcomm->subcommprefix));CHKERRQ(ierr);
90   }
91   PetscFunctionReturn(0);
92 }
93 
94 #undef __FUNCT__
95 #define __FUNCT__ "PetscSubcommView"
96 /*@C
97    PetsSubcommcView - Views a PetscSubcomm of values as either ASCII text or a binary file
98 
99    Collective on PetscSubcomm
100 
101    Input Parameter:
102 +  psubcomm - PetscSubcomm context
103 -  viewer - location to view the values
104 
105    Level: beginner
106 @*/
107 PetscErrorCode PetscSubcommView(PetscSubcomm psubcomm,PetscViewer viewer)
108 {
109   PetscErrorCode    ierr;
110   PetscBool         iascii;
111   PetscViewerFormat format;
112 
113   PetscFunctionBegin;
114   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
115   if (iascii) {
116     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
117     if (format == PETSC_VIEWER_DEFAULT) {
118       MPI_Comm    comm=psubcomm->parent;
119       PetscMPIInt rank,size,subsize,subrank,duprank;
120 
121       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
122       ierr = PetscViewerASCIIPrintf(viewer,"PetscSubcomm type %s with total %d MPI processes:\n",PetscSubcommTypes[psubcomm->type],size);CHKERRQ(ierr);
123       ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
124       ierr = MPI_Comm_size(psubcomm->child,&subsize);CHKERRQ(ierr);
125       ierr = MPI_Comm_rank(psubcomm->child,&subrank);CHKERRQ(ierr);
126       ierr = MPI_Comm_rank(psubcomm->dupparent,&duprank);CHKERRQ(ierr);
127       ierr = PetscSynchronizedPrintf(comm,"  [%d], color %d, sub-size %d, sub-rank %d, duprank %d\n",rank,psubcomm->color,subsize,subrank,duprank);CHKERRQ(ierr);
128       ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
129     }
130   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported yet");
131   PetscFunctionReturn(0);
132 }
133 
134 #undef __FUNCT__
135 #define __FUNCT__ "PetscSubcommSetNumber"
136 /*@C
137   PetscSubcommSetNumber - Set total number of subcommunicators.
138 
139    Collective on MPI_Comm
140 
141    Input Parameter:
142 +  psubcomm - PetscSubcomm context
143 -  nsubcomm - the total number of subcommunicators in psubcomm
144 
145    Level: advanced
146 
147 .keywords: communicator
148 
149 .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetType(),PetscSubcommSetTypeGeneral()
150 @*/
151 PetscErrorCode  PetscSubcommSetNumber(PetscSubcomm psubcomm,PetscInt nsubcomm)
152 {
153   PetscErrorCode ierr;
154   MPI_Comm       comm=psubcomm->parent;
155   PetscMPIInt    msub,size;
156 
157   PetscFunctionBegin;
158   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate() first");
159   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
160   ierr = PetscMPIIntCast(nsubcomm,&msub);CHKERRQ(ierr);
161   if (msub < 1 || msub > size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Num of subcommunicators %d cannot be < 1 or > input comm size %d",msub,size);
162 
163   psubcomm->n = msub;
164   PetscFunctionReturn(0);
165 }
166 
167 #undef __FUNCT__
168 #define __FUNCT__ "PetscSubcommSetType"
169 /*@C
170   PetscSubcommSetType - Set type of subcommunicators.
171 
172    Collective on MPI_Comm
173 
174    Input Parameter:
175 +  psubcomm - PetscSubcomm context
176 -  subcommtype - subcommunicator type, PETSC_SUBCOMM_CONTIGUOUS,PETSC_SUBCOMM_INTERLACED
177 
178    Level: advanced
179 
180 .keywords: communicator
181 
182 .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetTypeGeneral()
183 @*/
184 PetscErrorCode  PetscSubcommSetType(PetscSubcomm psubcomm,PetscSubcommType subcommtype)
185 {
186   PetscErrorCode ierr;
187 
188   PetscFunctionBegin;
189   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
190   if (psubcomm->n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()",psubcomm->n);
191 
192   if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS) {
193     ierr = PetscSubcommCreate_contiguous(psubcomm);CHKERRQ(ierr);
194   } else if (subcommtype == PETSC_SUBCOMM_INTERLACED) {
195     ierr = PetscSubcommCreate_interlaced(psubcomm);CHKERRQ(ierr);
196   } else SETERRQ1(psubcomm->parent,PETSC_ERR_SUP,"PetscSubcommType %s is not supported yet",PetscSubcommTypes[subcommtype]);
197   PetscFunctionReturn(0);
198 }
199 
200 #undef __FUNCT__
201 #define __FUNCT__ "PetscSubcommSetTypeGeneral"
202 /*@C
203   PetscSubcommSetTypeGeneral - Set a PetscSubcomm from user's specifications
204 
205    Collective on MPI_Comm
206 
207    Input Parameter:
208 +  psubcomm - PetscSubcomm context
209 .  color   - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator.
210 -  subrank - rank in the subcommunicator
211 
212    Level: advanced
213 
214 .keywords: communicator, create
215 
216 .seealso: PetscSubcommCreate(),PetscSubcommDestroy(),PetscSubcommSetNumber(),PetscSubcommSetType()
217 @*/
218 PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm,PetscMPIInt color,PetscMPIInt subrank)
219 {
220   PetscErrorCode ierr;
221   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
222   PetscMPIInt    size,icolor,duprank,*recvbuf,sendbuf[3],mysubsize,rank,*subsize;
223   PetscMPIInt    i,nsubcomm=psubcomm->n;
224 
225   PetscFunctionBegin;
226   if (!psubcomm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"PetscSubcomm is not created. Call PetscSubcommCreate()");
227   if (nsubcomm < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()",nsubcomm);
228 
229   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
230 
231   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
232   /* TODO: this can be done in an ostensibly scalale way (i.e., without allocating an array of size 'size') as is done in PetscObjectsCreateGlobalOrdering(). */
233   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
234   ierr = PetscMalloc1(2*size,&recvbuf);CHKERRQ(ierr);
235 
236   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
237   ierr = MPI_Comm_size(subcomm,&mysubsize);CHKERRQ(ierr);
238 
239   sendbuf[0] = color;
240   sendbuf[1] = mysubsize;
241   ierr = MPI_Allgather(sendbuf,2,MPI_INT,recvbuf,2,MPI_INT,comm);CHKERRQ(ierr);
242 
243   ierr = PetscCalloc1(nsubcomm,&subsize);CHKERRQ(ierr);
244   for (i=0; i<2*size; i+=2) {
245     subsize[recvbuf[i]] = recvbuf[i+1];
246   }
247   ierr = PetscFree(recvbuf);CHKERRQ(ierr);
248 
249   duprank = 0;
250   for (icolor=0; icolor<nsubcomm; icolor++) {
251     if (icolor != color) { /* not color of this process */
252       duprank += subsize[icolor];
253     } else {
254       duprank += subrank;
255       break;
256     }
257   }
258   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
259 
260   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
261   ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr);
262   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
263   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
264 
265   psubcomm->color   = color;
266   psubcomm->subsize = subsize;
267   psubcomm->type    = PETSC_SUBCOMM_GENERAL;
268   PetscFunctionReturn(0);
269 }
270 
271 #undef __FUNCT__
272 #define __FUNCT__ "PetscSubcommDestroy"
273 PetscErrorCode  PetscSubcommDestroy(PetscSubcomm *psubcomm)
274 {
275   PetscErrorCode ierr;
276 
277   PetscFunctionBegin;
278   if (!*psubcomm) PetscFunctionReturn(0);
279   ierr = PetscCommDestroy(&(*psubcomm)->dupparent);CHKERRQ(ierr);
280   ierr = PetscCommDestroy(&(*psubcomm)->child);CHKERRQ(ierr);
281   ierr = PetscFree((*psubcomm)->subsize);CHKERRQ(ierr);
282   if ((*psubcomm)->subcommprefix) { ierr = PetscFree((*psubcomm)->subcommprefix);CHKERRQ(ierr); }
283   ierr = PetscFree((*psubcomm));CHKERRQ(ierr);
284   PetscFunctionReturn(0);
285 }
286 
287 #undef __FUNCT__
288 #define __FUNCT__ "PetscSubcommCreate"
289 /*@C
290   PetscSubcommCreate - Create a PetscSubcomm context.
291 
292    Collective on MPI_Comm
293 
294    Input Parameter:
295 .  comm - MPI communicator
296 
297    Output Parameter:
298 .  psubcomm - location to store the PetscSubcomm context
299 
300    Level: advanced
301 
302 .keywords: communicator, create
303 
304 .seealso: PetscSubcommDestroy()
305 @*/
306 PetscErrorCode  PetscSubcommCreate(MPI_Comm comm,PetscSubcomm *psubcomm)
307 {
308   PetscErrorCode ierr;
309   PetscMPIInt    rank,size;
310 
311   PetscFunctionBegin;
312   ierr = PetscNew(psubcomm);CHKERRQ(ierr);
313 
314   /* set defaults */
315   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
316   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
317 
318   (*psubcomm)->parent    = comm;
319   (*psubcomm)->dupparent = comm;
320   (*psubcomm)->child     = PETSC_COMM_SELF;
321   (*psubcomm)->n         = size;
322   (*psubcomm)->color     = rank;
323   (*psubcomm)->subsize   = NULL;
324   (*psubcomm)->type      = PETSC_SUBCOMM_INTERLACED;
325   PetscFunctionReturn(0);
326 }
327 
328 #undef __FUNCT__
329 #define __FUNCT__ "PetscSubcommCreate_contiguous"
330 PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm)
331 {
332   PetscErrorCode ierr;
333   PetscMPIInt    rank,size,*subsize,duprank=-1,subrank=-1;
334   PetscMPIInt    np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n;
335   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
336 
337   PetscFunctionBegin;
338   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
339   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
340 
341   /* get size of each subcommunicator */
342   ierr = PetscMalloc1(1+nsubcomm,&subsize);CHKERRQ(ierr);
343 
344   np_subcomm = size/nsubcomm;
345   nleftover  = size - nsubcomm*np_subcomm;
346   for (i=0; i<nsubcomm; i++) {
347     subsize[i] = np_subcomm;
348     if (i<nleftover) subsize[i]++;
349   }
350 
351   /* get color and subrank of this proc */
352   rankstart = 0;
353   for (i=0; i<nsubcomm; i++) {
354     if (rank >= rankstart && rank < rankstart+subsize[i]) {
355       color   = i;
356       subrank = rank - rankstart;
357       duprank = rank;
358       break;
359     } else rankstart += subsize[i];
360   }
361 
362   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
363 
364   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
365   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
366   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
367   ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr);
368   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
369   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
370 
371   psubcomm->color   = color;
372   psubcomm->subsize = subsize;
373   psubcomm->type    = PETSC_SUBCOMM_CONTIGUOUS;
374   PetscFunctionReturn(0);
375 }
376 
377 #undef __FUNCT__
378 #define __FUNCT__ "PetscSubcommCreate_interlaced"
379 /*
380    Note:
381    In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators
382    by iteratively taking a process into a subcommunicator.
383    Example: size=4, nsubcomm=(*psubcomm)->n=3
384      comm=(*psubcomm)->parent:
385       rank:     [0]  [1]  [2]  [3]
386       color:     0    1    2    0
387 
388      subcomm=(*psubcomm)->comm:
389       subrank:  [0]  [0]  [0]  [1]
390 
391      dupcomm=(*psubcomm)->dupparent:
392       duprank:  [0]  [2]  [3]  [1]
393 
394      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
395            subcomm[color = 1] has subsize=1, owns process [1]
396            subcomm[color = 2] has subsize=1, owns process [2]
397            dupcomm has same number of processes as comm, and its duprank maps
398            processes in subcomm contiguously into a 1d array:
399             duprank: [0] [1]      [2]         [3]
400             rank:    [0] [3]      [1]         [2]
401                     subcomm[0] subcomm[1]  subcomm[2]
402 */
403 
404 PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm)
405 {
406   PetscErrorCode ierr;
407   PetscMPIInt    rank,size,*subsize,duprank,subrank;
408   PetscMPIInt    np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n;
409   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
410 
411   PetscFunctionBegin;
412   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
413   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
414 
415   /* get size of each subcommunicator */
416   ierr = PetscMalloc1(1+nsubcomm,&subsize);CHKERRQ(ierr);
417 
418   np_subcomm = size/nsubcomm;
419   nleftover  = size - nsubcomm*np_subcomm;
420   for (i=0; i<nsubcomm; i++) {
421     subsize[i] = np_subcomm;
422     if (i<nleftover) subsize[i]++;
423   }
424 
425   /* find color for this proc */
426   color   = rank%nsubcomm;
427   subrank = rank/nsubcomm;
428 
429   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
430 
431   j = 0; duprank = 0;
432   for (i=0; i<nsubcomm; i++) {
433     if (j == color) {
434       duprank += subrank;
435       break;
436     }
437     duprank += subsize[i]; j++;
438   }
439 
440   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
441   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
442   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,NULL);CHKERRQ(ierr);
443   ierr = PetscCommDuplicate(subcomm,&psubcomm->child,NULL);CHKERRQ(ierr);
444   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
445   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
446 
447   psubcomm->color   = color;
448   psubcomm->subsize = subsize;
449   psubcomm->type    = PETSC_SUBCOMM_INTERLACED;
450   PetscFunctionReturn(0);
451 }
452 
453 
454