Actual source code: sfallgatherv.c

  1: #include <../src/vec/is/sf/impls/basic/allgatherv/sfallgatherv.h>

  3: /* PetscSFGetGraph is non-collective. An implementation should not have collective calls */
  4: PETSC_INTERN PetscErrorCode PetscSFGetGraph_Allgatherv(PetscSF sf, PetscInt *nroots, PetscInt *nleaves, const PetscInt **ilocal, const PetscSFNode **iremote)
  5: {
  6:   PetscInt        j, k;
  7:   const PetscInt *range;
  8:   PetscMPIInt     size;

 10:   PetscFunctionBegin;
 11:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf), &size));
 12:   if (nroots) *nroots = sf->nroots;
 13:   if (nleaves) *nleaves = sf->nleaves;
 14:   if (ilocal) *ilocal = NULL; /* Contiguous leaves */
 15:   if (iremote) {
 16:     if (!sf->remote && sf->nleaves) { /* The && sf->nleaves makes sfgatherv able to inherit this routine */
 17:       PetscCall(PetscLayoutGetRanges(sf->map, &range));
 18:       PetscCall(PetscMalloc1(sf->nleaves, &sf->remote));
 19:       sf->remote_alloc = sf->remote;
 20:       for (PetscMPIInt i = 0; i < size; i++) {
 21:         for (j = range[i], k = 0; j < range[i + 1]; j++, k++) {
 22:           sf->remote[j].rank  = i;
 23:           sf->remote[j].index = k;
 24:         }
 25:       }
 26:     }
 27:     *iremote = sf->remote;
 28:   }
 29:   PetscFunctionReturn(PETSC_SUCCESS);
 30: }

 32: PETSC_INTERN PetscErrorCode PetscSFSetUp_Allgatherv(PetscSF sf)
 33: {
 34:   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
 35:   PetscMPIInt         size;
 36:   PetscInt            i;
 37:   const PetscInt     *range;
 38:   MPI_Comm            comm;

 40:   PetscFunctionBegin;
 41:   PetscCall(PetscSFSetUp_Allgather(sf));
 42:   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
 43:   PetscCallMPI(MPI_Comm_size(comm, &size));
 44:   if (sf->nleaves) { /* This if (sf->nleaves) test makes sfgatherv able to inherit this routine */
 45:     PetscBool isallgatherv = PETSC_FALSE;

 47:     PetscCall(PetscMalloc1(size, &dat->recvcounts));
 48:     PetscCall(PetscMalloc1(size, &dat->displs));
 49:     PetscCall(PetscLayoutGetRanges(sf->map, &range));

 51:     for (i = 0; i < size; i++) {
 52:       PetscCall(PetscMPIIntCast(range[i], &dat->displs[i]));
 53:       PetscCall(PetscMPIIntCast(range[i + 1] - range[i], &dat->recvcounts[i]));
 54:     }

 56:     /* check if we actually have a one-to-all pattern */
 57:     PetscCall(PetscObjectTypeCompare((PetscObject)sf, PETSCSFALLGATHERV, &isallgatherv));
 58:     if (isallgatherv) {
 59:       PetscMPIInt rank, nRanksWithZeroRoots;

 61:       nRanksWithZeroRoots = (sf->nroots == 0) ? 1 : 0; /* I have no roots */
 62:       PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &nRanksWithZeroRoots, 1, MPI_INT, MPI_SUM, comm));
 63:       if (nRanksWithZeroRoots == size - 1) { /* Only one rank has roots, which indicates a bcast pattern */
 64:         dat->bcast_pattern = PETSC_TRUE;
 65:         PetscCallMPI(MPI_Comm_rank(comm, &rank));
 66:         dat->bcast_root = sf->nroots > 0 ? rank : -1;
 67:         PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &dat->bcast_root, 1, MPI_INT, MPI_MAX, comm));
 68:       }
 69:     }
 70:   }
 71:   PetscFunctionReturn(PETSC_SUCCESS);
 72: }

 74: PETSC_INTERN PetscErrorCode PetscSFReset_Allgatherv(PetscSF sf)
 75: {
 76:   PetscSF_Allgatherv *dat  = (PetscSF_Allgatherv *)sf->data;
 77:   PetscSFLink         link = dat->avail, next;

 79:   PetscFunctionBegin;
 80:   PetscCall(PetscFree(dat->iranks));
 81:   PetscCall(PetscFree(dat->ioffset));
 82:   PetscCall(PetscFree(dat->irootloc));
 83:   PetscCall(PetscFree(dat->recvcounts));
 84:   PetscCall(PetscFree(dat->displs));
 85:   PetscCheck(!dat->inuse, PetscObjectComm((PetscObject)sf), PETSC_ERR_ARG_WRONGSTATE, "Outstanding operation has not been completed");
 86:   for (; link; link = next) {
 87:     next = link->next;
 88:     PetscCall(PetscSFLinkDestroy(sf, link));
 89:   }
 90:   dat->avail = NULL;
 91:   PetscFunctionReturn(PETSC_SUCCESS);
 92: }

 94: PETSC_INTERN PetscErrorCode PetscSFDestroy_Allgatherv(PetscSF sf)
 95: {
 96:   PetscFunctionBegin;
 97:   PetscCall(PetscSFReset_Allgatherv(sf));
 98:   PetscCall(PetscFree(sf->data));
 99:   PetscFunctionReturn(PETSC_SUCCESS);
100: }

