xref: /petsc/src/vec/is/sf/impls/basic/allgatherv/sfallgatherv.c (revision b23bfdefca792e3d9f47034521b8d4c437693123)
1 #include <../src/vec/is/sf/impls/basic/allgatherv/sfallgatherv.h>
2 
3 PETSC_INTERN PetscErrorCode PetscSFBcastAndOpBegin_Gatherv(PetscSF,MPI_Datatype,const void*,void*,MPI_Op);
4 
5 /*===================================================================================*/
6 /*              Internal routines for PetscSFPack_Allgatherv                         */
7 /*===================================================================================*/
8 PETSC_INTERN PetscErrorCode PetscSFPackGet_Allgatherv(PetscSF sf,MPI_Datatype unit,const void *rkey,const void *lkey,PetscSFPack_Allgatherv *mylink)
9 {
10   PetscSF_Allgatherv     *dat = (PetscSF_Allgatherv*)sf->data;
11   PetscSFPack_Allgatherv link;
12   PetscSFPack            *p;
13   PetscErrorCode         ierr;
14 
15   PetscFunctionBegin;
16   ierr = PetscSFPackSetErrorOnUnsupportedOverlap(sf,unit,rkey,lkey);CHKERRQ(ierr);
17   /* Look for types in cache */
18   for (p=&dat->avail; (link=(PetscSFPack_Allgatherv)*p); p=&link->next) {
19     PetscBool match;
20     ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
21     if (match) {
22       *p = link->next; /* Remove from available list */
23       goto found;
24     }
25   }
26 
27   ierr = PetscNew(&link);CHKERRQ(ierr);
28   ierr = PetscSFPackSetupType((PetscSFPack)link,unit);CHKERRQ(ierr);
29   /* Must be init'ed to NULL so that we can safely wait on them even they are not used */
30   link->request = MPI_REQUEST_NULL;
31   /* One tag per link */
32   ierr = PetscCommGetNewTag(PetscObjectComm((PetscObject)sf),&link->tag);CHKERRQ(ierr);
33 
34   /* DO NOT allocate link->root/leaf. We use lazy allocation since these buffers are likely not needed */
35 found:
36   link->rkey = rkey;
37   link->lkey = lkey;
38   link->next = dat->inuse;
39   dat->inuse = (PetscSFPack)link;
40 
41   *mylink     = link;
42   PetscFunctionReturn(0);
43 }
44 
45 /*===================================================================================*/
46 /*              Implementations of SF public APIs                                    */
47 /*===================================================================================*/
48 
49 /* PetscSFGetGraph is non-collective. An implementation should not have collective calls */
50 PETSC_INTERN PetscErrorCode PetscSFGetGraph_Allgatherv(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
51 {
52   PetscErrorCode ierr;
53   PetscInt       i,j,k;
54   const PetscInt *range;
55   PetscMPIInt    size;
56 
57   PetscFunctionBegin;
58   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
59   if (nroots)  *nroots  = sf->nroots;
60   if (nleaves) *nleaves = sf->nleaves;
61   if (ilocal)  *ilocal  = NULL; /* Contiguous leaves */
62   if (iremote) {
63     if (!sf->remote && sf->nleaves) { /* The && sf->nleaves makes sfgatherv able to inherit this routine */
64       ierr = PetscLayoutGetRanges(sf->map,&range);CHKERRQ(ierr);
65       ierr = PetscMalloc1(sf->nleaves,&sf->remote);CHKERRQ(ierr);
66       sf->remote_alloc = sf->remote;
67       for (i=0; i<size; i++) {
68         for (j=range[i],k=0; j<range[i+1]; j++,k++) {
69           sf->remote[j].rank  = i;
70           sf->remote[j].index = k;
71         }
72       }
73     }
74     *iremote = sf->remote;
75   }
76   PetscFunctionReturn(0);
77 }
78 
79 PETSC_INTERN PetscErrorCode PetscSFSetUp_Allgatherv(PetscSF sf)
80 {
81   PetscErrorCode     ierr;
82   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv*)sf->data;
83   PetscMPIInt        size;
84   PetscInt           i;
85   const PetscInt     *range;
86 
87   PetscFunctionBegin;
88   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
89   if (sf->nleaves) { /* This if (sf->nleaves) test makes sfgatherv able to inherit this routine */
90     ierr = PetscMalloc1(size,&dat->recvcounts);CHKERRQ(ierr);
91     ierr = PetscMalloc1(size,&dat->displs);CHKERRQ(ierr);
92     ierr = PetscLayoutGetRanges(sf->map,&range);CHKERRQ(ierr);
93 
94     for (i=0; i<size; i++) {
95       ierr = PetscMPIIntCast(range[i],&dat->displs[i]);CHKERRQ(ierr);
96       ierr = PetscMPIIntCast(range[i+1]-range[i],&dat->recvcounts[i]);CHKERRQ(ierr);
97     }
98   }
99   PetscFunctionReturn(0);
100 }
101 
102 PETSC_INTERN PetscErrorCode PetscSFReset_Allgatherv(PetscSF sf)
103 {
104   PetscSF_Allgatherv     *dat = (PetscSF_Allgatherv*)sf->data;
105   PetscSFPack_Allgatherv link,next;
106   PetscErrorCode         ierr;
107 
108   PetscFunctionBegin;
109   ierr = PetscFree(dat->iranks);CHKERRQ(ierr);
110   ierr = PetscFree(dat->ioffset);CHKERRQ(ierr);
111   ierr = PetscFree(dat->irootloc);CHKERRQ(ierr);
112   ierr = PetscFree(dat->recvcounts);CHKERRQ(ierr);
113   ierr = PetscFree(dat->displs);CHKERRQ(ierr);
114   if (dat->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed");
115   for (link=(PetscSFPack_Allgatherv)dat->avail; link; link=next) {
116     next = (PetscSFPack_Allgatherv)link->next;
117     if (!link->isbuiltin) {ierr = MPI_Type_free(&link->unit);CHKERRQ(ierr);}
118     ierr = PetscFree(link->rootbuf);CHKERRQ(ierr);
119     ierr = PetscFree(link->leafbuf);CHKERRQ(ierr);
120     ierr = PetscFree(link);CHKERRQ(ierr);
121   }
122   dat->avail = NULL;
123   PetscFunctionReturn(0);
124 }
125 
126 PETSC_INTERN PetscErrorCode PetscSFDestroy_Allgatherv(PetscSF sf)
127 {
128   PetscErrorCode ierr;
129 
130   PetscFunctionBegin;
131   ierr = PetscSFReset_Allgatherv(sf);CHKERRQ(ierr);
132   ierr = PetscFree(sf->data);CHKERRQ(ierr);
133   PetscFunctionReturn(0);
134 }
135 
136 static PetscErrorCode PetscSFBcastAndOpBegin_Allgatherv(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
137 {
138   PetscErrorCode         ierr;
139   PetscSFPack_Allgatherv link;
140   PetscMPIInt            sendcount;
141   MPI_Comm               comm;
142   void                   *recvbuf;
143   PetscSF_Allgatherv     *dat = (PetscSF_Allgatherv*)sf->data;
144 
145   PetscFunctionBegin;
146   ierr = PetscSFPackGet_Allgatherv(sf,unit,rootdata,leafdata,&link);CHKERRQ(ierr);
147   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
148   ierr = PetscMPIIntCast(sf->nroots,&sendcount);CHKERRQ(ierr);
149 
150   if (op == MPIU_REPLACE) {
151     recvbuf = leafdata;
152   } else {
153     if (!link->leafbuf) {ierr = PetscMalloc(sf->nleaves*link->unitbytes,&link->leafbuf);CHKERRQ(ierr);}
154     recvbuf = link->leafbuf;
155   }
156   ierr = MPIU_Iallgatherv(rootdata,sendcount,unit,recvbuf,dat->recvcounts,dat->displs,unit,comm,&link->request);CHKERRQ(ierr);
157   PetscFunctionReturn(0);
158 }
159 
160 PETSC_INTERN PetscErrorCode PetscSFBcastAndOpEnd_Allgatherv(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
161 {
162   PetscErrorCode         ierr,(*UnpackAndOp)(PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt,void*,const void*);
163   PetscSFPack_Allgatherv link;
164   PetscMPIInt            count;
165 
166   PetscFunctionBegin;
167   ierr = PetscSFPackGetInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,(PetscSFPack*)&link);CHKERRQ(ierr);
168   ierr = MPI_Wait(&link->request,MPI_STATUS_IGNORE);CHKERRQ(ierr);
169 
170   if (op != MPIU_REPLACE) {
171     ierr = PetscSFPackGetUnpackAndOp(sf,(PetscSFPack)link,op,&UnpackAndOp);CHKERRQ(ierr);
172     if (UnpackAndOp) {ierr = (*UnpackAndOp)(sf->nleaves,NULL,link->bs,NULL,leafdata,link->leafbuf);CHKERRQ(ierr);}
173     else { /* op might be user-defined */
174       ierr = PetscMPIIntCast(sf->nleaves,&count);CHKERRQ(ierr);
175       ierr = MPI_Reduce_local(link->leafbuf,leafdata,count,unit,op);CHKERRQ(ierr);
176     }
177   }
178   ierr = PetscSFPackReclaim(sf,(PetscSFPack*)&link);CHKERRQ(ierr);
179   PetscFunctionReturn(0);
180 }
181 
182 static PetscErrorCode PetscSFReduceBegin_Allgatherv(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
183 {
184   PetscErrorCode         ierr;
185   PetscSFPack_Allgatherv link;
186   PetscSF_Allgatherv     *dat = (PetscSF_Allgatherv*)sf->data;
187   PetscInt               rstart;
188   PetscMPIInt            rank,count,recvcount;
189   MPI_Comm               comm;
190 
191   PetscFunctionBegin;
192   ierr = PetscSFPackGet_Allgatherv(sf,unit,rootdata,leafdata,&link);CHKERRQ(ierr);
193   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
194   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
195 
196   if (op == MPIU_REPLACE) {
197     /* REPLACE is only meaningful when all processes have the same leafdata to readue. Therefore copy from local leafdata is fine */
198     ierr = PetscLayoutGetRange(sf->map,&rstart,NULL);CHKERRQ(ierr);
199     ierr = PetscMemcpy(rootdata,(const char*)leafdata+(size_t)rstart*link->unitbytes,(size_t)sf->nroots*link->unitbytes);CHKERRQ(ierr);
200   } else {
201     /* Reduce all leafdata on rank 0, then scatter the result to root buffer, then reduce root buffer to leafdata */
202     if (!rank && !link->leafbuf) {ierr = PetscMalloc(sf->nleaves*link->unitbytes,&link->leafbuf);CHKERRQ(ierr);}
203     ierr = PetscMPIIntCast(sf->nleaves*link->bs,&count);CHKERRQ(ierr);
204     ierr = PetscMPIIntCast(sf->nroots,&recvcount);CHKERRQ(ierr);
205     ierr = MPI_Reduce(leafdata,link->leafbuf,count,link->basicunit,op,0,comm);CHKERRQ(ierr); /* Must do reduce with MPI builltin datatype basicunit */
206     if (!link->rootbuf) {ierr = PetscMalloc(sf->nroots*link->unitbytes,&link->rootbuf);CHKERRQ(ierr);}
207     ierr = MPIU_Iscatterv(link->leafbuf,dat->recvcounts,dat->displs,unit,link->rootbuf,recvcount,unit,0,comm,&link->request);CHKERRQ(ierr);
208   }
209   PetscFunctionReturn(0);
210 }
211 
212 PETSC_INTERN PetscErrorCode PetscSFReduceEnd_Allgatherv(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
213 {
214   PetscErrorCode         ierr,(*UnpackAndOp)(PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt,void*,const void*);
215   PetscSFPack_Allgatherv link;
216   PetscMPIInt            count;
217 
218   PetscFunctionBegin;
219   ierr = PetscSFPackGetInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,(PetscSFPack*)&link);CHKERRQ(ierr);
220   ierr = MPI_Wait(&link->request,MPI_STATUS_IGNORE);CHKERRQ(ierr);
221   if (op != MPIU_REPLACE) {
222     ierr = PetscSFPackGetUnpackAndOp(sf,(PetscSFPack)link,op,&UnpackAndOp);CHKERRQ(ierr);
223     if (UnpackAndOp) {ierr = (*UnpackAndOp)(sf->nroots,NULL,link->bs,NULL,rootdata,link->rootbuf);CHKERRQ(ierr);}
224     else if (sf->nroots) {
225       ierr = PetscMPIIntCast(sf->nroots,&count);CHKERRQ(ierr);
226       ierr = MPI_Reduce_local(link->rootbuf,rootdata,count,unit,op);CHKERRQ(ierr);
227     }
228   }
229   ierr = PetscSFPackReclaim(sf,(PetscSFPack*)&link);CHKERRQ(ierr);
230   PetscFunctionReturn(0);
231 }
232 
233 static PetscErrorCode PetscSFBcastToZero_Allgatherv(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
234 {
235   PetscErrorCode         ierr;
236   PetscSFPack_Allgatherv link;
237 
238   PetscFunctionBegin;
239   ierr = PetscSFBcastAndOpBegin_Gatherv(sf,unit,rootdata,leafdata,MPIU_REPLACE);CHKERRQ(ierr);
240   /* A simplified PetscSFBcastAndOpEnd_Allgatherv */
241   ierr = PetscSFPackGetInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,(PetscSFPack*)&link);CHKERRQ(ierr);
242   ierr = MPI_Wait(&link->request,MPI_STATUS_IGNORE);CHKERRQ(ierr);
243   ierr = PetscSFPackReclaim(sf,(PetscSFPack*)&link);CHKERRQ(ierr);
244   PetscFunctionReturn(0);
245 }
246 
247 /* This routine is very tricky (I believe it is rarely used with this kind of graph so just provide a simple but not-optimal implementation).
248 
249    Suppose we have three ranks. Rank 0 has a root with value 1. Rank 0,1,2 has a leaf with value 2,3,4 respectively. The leaves are connected
250    to the root on rank 0. Suppose op=MPI_SUM and rank 0,1,2 gets root state in their rank order. By definition of this routine, rank 0 sees 1
251    in root, fetches it into its leafupate, then updates root to 1 + 2 = 3; rank 1 sees 3 in root, fetches it into its leafupate, then updates
252    root to 3 + 3 = 6; rank 2 sees 6 in root, fetches it into its leafupdate, then updates root to 6 + 4 = 10.  At the end, leafupdate on rank
253    0,1,2 is 1,3,6 respectively. root is 10.
254 
255    One optimized implementation could be: starting from the initial state:
256              rank-0   rank-1    rank-2
257         Root     1
258         Leaf     2       3         4
259 
260    Shift leaves rightwards to leafupdate. Rank 0 gathers the root value and puts it in leafupdate. We have:
261              rank-0   rank-1    rank-2
262         Root     1
263         Leaf     2       3         4
264      Leafupdate  1       2         3
265 
266    Then, do MPI_Scan on leafupdate and get:
267              rank-0   rank-1    rank-2
268         Root     1
269         Leaf     2       3         4
270      Leafupdate  1       3         6
271 
272    Rank 2 sums its leaf and leafupdate, scatters the result to the root, and gets
273              rank-0   rank-1    rank-2
274         Root     10
275         Leaf     2       3         4
276      Leafupdate  1       3         6
277 
278    We use a simpler implementation. From the same initial state, we copy leafdata to leafupdate
279              rank-0   rank-1    rank-2
280         Root     1
281         Leaf     2       3         4
282      Leafupdate  2       3         4
283 
284    Do MPI_Exscan on leafupdate,
285              rank-0   rank-1    rank-2
286         Root     1
287         Leaf     2       3         4
288      Leafupdate  2       2         5
289 
290    BcastAndOp from root to leafupdate,
291              rank-0   rank-1    rank-2
292         Root     1
293         Leaf     2       3         4
294      Leafupdate  3       3         6
295 
296    Copy root to leafupdate on rank-0
297              rank-0   rank-1    rank-2
298         Root     1
299         Leaf     2       3         4
300      Leafupdate  1       3         6
301 
302    Reduce from leaf to root,
303              rank-0   rank-1    rank-2
304         Root     10
305         Leaf     2       3         4
306      Leafupdate  1       3         6
307 */
308 PETSC_INTERN PetscErrorCode PetscSFFetchAndOpBegin_Allgatherv(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
309 {
310   PetscErrorCode         ierr;
311   PetscSFPack_Allgatherv link;
312   MPI_Comm               comm;
313   PetscMPIInt            count;
314 
315   PetscFunctionBegin;
316   /* Copy leafdata to leafupdate */
317   ierr = PetscSFPackGet_Allgatherv(sf,unit,rootdata,leafdata,&link);CHKERRQ(ierr);
318   ierr = PetscMemcpy(leafupdate,leafdata,sf->nleaves*link->unitbytes);CHKERRQ(ierr);
319   ierr = PetscSFPackGetInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,(PetscSFPack*)&link);CHKERRQ(ierr);
320 
321   /* Exscan on leafupdate and then BcastAndOp rootdata to leafupdate */
322   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
323   ierr = PetscMPIIntCast(sf->nleaves,&count);CHKERRQ(ierr);
324   if (op == MPIU_REPLACE) {
325     PetscMPIInt size,rank,prev,next;
326     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
327     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
328     prev = rank ?            rank-1 : MPI_PROC_NULL;
329     next = (rank < size-1) ? rank+1 : MPI_PROC_NULL;
330     ierr = MPI_Sendrecv_replace(leafupdate,count,unit,next,link->tag,prev,link->tag,comm,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
331   } else {ierr = MPI_Exscan(MPI_IN_PLACE,leafupdate,count,unit,op,comm);CHKERRQ(ierr);}
332   ierr = PetscSFPackReclaim(sf,(PetscSFPack*)&link);CHKERRQ(ierr);
333   ierr = PetscSFBcastAndOpBegin(sf,unit,rootdata,leafupdate,op);CHKERRQ(ierr);
334   ierr = PetscSFBcastAndOpEnd(sf,unit,rootdata,leafupdate,op);CHKERRQ(ierr);
335 
336   /* Bcast roots to rank 0's leafupdate */
337   ierr = PetscSFBcastToZero_Private(sf,unit,rootdata,leafupdate);CHKERRQ(ierr); /* Using this line makes Allgather SFs able to inherit this routine */
338 
339   /* Reduce leafdata to rootdata */
340   ierr = PetscSFReduceBegin(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
341   PetscFunctionReturn(0);
342 }
343 
344 PETSC_INTERN PetscErrorCode PetscSFFetchAndOpEnd_Allgatherv(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
345 {
346   PetscErrorCode         ierr;
347 
348   PetscFunctionBegin;
349   ierr = PetscSFReduceEnd(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
350   PetscFunctionReturn(0);
351 }
352 
353 /* Get root ranks accessing my leaves */
354 PETSC_INTERN PetscErrorCode PetscSFGetRootRanks_Allgatherv(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
355 {
356   PetscErrorCode ierr;
357   PetscInt       i,j,k,size;
358   const PetscInt *range;
359 
360   PetscFunctionBegin;
361   /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */
362   if (sf->nranks && !sf->ranks) { /* On rank!=0, sf->nranks=0. The sf->nranks test makes this routine also works for sfgatherv */
363     size = sf->nranks;
364     ierr = PetscLayoutGetRanges(sf->map,&range);CHKERRQ(ierr);
365     ierr = PetscMalloc4(size,&sf->ranks,size+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr);
366     for (i=0; i<size; i++) sf->ranks[i] = i;
367     ierr = PetscArraycpy(sf->roffset,range,size+1);CHKERRQ(ierr);
368     for (i=0; i<sf->nleaves; i++) sf->rmine[i] = i; /*rmine are never NULL even for contiguous leaves */
369     for (i=0; i<size; i++) {
370       for (j=range[i],k=0; j<range[i+1]; j++,k++) sf->rremote[j] = k;
371     }
372   }
373 
374   if (nranks)  *nranks  = sf->nranks;
375   if (ranks)   *ranks   = sf->ranks;
376   if (roffset) *roffset = sf->roffset;
377   if (rmine)   *rmine   = sf->rmine;
378   if (rremote) *rremote = sf->rremote;
379   PetscFunctionReturn(0);
380 }
381 
382 /* Get leaf ranks accessing my roots */
383 PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Allgatherv(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
384 {
385   PetscErrorCode     ierr;
386   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv*)sf->data;
387   MPI_Comm           comm;
388   PetscMPIInt        size,rank;
389   PetscInt           i,j;
390 
391   PetscFunctionBegin;
392   /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */
393   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
394   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
395   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
396   if (niranks) *niranks = size;
397 
398   /* PetscSF_Basic has distinguished incoming ranks. Here we do not need that. But we must put self as the first and
399      sort other ranks. See comments in PetscSFSetUp_Basic about MatGetBrowsOfAoCols_MPIAIJ on why.
400    */
401   if (iranks) {
402     if (!dat->iranks) {
403       ierr = PetscMalloc1(size,&dat->iranks);CHKERRQ(ierr);
404       dat->iranks[0] = rank;
405       for (i=0,j=1; i<size; i++) {if (i == rank) continue; dat->iranks[j++] = i;}
406     }
407     *iranks = dat->iranks; /* dat->iranks was init'ed to NULL by PetscNewLog */
408   }
409 
410   if (ioffset) {
411     if (!dat->ioffset) {
412       ierr = PetscMalloc1(size+1,&dat->ioffset);CHKERRQ(ierr);
413       for (i=0; i<=size; i++) dat->ioffset[i] = i*sf->nroots;
414     }
415     *ioffset = dat->ioffset;
416   }
417 
418   if (irootloc) {
419     if (!dat->irootloc) {
420       ierr = PetscMalloc1(sf->nleaves,&dat->irootloc);CHKERRQ(ierr);
421       for (i=0; i<size; i++) {
422         for (j=0; j<sf->nroots; j++) dat->irootloc[i*sf->nroots+j] = j;
423       }
424     }
425     *irootloc = dat->irootloc;
426   }
427   PetscFunctionReturn(0);
428 }
429 
430 PETSC_INTERN PetscErrorCode PetscSFCreateLocalSF_Allgatherv(PetscSF sf,PetscSF *out)
431 {
432   PetscInt       i,nroots,nleaves,rstart,*ilocal;
433   PetscSFNode    *iremote;
434   PetscSF        lsf;
435   PetscErrorCode ierr;
436 
437   PetscFunctionBegin;
438   nroots  = sf->nroots;
439   nleaves = sf->nleaves ? sf->nroots : 0;
440   ierr    = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr);
441   ierr    = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr);
442   ierr    = PetscLayoutGetRange(sf->map,&rstart,NULL);CHKERRQ(ierr);
443 
444   for (i=0; i<nleaves; i++) {
445     ilocal[i]        = rstart + i; /* lsf does not change leave indices */
446     iremote[i].rank  = 0;          /* rank in PETSC_COMM_SELF */
447     iremote[i].index = i;          /* root index */
448   }
449 
450   ierr = PetscSFCreate(PETSC_COMM_SELF,&lsf);CHKERRQ(ierr);
451   ierr = PetscSFSetGraph(lsf,nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
452   ierr = PetscSFSetUp(lsf);CHKERRQ(ierr);
453   *out = lsf;
454   PetscFunctionReturn(0);
455 }
456 
457 PETSC_INTERN PetscErrorCode PetscSFCreate_Allgatherv(PetscSF sf)
458 {
459   PetscErrorCode     ierr;
460   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv*)sf->data;
461 
462   PetscFunctionBegin;
463   sf->ops->SetUp           = PetscSFSetUp_Allgatherv;
464   sf->ops->Reset           = PetscSFReset_Allgatherv;
465   sf->ops->Destroy         = PetscSFDestroy_Allgatherv;
466   sf->ops->GetRootRanks    = PetscSFGetRootRanks_Allgatherv;
467   sf->ops->GetLeafRanks    = PetscSFGetLeafRanks_Allgatherv;
468   sf->ops->GetGraph        = PetscSFGetGraph_Allgatherv;
469   sf->ops->BcastAndOpBegin = PetscSFBcastAndOpBegin_Allgatherv;
470   sf->ops->BcastAndOpEnd   = PetscSFBcastAndOpEnd_Allgatherv;
471   sf->ops->ReduceBegin     = PetscSFReduceBegin_Allgatherv;
472   sf->ops->ReduceEnd       = PetscSFReduceEnd_Allgatherv;
473   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Allgatherv;
474   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Allgatherv;
475   sf->ops->CreateLocalSF   = PetscSFCreateLocalSF_Allgatherv;
476   sf->ops->BcastToZero     = PetscSFBcastToZero_Allgatherv;
477 
478   ierr = PetscNewLog(sf,&dat);CHKERRQ(ierr);
479   sf->data = (void*)dat;
480   PetscFunctionReturn(0);
481 }
482