xref: /petsc/src/sys/objects/subcomm.c (revision 17f8689b568b4cff545153ceaf0f1afe56c73cbb)
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_interlaced"
66 PetscErrorCode PetscSubcommCreate_contiguous(MPI_Comm comm,PetscInt nsubcomm,PetscSubcomm *psubcomm)
67 {
68   PetscErrorCode ierr;
69   PetscMPIInt    rank,size,*subsize,duprank,subrank;
70   PetscInt       np_subcomm,nleftover,i,color,rankstart;
71   MPI_Comm       subcomm=0,dupcomm=0;
72   PetscSubcomm   psubcomm_tmp;
73 
74   PetscFunctionBegin;
75   /* get size of each subcommunicator */
76   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
77   np_subcomm = size/nsubcomm;
78   nleftover  = size - nsubcomm*np_subcomm;
79   for (i=0; i<nsubcomm; i++){
80     subsize[i] = np_subcomm;
81     if (i<nleftover) subsize[i]++;
82   }
83 
84   /* get color and subrank of this proc */
85   rankstart = 0;
86   for (i=0; i<nsubcomm; i++){
87     if ( rank >= rankstart && rank < rankstart+subsize[i]) {
88       color   = i;
89       subrank = rank - rankstart;
90       duprank = rank;
91       break;
92     } else {
93       rankstart += subsize[i];
94     }
95   }
96   ierr = PetscSynchronizedPrintf(comm,"[%d] color %d, sub-rank %d\n",rank,color,subrank);
97   ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr);
98   ierr = PetscFree(subsize);CHKERRQ(ierr);
99 
100   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
101 
102   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
103   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
104 
105   ierr = PetscNew(struct _n_PetscSubcomm,&psubcomm_tmp);CHKERRQ(ierr);
106   psubcomm_tmp->parent    = comm;
107   psubcomm_tmp->dupparent = dupcomm;
108   psubcomm_tmp->comm      = subcomm;
109   psubcomm_tmp->n         = nsubcomm;
110   psubcomm_tmp->color     = color;
111   *psubcomm = psubcomm_tmp;
112   PetscFunctionReturn(0);
113 }
114 
115 #undef __FUNCT__
116 #define __FUNCT__ "PetscSubcommCreate_interlaced"
117 /*
118    Note:
119    In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators
120    by iteratively taking a process into a subcommunicator.
121    Example: size=4, nsubcomm=(*psubcomm)->n=3
122      comm=(*psubcomm)->parent:
123       rank:     [0]  [1]  [2]  [3]
124       color:     0    1    2    0
125 
126      subcomm=(*psubcomm)->comm:
127       subrank:  [0]  [0]  [0]  [1]
128 
129      dupcomm=(*psubcomm)->dupparent:
130       duprank:  [0]  [2]  [3]  [1]
131 
132      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
133            subcomm[color = 1] has subsize=1, owns process [1]
134            subcomm[color = 2] has subsize=1, owns process [2]
135            dupcomm has same number of processes as comm, and its duprank maps
136            processes in subcomm contiguously into a 1d array:
137             duprank: [0] [1]      [2]         [3]
138             rank:    [0] [3]      [1]         [2]
139                     subcomm[0] subcomm[1]  subcomm[2]
140 */
141 
142 PetscErrorCode PetscSubcommCreate_interlaced(MPI_Comm comm,PetscInt nsubcomm,PetscSubcomm *psubcomm)
143 {
144   PetscErrorCode ierr;
145   PetscMPIInt    rank,size,*subsize,duprank,subrank;
146   PetscInt       np_subcomm,nleftover,i,j,color;
147   MPI_Comm       subcomm=0,dupcomm=0;
148   PetscSubcomm   psubcomm_tmp;
149 
150   PetscFunctionBegin;
151   /* get size of each subcommunicator */
152   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
153   np_subcomm = size/nsubcomm;
154   nleftover  = size - nsubcomm*np_subcomm;
155   for (i=0; i<nsubcomm; i++){
156     subsize[i] = np_subcomm;
157     if (i<nleftover) subsize[i]++;
158   }
159 
160   /* find color for this proc */
161   color   = rank%nsubcomm;
162   subrank = rank/nsubcomm;
163 
164   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
165 
166   j = 0; duprank = 0;
167   for (i=0; i<nsubcomm; i++){
168     if (j == color){
169       duprank += subrank;
170       break;
171     }
172     duprank += subsize[i]; j++;
173   }
174   ierr = PetscFree(subsize);CHKERRQ(ierr);
175 
176   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
177   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
178 
179   ierr = PetscNew(struct _n_PetscSubcomm,&psubcomm_tmp);CHKERRQ(ierr);
180   psubcomm_tmp->parent    = comm;
181   psubcomm_tmp->dupparent = dupcomm;
182   psubcomm_tmp->comm      = subcomm;
183   psubcomm_tmp->n         = nsubcomm;
184   psubcomm_tmp->color     = color;
185   *psubcomm = psubcomm_tmp;
186   PetscFunctionReturn(0);
187 }
188 
189 
190