xref: /petsc/src/sys/objects/subcomm.c (revision cd05a4c019503193975def50f0ad7dec844973eb)
1*cd05a4c0SHong Zhang #define PETSC_DLL
2*cd05a4c0SHong Zhang /*
3*cd05a4c0SHong Zhang      Provides utility routines for split MPI communicator.
4*cd05a4c0SHong Zhang */
5*cd05a4c0SHong Zhang #include "petsc.h"  /*I   "petsc.h"    I*/
6*cd05a4c0SHong Zhang #include "petscsys.h"
7*cd05a4c0SHong Zhang 
8*cd05a4c0SHong Zhang #undef __FUNCT__
9*cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommDestroy"
10*cd05a4c0SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommDestroy(PetscSubcomm *psubcomm)
11*cd05a4c0SHong Zhang {
12*cd05a4c0SHong Zhang   PetscErrorCode ierr;
13*cd05a4c0SHong Zhang 
14*cd05a4c0SHong Zhang   PetscFunctionBegin;
15*cd05a4c0SHong Zhang   ierr = PetscFree(psubcomm);CHKERRQ(ierr);
16*cd05a4c0SHong Zhang   PetscFunctionReturn(0);
17*cd05a4c0SHong Zhang }
18*cd05a4c0SHong Zhang 
19*cd05a4c0SHong Zhang #undef __FUNCT__
20*cd05a4c0SHong Zhang #define __FUNCT__ "PetscSubcommCreate"
21*cd05a4c0SHong Zhang /*@
22*cd05a4c0SHong Zhang   PetscSubcommCreate - Create a PetscSubcomm context.
23*cd05a4c0SHong Zhang 
24*cd05a4c0SHong Zhang    Collective on MPI_Comm
25*cd05a4c0SHong Zhang 
26*cd05a4c0SHong Zhang    Input Parameter:
27*cd05a4c0SHong Zhang +  comm - MPI communicator
28*cd05a4c0SHong Zhang -  nsubcomm - the number of subcommunicators to be created
29*cd05a4c0SHong Zhang 
30*cd05a4c0SHong Zhang    Output Parameter:
31*cd05a4c0SHong Zhang .  psubcomm - location to store the PetscSubcomm context
32*cd05a4c0SHong Zhang 
33*cd05a4c0SHong Zhang    Options Database Keys:
34*cd05a4c0SHong Zhang 
35*cd05a4c0SHong Zhang    Notes:
36*cd05a4c0SHong Zhang    To avoid data scattering from subcomm back to original comm, we create subcommunicators
37*cd05a4c0SHong Zhang    by iteratively taking a processe into a subcommunicator.
38*cd05a4c0SHong Zhang    Example: size=4, nsubcomm=(*psubcomm)->n=3
39*cd05a4c0SHong Zhang      comm=(*psubcomm)->parent:
40*cd05a4c0SHong Zhang       rank:     [0]  [1]  [2]  [3]
41*cd05a4c0SHong Zhang       color:     0    1    2    0
42*cd05a4c0SHong Zhang 
43*cd05a4c0SHong Zhang      subcomm=(*psubcomm)->comm:
44*cd05a4c0SHong Zhang       subrank:  [0]  [0]  [0]  [1]
45*cd05a4c0SHong Zhang 
46*cd05a4c0SHong Zhang      dupcomm=(*psubcomm)->dupparent:
47*cd05a4c0SHong Zhang       duprank:  [0]  [2]  [3]  [1]
48*cd05a4c0SHong Zhang 
49*cd05a4c0SHong Zhang      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
50*cd05a4c0SHong Zhang            subcomm[color = 1] has subsize=1, owns process [1]
51*cd05a4c0SHong Zhang            subcomm[color = 2] has subsize=1, owns process [2]
52*cd05a4c0SHong Zhang            dupcomm has same number of processes as comm, and its duprank maps
53*cd05a4c0SHong Zhang            processes in subcomm contiguously into a 1d array:
54*cd05a4c0SHong Zhang             duprank: [0] [1]      [2]         [3]
55*cd05a4c0SHong Zhang             rank:    [0] [3]      [1]         [2]
56*cd05a4c0SHong Zhang                     subcomm[0] subcomm[1]  subcomm[2]
57*cd05a4c0SHong Zhang 
58*cd05a4c0SHong Zhang    Level: advanced
59*cd05a4c0SHong Zhang 
60*cd05a4c0SHong Zhang .keywords: communicator, create
61*cd05a4c0SHong Zhang 
62*cd05a4c0SHong Zhang .seealso: PetscSubcommDestroy()
63*cd05a4c0SHong Zhang @*/
64*cd05a4c0SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT PetscSubcommCreate(MPI_Comm comm,PetscInt nsubcomm,PetscSubcomm **psubcomm)
65*cd05a4c0SHong Zhang {
66*cd05a4c0SHong Zhang   PetscErrorCode ierr;
67*cd05a4c0SHong Zhang   PetscMPIInt    rank,size,*subsize,duprank,subrank;
68*cd05a4c0SHong Zhang   PetscInt       np_subcomm,nleftover,i,j,color;
69*cd05a4c0SHong Zhang   MPI_Comm       subcomm=0,dupcomm=0;
70*cd05a4c0SHong Zhang   const char     *prefix;
71*cd05a4c0SHong Zhang   PetscSubcomm   *psubcomm_tmp;
72*cd05a4c0SHong Zhang 
73*cd05a4c0SHong Zhang   PetscFunctionBegin;
74*cd05a4c0SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
75*cd05a4c0SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
76*cd05a4c0SHong Zhang   if (nsubcomm < 1 || nsubcomm > size) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Num of subcommunicators %D cannot be < 1 or > input comm size %D",nsubcomm,size);
77*cd05a4c0SHong Zhang 
78*cd05a4c0SHong Zhang   /* get size of each subcommunicator */
79*cd05a4c0SHong Zhang   ierr = PetscMalloc((1+nsubcomm)*sizeof(PetscMPIInt),&subsize);CHKERRQ(ierr);
80*cd05a4c0SHong Zhang   np_subcomm = size/nsubcomm;
81*cd05a4c0SHong Zhang   nleftover  = size - nsubcomm*np_subcomm;
82*cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++){
83*cd05a4c0SHong Zhang     subsize[i] = np_subcomm;
84*cd05a4c0SHong Zhang     if (i<nleftover) subsize[i]++;
85*cd05a4c0SHong Zhang   }
86*cd05a4c0SHong Zhang 
87*cd05a4c0SHong Zhang   /* find color for this proc */
88*cd05a4c0SHong Zhang   color   = rank%nsubcomm;
89*cd05a4c0SHong Zhang   subrank = rank/nsubcomm;
90*cd05a4c0SHong Zhang 
91*cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,color,subrank,&subcomm);CHKERRQ(ierr);
92*cd05a4c0SHong Zhang 
93*cd05a4c0SHong Zhang   j = 0; duprank = 0;
94*cd05a4c0SHong Zhang   for (i=0; i<nsubcomm; i++){
95*cd05a4c0SHong Zhang     if (j == color){
96*cd05a4c0SHong Zhang       duprank += subrank;
97*cd05a4c0SHong Zhang       break;
98*cd05a4c0SHong Zhang     }
99*cd05a4c0SHong Zhang     duprank += subsize[i]; j++;
100*cd05a4c0SHong Zhang   }
101*cd05a4c0SHong Zhang 
102*cd05a4c0SHong Zhang   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
103*cd05a4c0SHong Zhang   ierr = MPI_Comm_split(comm,0,duprank,&dupcomm);CHKERRQ(ierr);
104*cd05a4c0SHong Zhang   ierr = PetscFree(subsize);CHKERRQ(ierr);
105*cd05a4c0SHong Zhang 
106*cd05a4c0SHong Zhang   ierr = PetscNew(PetscSubcomm,&psubcomm_tmp);CHKERRQ(ierr);
107*cd05a4c0SHong Zhang   psubcomm_tmp->parent    = comm;
108*cd05a4c0SHong Zhang   psubcomm_tmp->dupparent = dupcomm;
109*cd05a4c0SHong Zhang   psubcomm_tmp->comm      = subcomm;
110*cd05a4c0SHong Zhang   psubcomm_tmp->n         = nsubcomm;
111*cd05a4c0SHong Zhang   psubcomm_tmp->color     = color;
112*cd05a4c0SHong Zhang   *psubcomm = psubcomm_tmp;
113*cd05a4c0SHong Zhang   PetscFunctionReturn(0);
114*cd05a4c0SHong Zhang }
115