xref: /petsc/src/sys/objects/subcomm.c (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
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   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 
164   ierr = PetscNew(struct _n_PetscSubcomm,psubcomm);CHKERRQ(ierr);
165   (*psubcomm)->parent = comm;
166   PetscFunctionReturn(0);
167 }
168 
169 #undef __FUNCT__
170 #define __FUNCT__ "PetscSubcommCreate_contiguous"
171 PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm)
172 {
173   PetscErrorCode ierr;
174   PetscMPIInt    rank,size,*subsize,duprank=-1,subrank=-1;
175   PetscInt       np_subcomm,nleftover,i,color=-1,rankstart,nsubcomm=psubcomm->n;
176   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
177 
178   PetscFunctionBegin;
179   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
180   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
181 
182   /* get size of each subcommunicator */
183   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
184   np_subcomm = size/nsubcomm;
185   nleftover  = size - nsubcomm*np_subcomm;
186   for (i=0; i<nsubcomm; i++){
187     subsize[i] = np_subcomm;
188     if (i<nleftover) subsize[i]++;
189   }
190 
191   /* get color and subrank of this proc */
192   rankstart = 0;
193   for (i=0; i<nsubcomm; i++){
194     if ( rank >= rankstart && rank < rankstart+subsize[i]) {
195       color   = i;
196       subrank = rank - rankstart;
197       duprank = rank;
198       break;
199     } else {
200       rankstart += subsize[i];
201     }
202   }
203   ierr = PetscFree(subsize);CHKERRQ(ierr);
204 
205   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
206 
207   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
208   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
209 
210   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr);
211   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr);
212   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
213   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
214   psubcomm->color     = color;
215   PetscFunctionReturn(0);
216 }
217 
218 #undef __FUNCT__
219 #define __FUNCT__ "PetscSubcommCreate_interlaced"
220 /*
221    Note:
222    In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators
223    by iteratively taking a process into a subcommunicator.
224    Example: size=4, nsubcomm=(*psubcomm)->n=3
225      comm=(*psubcomm)->parent:
226       rank:     [0]  [1]  [2]  [3]
227       color:     0    1    2    0
228 
229      subcomm=(*psubcomm)->comm:
230       subrank:  [0]  [0]  [0]  [1]
231 
232      dupcomm=(*psubcomm)->dupparent:
233       duprank:  [0]  [2]  [3]  [1]
234 
235      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
236            subcomm[color = 1] has subsize=1, owns process [1]
237            subcomm[color = 2] has subsize=1, owns process [2]
238            dupcomm has same number of processes as comm, and its duprank maps
239            processes in subcomm contiguously into a 1d array:
240             duprank: [0] [1]      [2]         [3]
241             rank:    [0] [3]      [1]         [2]
242                     subcomm[0] subcomm[1]  subcomm[2]
243 */
244 
245 PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm)
246 {
247   PetscErrorCode ierr;
248   PetscMPIInt    rank,size,*subsize,duprank,subrank;
249   PetscInt       np_subcomm,nleftover,i,j,color,nsubcomm=psubcomm->n;
250   MPI_Comm       subcomm=0,dupcomm=0,comm=psubcomm->parent;
251 
252   PetscFunctionBegin;
253   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
254   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
255 
256   /* get size of each subcommunicator */
257   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
258   np_subcomm = size/nsubcomm;
259   nleftover  = size - nsubcomm*np_subcomm;
260   for (i=0; i<nsubcomm; i++){
261     subsize[i] = np_subcomm;
262     if (i<nleftover) subsize[i]++;
263   }
264 
265   /* find color for this proc */
266   color   = rank%nsubcomm;
267   subrank = rank/nsubcomm;
268 
269   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
270 
271   j = 0; duprank = 0;
272   for (i=0; i<nsubcomm; i++){
273     if (j == color){
274       duprank += subrank;
275       break;
276     }
277     duprank += subsize[i]; j++;
278   }
279   ierr = PetscFree(subsize);CHKERRQ(ierr);
280 
281   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
282   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
283 
284   ierr = PetscCommDuplicate(dupcomm,&psubcomm->dupparent,PETSC_NULL);CHKERRQ(ierr);
285   ierr = PetscCommDuplicate(subcomm,&psubcomm->comm,PETSC_NULL);CHKERRQ(ierr);
286   ierr = MPI_Comm_free(&dupcomm);CHKERRQ(ierr);
287   ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr);
288   psubcomm->color     = color;
289   PetscFunctionReturn(0);
290 }
291 
292 
293