102: static PetscErrorCode PetscSFBcastBegin_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, const void *rootdata, PetscMemType leafmtype, void *leafdata, MPI_Op op)
103: {
104:   PetscSFLink         link;
105:   PetscMPIInt         sendcount, rank, nleaves;
106:   MPI_Comm            comm;
107:   void               *rootbuf = NULL, *leafbuf = NULL;
108:   MPI_Request        *req = NULL;
109:   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;

111:   PetscFunctionBegin;
112:   PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_BCAST, &link));
113:   PetscCall(PetscSFLinkPackRootData(sf, link, PETSCSF_REMOTE, rootdata));
114:   PetscCall(PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE /* device2host before sending */));
115:   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
116:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
117:   PetscCall(PetscMPIIntCast(sf->nroots, &sendcount));
118:   PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, PETSCSF_ROOT2LEAF, &rootbuf, &leafbuf, &req, NULL));

120:   if (dat->bcast_pattern && rank == dat->bcast_root) PetscCall((*link->Memcpy)(link, link->leafmtype_mpi, leafbuf, link->rootmtype_mpi, rootbuf, (size_t)sendcount * link->unitbytes));
121:   /* Ready the buffers for MPI */
122:   PetscCall(PetscSFLinkSyncStreamBeforeCallMPI(sf, link));
123:   PetscCall(PetscMPIIntCast(sf->nleaves, &nleaves));
124:   if (dat->bcast_pattern) PetscCallMPI(MPIU_Ibcast(leafbuf, nleaves, unit, dat->bcast_root, comm, req));
125:   else PetscCallMPI(MPIU_Iallgatherv(rootbuf, sendcount, unit, leafbuf, dat->recvcounts, dat->displs, unit, comm, req));
126:   PetscFunctionReturn(PETSC_SUCCESS);
127: }

