xref: /petsc/src/sys/objects/subcomm.c (revision ace3abfce6818a63371c0ee8fdb525d96915299c)
1 #define PETSC_DLL
2 /*
3      Provides utility routines for split MPI communicator.
4 */
5 #include "petscsys.h"    /*I   "petscsys.h"    I*/
6 
7 extern PetscErrorCode PetscSubcommCreate_contiguous(MPI_Comm,PetscInt,PetscSubcomm*);
8 extern PetscErrorCode PetscSubcommCreate_interlaced(MPI_Comm,PetscInt,PetscSubcomm*);
9 
10 #undef __FUNCT__
11 #define __FUNCT__ "PetscSubcommDestroy"
12 PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommDestroy(PetscSubcomm psubcomm)
13 {
14   PetscErrorCode ierr;
15 
16   PetscFunctionBegin;
17   ierr = MPI_Comm_free(&psubcomm->dupparent);CHKERRQ(ierr);
18   ierr = MPI_Comm_free(&psubcomm->comm);CHKERRQ(ierr);
19   ierr = PetscFree(psubcomm);CHKERRQ(ierr);
20   PetscFunctionReturn(0);
21 }
22 
23 #undef __FUNCT__
24 #define __FUNCT__ "PetscSubcommCreate"
25 /*@C
26   PetscSubcommCreate - Create a PetscSubcomm context.
27 
28    Collective on MPI_Comm
29 
30    Input Parameter:
31 +  comm - MPI communicator
32 .  nsubcomm - the number of subcommunicators to be created
33 -  subcommtype - subcommunicator type
34 
35    Output Parameter:
36 .  psubcomm - location to store the PetscSubcomm context
37 
38    Level: advanced
39 
40 .keywords: communicator, create
41 
42 .seealso: PetscSubcommDestroy()
43 @*/
44 PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommCreate(MPI_Comm comm,PetscInt nsubcomm,PetscSubcommType subcommtype,PetscSubcomm *psubcomm)
45 {
46   PetscErrorCode ierr;
47   PetscMPIInt    rank,size;
48 
49   PetscFunctionBegin;
50   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
51   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
52   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);
53 
54   if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS){
55     ierr = PetscSubcommCreate_contiguous(comm,nsubcomm,psubcomm);CHKERRQ(ierr);
56   } else if (subcommtype == PETSC_SUBCOMM_INTERLACED){
57     ierr = PetscSubcommCreate_interlaced(comm,nsubcomm,psubcomm);CHKERRQ(ierr);
58   } else {
59     SETERRQ1(comm,PETSC_ERR_SUP,"PetscSubcommType %D is not supported yet",subcommtype);
60   }
61   PetscFunctionReturn(0);
62 }
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "PetscSubcommCreate_contiguous"
66 PetscErrorCode PetscSubcommCreate_contiguous(MPI_Comm comm,PetscInt nsubcomm,PetscSubcomm *psubcomm)
67 {
68   PetscErrorCode ierr;
69   PetscMPIInt    rank,size,*subsize,duprank=-1,subrank=-1;
70   PetscInt       np_subcomm,nleftover,i,color=-1,rankstart;
71   MPI_Comm       subcomm=0,dupcomm=0;
72   PetscSubcomm   psubcomm_tmp;
73 
74   PetscFunctionBegin;
75   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
76   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
77 
78   /* get size of each subcommunicator */
79   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
80   np_subcomm = size/nsubcomm;
81   nleftover  = size - nsubcomm*np_subcomm;
82   for (i=0; i<nsubcomm; i++){
83     subsize[i] = np_subcomm;
84     if (i<nleftover) subsize[i]++;
85   }
86 
87   /* get color and subrank of this proc */
88   rankstart = 0;
89   for (i=0; i<nsubcomm; i++){
90     if ( rank >= rankstart && rank < rankstart+subsize[i]) {
91       color   = i;
92       subrank = rank - rankstart;
93       duprank = rank;
94       break;
95     } else {
96       rankstart += subsize[i];
97     }
98   }
99   ierr = PetscFree(subsize);CHKERRQ(ierr);
100 
101   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
102 
103   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
104   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
105 
106   ierr = PetscNew(struct _n_PetscSubcomm,&psubcomm_tmp);CHKERRQ(ierr);
107   psubcomm_tmp->parent    = comm;
108   psubcomm_tmp->dupparent = dupcomm;
109   psubcomm_tmp->comm      = subcomm;
110   psubcomm_tmp->n         = nsubcomm;
111   psubcomm_tmp->color     = color;
112   *psubcomm = psubcomm_tmp;
113   PetscFunctionReturn(0);
114 }
115 
116 #undef __FUNCT__
117 #define __FUNCT__ "PetscSubcommCreate_interlaced"
118 /*
119    Note:
120    In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators
121    by iteratively taking a process into a subcommunicator.
122    Example: size=4, nsubcomm=(*psubcomm)->n=3
123      comm=(*psubcomm)->parent:
124       rank:     [0]  [1]  [2]  [3]
125       color:     0    1    2    0
126 
127      subcomm=(*psubcomm)->comm:
128       subrank:  [0]  [0]  [0]  [1]
129 
130      dupcomm=(*psubcomm)->dupparent:
131       duprank:  [0]  [2]  [3]  [1]
132 
133      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
134            subcomm[color = 1] has subsize=1, owns process [1]
135            subcomm[color = 2] has subsize=1, owns process [2]
136            dupcomm has same number of processes as comm, and its duprank maps
137            processes in subcomm contiguously into a 1d array:
138             duprank: [0] [1]      [2]         [3]
139             rank:    [0] [3]      [1]         [2]
140                     subcomm[0] subcomm[1]  subcomm[2]
141 */
142 
143 PetscErrorCode PetscSubcommCreate_interlaced(MPI_Comm comm,PetscInt nsubcomm,PetscSubcomm *psubcomm)
144 {
145   PetscErrorCode ierr;
146   PetscMPIInt    rank,size,*subsize,duprank,subrank;
147   PetscInt       np_subcomm,nleftover,i,j,color;
148   MPI_Comm       subcomm=0,dupcomm=0;
149   PetscSubcomm   psubcomm_tmp;
150 
151   PetscFunctionBegin;
152   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
153   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
154 
155   /* get size of each subcommunicator */
156   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
157   np_subcomm = size/nsubcomm;
158   nleftover  = size - nsubcomm*np_subcomm;
159   for (i=0; i<nsubcomm; i++){
160     subsize[i] = np_subcomm;
161     if (i<nleftover) subsize[i]++;
162   }
163 
164   /* find color for this proc */
165   color   = rank%nsubcomm;
166   subrank = rank/nsubcomm;
167 
168   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
169 
170   j = 0; duprank = 0;
171   for (i=0; i<nsubcomm; i++){
172     if (j == color){
173       duprank += subrank;
174       break;
175     }
176     duprank += subsize[i]; j++;
177   }
178   ierr = PetscFree(subsize);CHKERRQ(ierr);
179 
180   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
181   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
182 
183   ierr = PetscNew(struct _n_PetscSubcomm,&psubcomm_tmp);CHKERRQ(ierr);
184   psubcomm_tmp->parent    = comm;
185   psubcomm_tmp->dupparent = dupcomm;
186   psubcomm_tmp->comm      = subcomm;
187   psubcomm_tmp->n         = nsubcomm;
188   psubcomm_tmp->color     = color;
189   *psubcomm = psubcomm_tmp;
190   PetscFunctionReturn(0);
191 }
192 
193 
194