xref: /petsc/src/vec/is/sf/impls/basic/alltoall/sfalltoall.c (revision ffc4695bcb29f4b022f59a5fd6bc99fc280ff6d8)
1dd5b3ca6SJunchao Zhang #include <../src/vec/is/sf/impls/basic/allgatherv/sfallgatherv.h>
2dd5b3ca6SJunchao Zhang #include <../src/vec/is/sf/impls/basic/allgather/sfallgather.h>
3dd5b3ca6SJunchao Zhang #include <../src/vec/is/sf/impls/basic/gatherv/sfgatherv.h>
4dd5b3ca6SJunchao Zhang 
5dd5b3ca6SJunchao Zhang /* Reuse the type. The difference is some fields (i.e., displs, recvcounts) are not used, which is not a big deal */
6dd5b3ca6SJunchao Zhang typedef PetscSF_Allgatherv PetscSF_Alltoall;
7dd5b3ca6SJunchao Zhang 
8dd5b3ca6SJunchao Zhang /*===================================================================================*/
9dd5b3ca6SJunchao Zhang /*              Implementations of SF public APIs                                    */
10dd5b3ca6SJunchao Zhang /*===================================================================================*/
11dd5b3ca6SJunchao Zhang static PetscErrorCode PetscSFGetGraph_Alltoall(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
12dd5b3ca6SJunchao Zhang {
13dd5b3ca6SJunchao Zhang   PetscErrorCode ierr;
14dd5b3ca6SJunchao Zhang   PetscInt       i;
15dd5b3ca6SJunchao Zhang 
16dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
17dd5b3ca6SJunchao Zhang   if (nroots)  *nroots  = sf->nroots;
18dd5b3ca6SJunchao Zhang   if (nleaves) *nleaves = sf->nleaves;
19dd5b3ca6SJunchao Zhang   if (ilocal)  *ilocal  = NULL; /* Contiguous local indices */
20dd5b3ca6SJunchao Zhang   if (iremote) {
21dd5b3ca6SJunchao Zhang     if (!sf->remote) {
22dd5b3ca6SJunchao Zhang       ierr = PetscMalloc1(sf->nleaves,&sf->remote);CHKERRQ(ierr);
23dd5b3ca6SJunchao Zhang       sf->remote_alloc = sf->remote;
24dd5b3ca6SJunchao Zhang       for (i=0; i<sf->nleaves; i++) {
25dd5b3ca6SJunchao Zhang         sf->remote[i].rank  = i;
26dd5b3ca6SJunchao Zhang         sf->remote[i].index = i;
27dd5b3ca6SJunchao Zhang       }
28dd5b3ca6SJunchao Zhang     }
29dd5b3ca6SJunchao Zhang     *iremote = sf->remote;
30dd5b3ca6SJunchao Zhang   }
31dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
32dd5b3ca6SJunchao Zhang }
33dd5b3ca6SJunchao Zhang 
34eb02082bSJunchao Zhang static PetscErrorCode PetscSFBcastAndOpBegin_Alltoall(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,const void *rootdata,PetscMemType leafmtype,void *leafdata,MPI_Op op)
35dd5b3ca6SJunchao Zhang {
36dd5b3ca6SJunchao Zhang   PetscErrorCode       ierr;
37cd620004SJunchao Zhang   PetscSFLink          link;
38dd5b3ca6SJunchao Zhang   MPI_Comm             comm;
39cd620004SJunchao Zhang   void                 *rootbuf = NULL,*leafbuf = NULL; /* buffer used by MPI */
40cd620004SJunchao Zhang   MPI_Request          *req;
41dd5b3ca6SJunchao Zhang 
42dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
43cd620004SJunchao Zhang   ierr = PetscSFLinkCreate(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op,PETSCSF_BCAST,&link);CHKERRQ(ierr);
44cd620004SJunchao Zhang   ierr = PetscSFLinkPackRootData(sf,link,PETSCSF_REMOTE,rootdata);CHKERRQ(ierr);
45dd5b3ca6SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
46cd620004SJunchao Zhang   ierr = PetscSFLinkGetMPIBuffersAndRequests(sf,link,PETSCSF_ROOT2LEAF,&rootbuf,&leafbuf,&req,NULL);CHKERRQ(ierr);
47cd620004SJunchao Zhang   ierr = MPIU_Ialltoall(rootbuf,1,unit,leafbuf,1,unit,comm,req);CHKERRQ(ierr);
48dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
49dd5b3ca6SJunchao Zhang }
50dd5b3ca6SJunchao Zhang 
51eb02082bSJunchao Zhang static PetscErrorCode PetscSFReduceBegin_Alltoall(PetscSF sf,MPI_Datatype unit,PetscMemType leafmtype,const void *leafdata,PetscMemType rootmtype,void *rootdata,MPI_Op op)
52dd5b3ca6SJunchao Zhang {
53dd5b3ca6SJunchao Zhang   PetscErrorCode       ierr;
54cd620004SJunchao Zhang   PetscSFLink          link;
55dd5b3ca6SJunchao Zhang   MPI_Comm             comm;
56cd620004SJunchao Zhang   void                 *rootbuf = NULL,*leafbuf = NULL; /* buffer used by MPI */
57cd620004SJunchao Zhang   MPI_Request          *req;
58dd5b3ca6SJunchao Zhang 
59dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
60cd620004SJunchao Zhang   ierr = PetscSFLinkCreate(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op,PETSCSF_REDUCE,&link);CHKERRQ(ierr);
61cd620004SJunchao Zhang   ierr = PetscSFLinkPackLeafData(sf,link,PETSCSF_REMOTE,leafdata);CHKERRQ(ierr);
62dd5b3ca6SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
63cd620004SJunchao Zhang   ierr = PetscSFLinkGetMPIBuffersAndRequests(sf,link,PETSCSF_LEAF2ROOT,&rootbuf,&leafbuf,&req,NULL);CHKERRQ(ierr);
64cd620004SJunchao Zhang   ierr = MPIU_Ialltoall(leafbuf,1,unit,rootbuf,1,unit,comm,req);CHKERRQ(ierr);
65dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
66dd5b3ca6SJunchao Zhang }
67dd5b3ca6SJunchao Zhang 
68dd5b3ca6SJunchao Zhang static PetscErrorCode PetscSFCreateLocalSF_Alltoall(PetscSF sf,PetscSF *out)
69dd5b3ca6SJunchao Zhang {
70dd5b3ca6SJunchao Zhang   PetscErrorCode ierr;
71dd5b3ca6SJunchao Zhang   PetscInt       nroots = 1,nleaves = 1,*ilocal;
72dd5b3ca6SJunchao Zhang   PetscSFNode    *iremote = NULL;
73dd5b3ca6SJunchao Zhang   PetscSF        lsf;
74dd5b3ca6SJunchao Zhang   PetscMPIInt    rank;
75dd5b3ca6SJunchao Zhang 
76dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
77dd5b3ca6SJunchao Zhang   nroots  = 1;
78dd5b3ca6SJunchao Zhang   nleaves = 1;
79*ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRMPI(ierr);
80dd5b3ca6SJunchao Zhang   ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr);
81dd5b3ca6SJunchao Zhang   ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr);
82dd5b3ca6SJunchao Zhang   ilocal[0]        = rank;
83dd5b3ca6SJunchao Zhang   iremote[0].rank  = 0;    /* rank in PETSC_COMM_SELF */
84dd5b3ca6SJunchao Zhang   iremote[0].index = rank; /* LocalSF is an embedded SF. Indices are not remapped */
85dd5b3ca6SJunchao Zhang 
86dd5b3ca6SJunchao Zhang   ierr = PetscSFCreate(PETSC_COMM_SELF,&lsf);CHKERRQ(ierr);
87dd5b3ca6SJunchao Zhang   ierr = PetscSFSetGraph(lsf,nroots,nleaves,NULL/*contiguous leaves*/,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
88dd5b3ca6SJunchao Zhang   ierr = PetscSFSetUp(lsf);CHKERRQ(ierr);
89dd5b3ca6SJunchao Zhang   *out = lsf;
90dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
91dd5b3ca6SJunchao Zhang }
92dd5b3ca6SJunchao Zhang 
93dd5b3ca6SJunchao Zhang static PetscErrorCode PetscSFCreateEmbeddedSF_Alltoall(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
94dd5b3ca6SJunchao Zhang {
95dd5b3ca6SJunchao Zhang   PetscErrorCode ierr;
96dd5b3ca6SJunchao Zhang   PetscInt       i,*tmproots,*ilocal,ndranks,ndiranks;
97dd5b3ca6SJunchao Zhang   PetscSFNode    *iremote;
98dd5b3ca6SJunchao Zhang   PetscMPIInt    nroots,*roots,nleaves,*leaves,rank;
99dd5b3ca6SJunchao Zhang   MPI_Comm       comm;
100dd5b3ca6SJunchao Zhang   PetscSF_Basic  *bas;
101dd5b3ca6SJunchao Zhang   PetscSF        esf;
102dd5b3ca6SJunchao Zhang 
103dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
104dd5b3ca6SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
105*ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
106dd5b3ca6SJunchao Zhang 
107dd5b3ca6SJunchao Zhang   /* Uniq selected[] and store the result in roots[] */
108dd5b3ca6SJunchao Zhang   ierr = PetscMalloc1(nselected,&tmproots);CHKERRQ(ierr);
109da2e4c71SJunchao Zhang   ierr = PetscArraycpy(tmproots,selected,nselected);CHKERRQ(ierr);
110dd5b3ca6SJunchao Zhang   ierr = PetscSortRemoveDupsInt(&nselected,tmproots);CHKERRQ(ierr); /* nselected might be changed */
111dd5b3ca6SJunchao Zhang   if (tmproots[0] < 0 || tmproots[nselected-1] >= sf->nroots) SETERRQ3(comm,PETSC_ERR_ARG_OUTOFRANGE,"Min/Max root indices %D/%D are not in [0,%D)",tmproots[0],tmproots[nselected-1],sf->nroots);
112dd5b3ca6SJunchao Zhang   nroots = nselected;   /* For Alltoall, we know root indices will not overflow MPI_INT */
113dd5b3ca6SJunchao Zhang   ierr   = PetscMalloc1(nselected,&roots);CHKERRQ(ierr);
114dd5b3ca6SJunchao Zhang   for (i=0; i<nselected; i++) roots[i] = tmproots[i];
115dd5b3ca6SJunchao Zhang   ierr   = PetscFree(tmproots);CHKERRQ(ierr);
116dd5b3ca6SJunchao Zhang 
117dd5b3ca6SJunchao Zhang   /* Find out which leaves are still connected to roots in the embedded sf. Expect PetscCommBuildTwoSided is more scalable than MPI_Alltoall */
118dd5b3ca6SJunchao Zhang   ierr   = PetscCommBuildTwoSided(comm,0/*empty msg*/,MPI_INT/*fake*/,nroots,roots,NULL/*todata*/,&nleaves,&leaves,NULL/*fromdata*/);CHKERRQ(ierr);
119dd5b3ca6SJunchao Zhang 
120dd5b3ca6SJunchao Zhang   /* Move myself ahead if rank is in leaves[], since I am a distinguished rank */
121dd5b3ca6SJunchao Zhang   ndranks = 0;
122dd5b3ca6SJunchao Zhang   for (i=0; i<nleaves; i++) {
123dd5b3ca6SJunchao Zhang     if (leaves[i] == rank) {leaves[i] = -rank; ndranks = 1; break;}
124dd5b3ca6SJunchao Zhang   }
125dd5b3ca6SJunchao Zhang   ierr = PetscSortMPIInt(nleaves,leaves);CHKERRQ(ierr);
126dd5b3ca6SJunchao Zhang   if (nleaves && leaves[0] < 0) leaves[0] = rank;
127dd5b3ca6SJunchao Zhang 
128dd5b3ca6SJunchao Zhang   /* Build esf and fill its fields manually (without calling PetscSFSetUp) */
129dd5b3ca6SJunchao Zhang   ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr);
130dd5b3ca6SJunchao Zhang   ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr);
131dd5b3ca6SJunchao Zhang   for (i=0; i<nleaves; i++) { /* 1:1 map from roots to leaves */
132dd5b3ca6SJunchao Zhang     ilocal[i]        = leaves[i];
133dd5b3ca6SJunchao Zhang     iremote[i].rank  = leaves[i];
134dd5b3ca6SJunchao Zhang     iremote[i].index = leaves[i];
135dd5b3ca6SJunchao Zhang   }
136dd5b3ca6SJunchao Zhang   ierr = PetscSFCreate(comm,&esf);CHKERRQ(ierr);
137dd5b3ca6SJunchao Zhang   ierr = PetscSFSetType(esf,PETSCSFBASIC);CHKERRQ(ierr); /* This optimized routine can only create a basic sf */
138dd5b3ca6SJunchao Zhang   ierr = PetscSFSetGraph(esf,sf->nleaves,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
139dd5b3ca6SJunchao Zhang 
140dd5b3ca6SJunchao Zhang   /* As if we are calling PetscSFSetUpRanks(esf,self's group) */
141dd5b3ca6SJunchao Zhang   ierr = PetscMalloc4(nleaves,&esf->ranks,nleaves+1,&esf->roffset,nleaves,&esf->rmine,nleaves,&esf->rremote);CHKERRQ(ierr);
142dd5b3ca6SJunchao Zhang   esf->nranks     = nleaves;
143dd5b3ca6SJunchao Zhang   esf->ndranks    = ndranks;
144dd5b3ca6SJunchao Zhang   esf->roffset[0] = 0;
145dd5b3ca6SJunchao Zhang   for (i=0; i<nleaves; i++) {
146dd5b3ca6SJunchao Zhang     esf->ranks[i]     = leaves[i];
147dd5b3ca6SJunchao Zhang     esf->roffset[i+1] = i+1;
148dd5b3ca6SJunchao Zhang     esf->rmine[i]     = leaves[i];
149dd5b3ca6SJunchao Zhang     esf->rremote[i]   = leaves[i];
150dd5b3ca6SJunchao Zhang   }
151dd5b3ca6SJunchao Zhang 
152dd5b3ca6SJunchao Zhang   /* Set up esf->data, the incoming communication (i.e., recv info), which is usually done by PetscSFSetUp_Basic */
153dd5b3ca6SJunchao Zhang   bas  = (PetscSF_Basic*)esf->data;
154dd5b3ca6SJunchao Zhang   ierr = PetscMalloc2(nroots,&bas->iranks,nroots+1,&bas->ioffset);CHKERRQ(ierr);
155dd5b3ca6SJunchao Zhang   ierr = PetscMalloc1(nroots,&bas->irootloc);CHKERRQ(ierr);
156dd5b3ca6SJunchao Zhang   /* Move myself ahead if rank is in roots[], since I am a distinguished irank */
157dd5b3ca6SJunchao Zhang   ndiranks = 0;
158dd5b3ca6SJunchao Zhang   for (i=0; i<nroots; i++) {
159dd5b3ca6SJunchao Zhang     if (roots[i] == rank) {roots[i] = -rank; ndiranks = 1; break;}
160dd5b3ca6SJunchao Zhang   }
161dd5b3ca6SJunchao Zhang   ierr = PetscSortMPIInt(nroots,roots);CHKERRQ(ierr);
162dd5b3ca6SJunchao Zhang   if (nroots && roots[0] < 0) roots[0] = rank;
163dd5b3ca6SJunchao Zhang 
164dd5b3ca6SJunchao Zhang   bas->niranks    = nroots;
165dd5b3ca6SJunchao Zhang   bas->ndiranks   = ndiranks;
166dd5b3ca6SJunchao Zhang   bas->ioffset[0] = 0;
167dd5b3ca6SJunchao Zhang   bas->itotal     = nroots;
168dd5b3ca6SJunchao Zhang   for (i=0; i<nroots; i++) {
169dd5b3ca6SJunchao Zhang     bas->iranks[i]    = roots[i];
170dd5b3ca6SJunchao Zhang     bas->ioffset[i+1] = i+1;
171dd5b3ca6SJunchao Zhang     bas->irootloc[i]  = roots[i];
172dd5b3ca6SJunchao Zhang   }
173dd5b3ca6SJunchao Zhang 
174cd620004SJunchao Zhang   /* See PetscSFCreateEmbeddedSF_Basic */
175cd620004SJunchao Zhang   esf->nleafreqs  = esf->nranks - esf->ndranks;
176cd620004SJunchao Zhang   bas->nrootreqs  = bas->niranks - bas->ndiranks;
177cd620004SJunchao Zhang   esf->persistent = PETSC_TRUE;
178cd620004SJunchao Zhang   /* Setup packing related fields */
179cd620004SJunchao Zhang   ierr = PetscSFSetUpPackFields(esf);CHKERRQ(ierr);
180cd620004SJunchao Zhang 
181cd620004SJunchao Zhang   esf->setupcalled = PETSC_TRUE; /* We have done setup ourselves! */
182dd5b3ca6SJunchao Zhang   *newsf = esf;
183dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
184dd5b3ca6SJunchao Zhang }
185dd5b3ca6SJunchao Zhang 
186dd5b3ca6SJunchao Zhang PETSC_INTERN PetscErrorCode PetscSFCreate_Alltoall(PetscSF sf)
187dd5b3ca6SJunchao Zhang {
188dd5b3ca6SJunchao Zhang   PetscErrorCode   ierr;
189dd5b3ca6SJunchao Zhang   PetscSF_Alltoall *dat = (PetscSF_Alltoall*)sf->data;
190dd5b3ca6SJunchao Zhang 
191dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
192cd620004SJunchao Zhang   sf->ops->BcastAndOpEnd   = PetscSFBcastAndOpEnd_Basic;
193cd620004SJunchao Zhang   sf->ops->ReduceEnd       = PetscSFReduceEnd_Basic;
194cd620004SJunchao Zhang 
195dd5b3ca6SJunchao Zhang   /* Inherit from Allgatherv. It is astonishing Alltoall can inherit so much from Allgather(v) */
196dd5b3ca6SJunchao Zhang   sf->ops->Destroy         = PetscSFDestroy_Allgatherv;
197dd5b3ca6SJunchao Zhang   sf->ops->Reset           = PetscSFReset_Allgatherv;
198dd5b3ca6SJunchao Zhang   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Allgatherv;
199dd5b3ca6SJunchao Zhang   sf->ops->GetRootRanks    = PetscSFGetRootRanks_Allgatherv;
200dd5b3ca6SJunchao Zhang 
201dd5b3ca6SJunchao Zhang   /* Inherit from Allgather. Every process gathers equal-sized data from others, which enables this inheritance. */
202dd5b3ca6SJunchao Zhang   sf->ops->GetLeafRanks    = PetscSFGetLeafRanks_Allgatherv;
203cd620004SJunchao Zhang   sf->ops->SetUp           = PetscSFSetUp_Allgather;
204dd5b3ca6SJunchao Zhang 
205dd5b3ca6SJunchao Zhang   /* Inherit from Gatherv. Each root has only one leaf connected, which enables this inheritance */
206dd5b3ca6SJunchao Zhang   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Gatherv;
207dd5b3ca6SJunchao Zhang 
208dd5b3ca6SJunchao Zhang   /* Alltoall stuff */
209dd5b3ca6SJunchao Zhang   sf->ops->GetGraph         = PetscSFGetGraph_Alltoall;
210dd5b3ca6SJunchao Zhang   sf->ops->BcastAndOpBegin  = PetscSFBcastAndOpBegin_Alltoall;
211dd5b3ca6SJunchao Zhang   sf->ops->ReduceBegin      = PetscSFReduceBegin_Alltoall;
212dd5b3ca6SJunchao Zhang   sf->ops->CreateLocalSF    = PetscSFCreateLocalSF_Alltoall;
213dd5b3ca6SJunchao Zhang   sf->ops->CreateEmbeddedSF = PetscSFCreateEmbeddedSF_Alltoall;
214dd5b3ca6SJunchao Zhang 
215dd5b3ca6SJunchao Zhang   ierr = PetscNewLog(sf,&dat);CHKERRQ(ierr);
216dd5b3ca6SJunchao Zhang   sf->data = (void*)dat;
217dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
218dd5b3ca6SJunchao Zhang }
219dd5b3ca6SJunchao Zhang 
220dd5b3ca6SJunchao Zhang 
221