129: static PetscErrorCode PetscSFReduceBegin_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType leafmtype, const void *leafdata, PetscMemType rootmtype, void *rootdata, MPI_Op op)
130: {
131:   PetscSFLink         link;
132:   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
133:   PetscInt            rstart;
134:   PetscMPIInt         rank, count, recvcount;
135:   MPI_Comm            comm;
136:   void               *rootbuf = NULL, *leafbuf = NULL;
137:   MPI_Request        *req = NULL;

139:   PetscFunctionBegin;
140:   PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_REDUCE, &link));
141:   if (op == MPI_REPLACE) {
142:     /* REPLACE is only meaningful when all processes have the same leafdata to reduce. Therefore copying from local leafdata is fine */
143:     PetscCall(PetscLayoutGetRange(sf->map, &rstart, NULL));
144:     PetscCall((*link->Memcpy)(link, rootmtype, rootdata, leafmtype, (const char *)leafdata + (size_t)rstart * link->unitbytes, (size_t)sf->nroots * link->unitbytes));
145:     if (PetscMemTypeDevice(leafmtype) && PetscMemTypeHost(rootmtype)) PetscCall((*link->SyncStream)(link));
146:   } else {
147:     PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
148:     PetscCall(PetscSFLinkPackLeafData(sf, link, PETSCSF_REMOTE, leafdata));
149:     PetscCall(PetscSFLinkCopyLeafBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE /* device2host before sending */));
150:     PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, PETSCSF_LEAF2ROOT, &rootbuf, &leafbuf, &req, NULL));
151:     PetscCall(PetscSFLinkSyncStreamBeforeCallMPI(sf, link));
152:     if (dat->bcast_pattern) {
153:       PetscMPIInt nleavesi;

155:       PetscCall(PetscMPIIntCast(sf->nleaves, &nleavesi));
156: #if defined(PETSC_HAVE_OPENMPI) /* Workaround: cuda-aware Open MPI 4.1.3 does not support MPI_Ireduce() with device buffers */
157:       *req = MPI_REQUEST_NULL;  /* Set NULL so that we can safely MPI_Wait(req) */
158:       PetscCallMPI(MPIU_Reduce(leafbuf, rootbuf, nleavesi, unit, op, dat->bcast_root, comm));
159: #else
160:       PetscCallMPI(MPIU_Ireduce(leafbuf, rootbuf, nleavesi, unit, op, dat->bcast_root, comm, req));
161: #endif
162:     } else { /* Reduce leafdata, then scatter to rootdata */
163:       PetscCallMPI(MPI_Comm_rank(comm, &rank));
164:       PetscCall(PetscMPIIntCast(dat->rootbuflen[PETSCSF_REMOTE], &recvcount));
165:       /* Allocate a separate leaf buffer on rank 0 */
166:       if (rank == 0 && !link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]) {
167:         PetscCall(PetscSFMalloc(sf, link->leafmtype_mpi, sf->leafbuflen[PETSCSF_REMOTE] * link->unitbytes, (void **)&link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]));
168:       }
169:       /* In case we already copied leafdata from device to host (i.e., no use_gpu_aware_mpi), we need to adjust leafbuf on rank 0 */
170:       if (rank == 0 && link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi] == leafbuf) leafbuf = MPI_IN_PLACE;
171:       PetscCall(PetscMPIIntCast(sf->nleaves * link->bs, &count));
172:       PetscCallMPI(MPI_Reduce(leafbuf, link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi], count, link->basicunit, op, 0, comm)); /* Must do reduce with MPI builtin datatype basicunit */
173:       PetscCallMPI(MPIU_Iscatterv(link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi], dat->recvcounts, dat->displs, unit, rootbuf, recvcount, unit, 0, comm, req));
174:     }
175:   }
176:   PetscFunctionReturn(PETSC_SUCCESS);
177: }

179: PETSC_INTERN PetscErrorCode PetscSFReduceEnd_Allgatherv(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *rootdata, MPI_Op op)
180: {
181:   PetscSFLink link;

183:   PetscFunctionBegin;
184:   if (op == MPI_REPLACE) {
185:     /* A rare case happens when op is MPI_REPLACE, using GPUs but no GPU aware MPI. In PetscSFReduceBegin_Allgather(v),
186:       we did a device to device copy and in effect finished the communication. But in PetscSFLinkFinishCommunication()
187:       of PetscSFReduceEnd_Basic(), it thinks since there is rootbuf, it calls PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI().
188:       It does a host to device memory copy on rootbuf, wrongly overwriting the results. So we don't overload
189:       PetscSFReduceEnd_Basic() in this case, and just reclaim the link.
190:      */
191:     PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
192:     PetscCall(PetscSFLinkReclaim(sf, &link));
193:   } else {
194:     PetscCall(PetscSFReduceEnd_Basic(sf, unit, leafdata, rootdata, op));
195:   }
196:   PetscFunctionReturn(PETSC_SUCCESS);
197: }

