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