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