199: static PetscErrorCode PetscSFBcastToZero_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, const void *rootdata, PetscMemType leafmtype, void *leafdata)
200: {
201:   PetscSFLink         link;
202:   PetscMPIInt         rank;
203:   PetscMPIInt         sendcount;
204:   MPI_Comm            comm;
205:   PetscSF_Allgatherv *dat     = (PetscSF_Allgatherv *)sf->data;
206:   void               *rootbuf = NULL, *leafbuf = NULL; /* buffer seen by MPI */
207:   MPI_Request        *req = NULL;

209:   PetscFunctionBegin;
210:   PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, MPI_REPLACE, PETSCSF_BCAST, &link));
211:   PetscCall(PetscSFLinkPackRootData(sf, link, PETSCSF_REMOTE, rootdata));
212:   PetscCall(PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE /* device2host before sending */));
213:   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
214:   PetscCall(PetscMPIIntCast(sf->nroots, &sendcount));
215:   PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, PETSCSF_ROOT2LEAF, &rootbuf, &leafbuf, &req, NULL));
216:   PetscCall(PetscSFLinkSyncStreamBeforeCallMPI(sf, link));
217:   PetscCallMPI(MPIU_Igatherv(rootbuf, sendcount, unit, leafbuf, dat->recvcounts, dat->displs, unit, 0 /*rank 0*/, comm, req));

219:   PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
220:   PetscCall(PetscSFLinkFinishCommunication(sf, link, PETSCSF_ROOT2LEAF));
221:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
222:   if (rank == 0 && PetscMemTypeDevice(leafmtype) && !sf->use_gpu_aware_mpi) PetscCall((*link->Memcpy)(link, PETSC_MEMTYPE_DEVICE, leafdata, PETSC_MEMTYPE_HOST, link->leafbuf[PETSC_MEMTYPE_HOST], sf->leafbuflen[PETSCSF_REMOTE] * link->unitbytes));
223:   PetscCall(PetscSFLinkReclaim(sf, &link));
224:   PetscFunctionReturn(PETSC_SUCCESS);
225: }

227: /* 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).

229:    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
230:    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
231:    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
232:    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
233:    0,1,2 is 1,3,6 respectively. root is 10.

235:    We use a simpler implementation. From the same initial state, we copy leafdata to leafupdate
236:              rank-0   rank-1    rank-2
237:         Root     1
238:         Leaf     2       3         4
239:      Leafupdate  2       3         4

241:    Do MPI_Exscan on leafupdate,
242:              rank-0   rank-1    rank-2
243:         Root     1
244:         Leaf     2       3         4
245:      Leafupdate  2       2         5

247:    BcastAndOp from root to leafupdate,
248:              rank-0   rank-1    rank-2
249:         Root     1
250:         Leaf     2       3         4
251:      Leafupdate  3       3         6

253:    Copy root to leafupdate on rank-0
254:              rank-0   rank-1    rank-2
255:         Root     1
256:         Leaf     2       3         4
257:      Leafupdate  1       3         6

259:    Reduce from leaf to root,
260:              rank-0   rank-1    rank-2
261:         Root     10
262:         Leaf     2       3         4
263:      Leafupdate  1       3         6
264: */
265: PETSC_INTERN PetscErrorCode PetscSFFetchAndOpBegin_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, void *rootdata, PetscMemType leafmtype, const void *leafdata, void *leafupdate, MPI_Op op)
266: {
267:   PetscSFLink link;
268:   MPI_Comm    comm;
269:   PetscMPIInt count;

271:   PetscFunctionBegin;
272:   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
273:   PetscCheck(!PetscMemTypeDevice(rootmtype) && !PetscMemTypeDevice(leafmtype), comm, PETSC_ERR_SUP, "Do FetchAndOp on device");
274:   /* Copy leafdata to leafupdate */
275:   PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_FETCH, &link));
276:   PetscCall(PetscSFLinkPackLeafData(sf, link, PETSCSF_REMOTE, leafdata)); /* Sync the device */
277:   PetscCall((*link->Memcpy)(link, leafmtype, leafupdate, leafmtype, leafdata, sf->nleaves * link->unitbytes));
278:   PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));

