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