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