280:   /* Exscan on leafupdate and then BcastAndOp rootdata to leafupdate */
281:   if (op == MPI_REPLACE) {
282:     PetscMPIInt size, rank, prev, next;
283:     PetscCallMPI(MPI_Comm_rank(comm, &rank));
284:     PetscCallMPI(MPI_Comm_size(comm, &size));
285:     prev = rank ? rank - 1 : MPI_PROC_NULL;
286:     next = (rank < size - 1) ? rank + 1 : MPI_PROC_NULL;
287:     PetscCall(PetscMPIIntCast(sf->nleaves, &count));
288:     PetscCallMPI(MPI_Sendrecv_replace(leafupdate, count, unit, next, link->tag, prev, link->tag, comm, MPI_STATUSES_IGNORE));
289:   } else {
290:     PetscCall(PetscMPIIntCast(sf->nleaves * link->bs, &count));
291:     PetscCallMPI(MPI_Exscan(MPI_IN_PLACE, leafupdate, count, link->basicunit, op, comm));
292:   }
293:   PetscCall(PetscSFLinkReclaim(sf, &link));
294:   PetscCall(PetscSFBcastBegin(sf, unit, rootdata, leafupdate, op));
295:   PetscCall(PetscSFBcastEnd(sf, unit, rootdata, leafupdate, op));

297:   /* Bcast roots to rank 0's leafupdate */
298:   PetscCall(PetscSFBcastToZero_Private(sf, unit, rootdata, leafupdate)); /* Using this line makes Allgather SFs able to inherit this routine */

300:   /* Reduce leafdata to rootdata */
301:   PetscCall(PetscSFReduceBegin(sf, unit, leafdata, rootdata, op));
302:   PetscFunctionReturn(PETSC_SUCCESS);
303: }

305: PETSC_INTERN PetscErrorCode PetscSFFetchAndOpEnd_Allgatherv(PetscSF sf, MPI_Datatype unit, void *rootdata, const void *leafdata, void *leafupdate, MPI_Op op)
306: {
307:   PetscFunctionBegin;
308:   PetscCall(PetscSFReduceEnd(sf, unit, leafdata, rootdata, op));
309:   PetscFunctionReturn(PETSC_SUCCESS);
310: }

312: /* Get root ranks accessing my leaves */
313: PETSC_INTERN PetscErrorCode PetscSFGetRootRanks_Allgatherv(PetscSF sf, PetscMPIInt *nranks, const PetscMPIInt **ranks, const PetscInt **roffset, const PetscInt **rmine, const PetscInt **rremote)
314: {
315:   PetscInt        j, k, size;
316:   const PetscInt *range;

318:   PetscFunctionBegin;
319:   /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */
320:   if (sf->nranks && !sf->ranks) { /* On rank!=0, sf->nranks=0. The sf->nranks test makes this routine also works for sfgatherv */
321:     size = sf->nranks;
322:     PetscCall(PetscLayoutGetRanges(sf->map, &range));
323:     PetscCall(PetscMalloc4(size, &sf->ranks, size + 1, &sf->roffset, sf->nleaves, &sf->rmine, sf->nleaves, &sf->rremote));
324:     for (PetscMPIInt i = 0; i < size; i++) sf->ranks[i] = i;
325:     PetscCall(PetscArraycpy(sf->roffset, range, size + 1));
326:     for (PetscInt i = 0; i < sf->nleaves; i++) sf->rmine[i] = i; /*rmine are never NULL even for contiguous leaves */
327:     for (PetscMPIInt i = 0; i < size; i++) {
328:       for (j = range[i], k = 0; j < range[i + 1]; j++, k++) sf->rremote[j] = k;
329:     }
330:   }

332:   if (nranks) *nranks = sf->nranks;
333:   if (ranks) *ranks = sf->ranks;
334:   if (roffset) *roffset = sf->roffset;
335:   if (rmine) *rmine = sf->rmine;
336:   if (rremote) *rremote = sf->rremote;
337:   PetscFunctionReturn(PETSC_SUCCESS);
338: }

340: /* Get leaf ranks accessing my roots */
341: PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Allgatherv(PetscSF sf, PetscMPIInt *niranks, const PetscMPIInt **iranks, const PetscInt **ioffset, const PetscInt **irootloc)
342: {
343:   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
344:   MPI_Comm            comm;
345:   PetscMPIInt         size, rank;

347:   PetscFunctionBegin;
348:   /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */
349:   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
350:   PetscCallMPI(MPI_Comm_size(comm, &size));
351:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
352:   if (niranks) *niranks = size;

354:   /* PetscSF_Basic has distinguished incoming ranks. Here we do not need that. But we must put self as the first and
355:      sort other ranks. See comments in PetscSFSetUp_Basic about MatGetBrowsOfAoCols_MPIAIJ on why.
356:    */
357:   if (iranks) {
358:     if (!dat->iranks) {
359:       PetscCall(PetscMalloc1(size, &dat->iranks));
360:       dat->iranks[0] = rank;
361:       for (PetscMPIInt i = 0, j = 1; i < size; i++) {
362:         if (i == rank) continue;
363:         dat->iranks[j++] = i;
364:       }
365:     }
366:     *iranks = dat->iranks; /* dat->iranks was init'ed to NULL by PetscNew */
367:   }

369:   if (ioffset) {
370:     if (!dat->ioffset) {
371:       PetscCall(PetscMalloc1(size + 1, &dat->ioffset));
372:       for (PetscMPIInt i = 0; i <= size; i++) dat->ioffset[i] = i * sf->nroots;
373:     }
374:     *ioffset = dat->ioffset;
375:   }

377:   if (irootloc) {
378:     if (!dat->irootloc) {
379:       PetscCall(PetscMalloc1(sf->nleaves, &dat->irootloc));
380:       for (PetscMPIInt i = 0; i < size; i++) {
381:         for (PetscInt j = 0; j < sf->nroots; j++) dat->irootloc[i * sf->nroots + j] = j;
382:       }
383:     }
384:     *irootloc = dat->irootloc;
385:   }
386:   PetscFunctionReturn(PETSC_SUCCESS);
387: }

389: PETSC_INTERN PetscErrorCode PetscSFCreateLocalSF_Allgatherv(PetscSF sf, PetscSF *out)
390: {
391:   PetscInt     i, nroots, nleaves, rstart, *ilocal;
392:   PetscSFNode *iremote;
393:   PetscSF      lsf;

395:   PetscFunctionBegin;
396:   nleaves = sf->nleaves ? sf->nroots : 0; /* sf->nleaves can be zero with SFGather(v) */
397:   nroots  = nleaves;
398:   PetscCall(PetscMalloc1(nleaves, &ilocal));
399:   PetscCall(PetscMalloc1(nleaves, &iremote));
400:   PetscCall(PetscLayoutGetRange(sf->map, &rstart, NULL));

402:   for (i = 0; i < nleaves; i++) {
403:     ilocal[i]        = rstart + i; /* lsf does not change leave indices */
404:     iremote[i].rank  = 0;          /* rank in PETSC_COMM_SELF */
405:     iremote[i].index = i;          /* root index */
406:   }

408:   PetscCall(PetscSFCreate(PETSC_COMM_SELF, &lsf));
409:   PetscCall(PetscSFSetGraph(lsf, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
410:   PetscCall(PetscSFSetUp(lsf));
411:   *out = lsf;
412:   PetscFunctionReturn(PETSC_SUCCESS);
413: }

415: PETSC_INTERN PetscErrorCode PetscSFCreate_Allgatherv(PetscSF sf)
416: {
417:   PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;

419:   PetscFunctionBegin;
420:   sf->ops->BcastEnd  = PetscSFBcastEnd_Basic;
421:   sf->ops->ReduceEnd = PetscSFReduceEnd_Allgatherv;

423:   sf->ops->SetUp           = PetscSFSetUp_Allgatherv;
424:   sf->ops->Reset           = PetscSFReset_Allgatherv;
425:   sf->ops->Destroy         = PetscSFDestroy_Allgatherv;
426:   sf->ops->GetRootRanks    = PetscSFGetRootRanks_Allgatherv;
427:   sf->ops->GetLeafRanks    = PetscSFGetLeafRanks_Allgatherv;
428:   sf->ops->GetGraph        = PetscSFGetGraph_Allgatherv;
429:   sf->ops->BcastBegin      = PetscSFBcastBegin_Allgatherv;
430:   sf->ops->ReduceBegin     = PetscSFReduceBegin_Allgatherv;
431:   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Allgatherv;
432:   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Allgatherv;
433:   sf->ops->CreateLocalSF   = PetscSFCreateLocalSF_Allgatherv;
434:   sf->ops->BcastToZero     = PetscSFBcastToZero_Allgatherv;

436:   sf->collective = PETSC_TRUE;

438:   PetscCall(PetscNew(&dat));
439:   dat->bcast_root = -1;
440:   sf->data        = (void *)dat;
441:   PetscFunctionReturn(PETSC_SUCCESS);
442: }