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: const PetscInt *range;
37: MPI_Comm comm;
39: PetscFunctionBegin;
40: PetscCall(PetscSFSetUp_Allgather(sf));
41: PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
42: PetscCallMPI(MPI_Comm_size(comm, &size));
43: if (sf->nleaves) { /* This if (sf->nleaves) test makes sfgatherv able to inherit this routine */
44: PetscBool isallgatherv = PETSC_FALSE;
46: PetscCall(PetscMalloc1(size, &dat->recvcounts));
47: PetscCall(PetscMalloc1(size, &dat->displs));
48: PetscCall(PetscLayoutGetRanges(sf->map, &range));
50: for (PetscInt i = 0; i < size; i++) {
51: PetscCall(PetscMPIIntCast(range[i], &dat->displs[i]));
52: PetscCall(PetscMPIIntCast(range[i + 1] - range[i], &dat->recvcounts[i]));
53: }
55: /* check if we actually have a one-to-all pattern */
56: PetscCall(PetscObjectTypeCompare((PetscObject)sf, PETSCSFALLGATHERV, &isallgatherv));
57: if (isallgatherv) {
58: PetscMPIInt rank, nRanksWithZeroRoots;
60: nRanksWithZeroRoots = (sf->nroots == 0) ? 1 : 0; /* I have no roots */
61: PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &nRanksWithZeroRoots, 1, MPI_INT, MPI_SUM, comm));
62: if (nRanksWithZeroRoots == size - 1) { /* Only one rank has roots, which indicates a bcast pattern */
63: dat->bcast_pattern = PETSC_TRUE;
64: PetscCallMPI(MPI_Comm_rank(comm, &rank));
65: dat->bcast_root = sf->nroots > 0 ? rank : -1;
66: PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &dat->bcast_root, 1, MPI_INT, MPI_MAX, comm));
67: }
68: }
69: }
70: PetscFunctionReturn(PETSC_SUCCESS);
71: }
73: PETSC_INTERN PetscErrorCode PetscSFReset_Allgatherv(PetscSF sf)
74: {
75: PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
76: PetscSFLink link = dat->avail, next;
78: PetscFunctionBegin;
79: PetscCall(PetscFree(dat->iranks));
80: PetscCall(PetscFree(dat->ioffset));
81: PetscCall(PetscFree(dat->irootloc));
82: PetscCall(PetscFree(dat->recvcounts));
83: PetscCall(PetscFree(dat->displs));
84: PetscCheck(!dat->inuse, PetscObjectComm((PetscObject)sf), PETSC_ERR_ARG_WRONGSTATE, "Outstanding operation has not been completed");
85: for (; link; link = next) {
86: next = link->next;
87: PetscCall(PetscSFLinkDestroy(sf, link));
88: }
89: dat->avail = NULL;
90: PetscFunctionReturn(PETSC_SUCCESS);
91: }
93: PETSC_INTERN PetscErrorCode PetscSFDestroy_Allgatherv(PetscSF sf)
94: {
95: PetscFunctionBegin;
96: PetscCall(PetscSFReset_Allgatherv(sf));
97: PetscCall(PetscFree(sf->data));
98: PetscFunctionReturn(PETSC_SUCCESS);
99: }
101: static PetscErrorCode PetscSFBcastBegin_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, const void *rootdata, PetscMemType leafmtype, void *leafdata, MPI_Op op)
102: {
103: PetscSFLink link;
104: PetscMPIInt sendcount, rank, nleaves;
105: MPI_Comm comm;
106: void *rootbuf = NULL, *leafbuf = NULL;
107: MPI_Request *req = NULL;
108: PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
110: PetscFunctionBegin;
111: PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_BCAST, &link));
112: PetscCall(PetscSFLinkPackRootData(sf, link, PETSCSF_REMOTE, rootdata));
113: PetscCall(PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE /* device2host before sending */));
114: PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
115: PetscCallMPI(MPI_Comm_rank(comm, &rank));
116: PetscCall(PetscMPIIntCast(sf->nroots, &sendcount));
117: PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, PETSCSF_ROOT2LEAF, &rootbuf, &leafbuf, &req, NULL));
119: 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));
120: /* Ready the buffers for MPI */
121: PetscCall(PetscSFLinkSyncStreamBeforeCallMPI(sf, link));
122: PetscCall(PetscMPIIntCast(sf->nleaves, &nleaves));
123: if (dat->bcast_pattern) PetscCallMPI(MPIU_Ibcast(leafbuf, nleaves, unit, dat->bcast_root, comm, req));
124: else PetscCallMPI(MPIU_Iallgatherv(rootbuf, sendcount, unit, leafbuf, dat->recvcounts, dat->displs, unit, comm, req));
125: PetscFunctionReturn(PETSC_SUCCESS);
126: }
128: static PetscErrorCode PetscSFReduceBegin_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType leafmtype, const void *leafdata, PetscMemType rootmtype, void *rootdata, MPI_Op op)
129: {
130: PetscSFLink link;
131: PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
132: PetscInt rstart;
133: PetscMPIInt rank, count, recvcount;
134: MPI_Comm comm;
135: void *rootbuf = NULL, *leafbuf = NULL;
136: MPI_Request *req = NULL;
138: PetscFunctionBegin;
139: PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_REDUCE, &link));
140: if (op == MPI_REPLACE) {
141: /* REPLACE is only meaningful when all processes have the same leafdata to reduce. Therefore copying from local leafdata is fine */
142: PetscCall(PetscLayoutGetRange(sf->map, &rstart, NULL));
143: PetscCall((*link->Memcpy)(link, rootmtype, rootdata, leafmtype, (const char *)leafdata + (size_t)rstart * link->unitbytes, (size_t)sf->nroots * link->unitbytes));
144: if (PetscMemTypeDevice(leafmtype) && PetscMemTypeHost(rootmtype)) PetscCall((*link->SyncStream)(link));
145: } else {
146: PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
147: PetscCall(PetscSFLinkPackLeafData(sf, link, PETSCSF_REMOTE, leafdata));
148: PetscCall(PetscSFLinkCopyLeafBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE /* device2host before sending */));
149: PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, PETSCSF_LEAF2ROOT, &rootbuf, &leafbuf, &req, NULL));
150: PetscCall(PetscSFLinkSyncStreamBeforeCallMPI(sf, link));
151: if (dat->bcast_pattern) {
152: PetscInt nleaves = sf->nleaves;
153: PetscInt nreal;
154: PetscMPIInt nleavesi;
155: MPI_Datatype baseunit = unit;
157: PetscCall(MPIPetsc_Type_compare_contig(unit, MPIU_REAL, &nreal));
158: if (nreal > 0) {
159: baseunit = MPIU_REAL;
160: nleaves *= nreal;
161: #if PetscDefined(HAVE_COMPLEX)
162: } else {
163: PetscInt ncomplex;
165: PetscCall(MPIPetsc_Type_compare_contig(unit, MPIU_COMPLEX, &ncomplex));
166: if (ncomplex > 0) {
167: baseunit = MPIU_COMPLEX;
168: nleaves *= ncomplex;
169: }
170: #endif
171: }
172: PetscCall(PetscMPIIntCast(nleaves, &nleavesi));
173: #if defined(PETSC_HAVE_OPENMPI) /* Workaround: cuda-aware Open MPI 4.1.3 does not support MPI_Ireduce() with device buffers */
174: *req = MPI_REQUEST_NULL; /* Set NULL so that we can safely MPI_Wait(req) */
175: PetscCallMPI(MPIU_Reduce(leafbuf, rootbuf, nleavesi, baseunit, op, dat->bcast_root, comm));
176: #else
177: PetscCallMPI(MPIU_Ireduce(leafbuf, rootbuf, nleavesi, baseunit, op, dat->bcast_root, comm, req));
178: #endif
179: } else { /* Reduce leafdata, then scatter to rootdata */
180: PetscCallMPI(MPI_Comm_rank(comm, &rank));
181: PetscCall(PetscMPIIntCast(dat->rootbuflen[PETSCSF_REMOTE], &recvcount));
182: /* Allocate a separate leaf buffer on rank 0 */
183: if (rank == 0 && !link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]) PetscCall(PetscSFMalloc(sf, link->leafmtype_mpi, sf->leafbuflen[PETSCSF_REMOTE] * link->unitbytes, (void **)&link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi]));
184: /* 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 */
185: if (rank == 0 && link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi] == leafbuf) leafbuf = MPI_IN_PLACE;
186: PetscCall(PetscMPIIntCast(sf->nleaves * link->bs, &count));
187: 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 */
188: PetscCallMPI(MPIU_Iscatterv(link->leafbuf_alloc[PETSCSF_REMOTE][link->leafmtype_mpi], dat->recvcounts, dat->displs, unit, rootbuf, recvcount, unit, 0, comm, req));
189: }
190: }
191: PetscFunctionReturn(PETSC_SUCCESS);
192: }
194: PETSC_INTERN PetscErrorCode PetscSFReduceEnd_Allgatherv(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *rootdata, MPI_Op op)
195: {
196: PetscSFLink link;
198: PetscFunctionBegin;
199: if (op == MPI_REPLACE) {
200: /* A rare case happens when op is MPI_REPLACE, using GPUs but no GPU aware MPI. In PetscSFReduceBegin_Allgather(v),
201: we did a device to device copy and in effect finished the communication. But in PetscSFLinkFinishCommunication()
202: of PetscSFReduceEnd_Basic(), it thinks since there is rootbuf, it calls PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI().
203: It does a host to device memory copy on rootbuf, wrongly overwriting the results. So we don't overload
204: PetscSFReduceEnd_Basic() in this case, and just reclaim the link.
205: */
206: PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
207: PetscCall(PetscSFLinkReclaim(sf, &link));
208: } else {
209: PetscCall(PetscSFReduceEnd_Basic(sf, unit, leafdata, rootdata, op));
210: }
211: PetscFunctionReturn(PETSC_SUCCESS);
212: }
214: static PetscErrorCode PetscSFBcastToZero_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, const void *rootdata, PetscMemType leafmtype, void *leafdata)
215: {
216: PetscSFLink link;
217: PetscMPIInt rank;
218: PetscMPIInt sendcount;
219: MPI_Comm comm;
220: PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
221: void *rootbuf = NULL, *leafbuf = NULL; /* buffer seen by MPI */
222: MPI_Request *req = NULL;
224: PetscFunctionBegin;
225: PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, MPI_REPLACE, PETSCSF_BCAST, &link));
226: PetscCall(PetscSFLinkPackRootData(sf, link, PETSCSF_REMOTE, rootdata));
227: PetscCall(PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE /* device2host before sending */));
228: PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
229: PetscCall(PetscMPIIntCast(sf->nroots, &sendcount));
230: PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, PETSCSF_ROOT2LEAF, &rootbuf, &leafbuf, &req, NULL));
231: PetscCall(PetscSFLinkSyncStreamBeforeCallMPI(sf, link));
232: PetscCallMPI(MPIU_Igatherv(rootbuf, sendcount, unit, leafbuf, dat->recvcounts, dat->displs, unit, 0 /*rank 0*/, comm, req));
234: PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
235: PetscCall(PetscSFLinkFinishCommunication(sf, link, PETSCSF_ROOT2LEAF));
236: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
237: 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));
238: PetscCall(PetscSFLinkReclaim(sf, &link));
239: PetscFunctionReturn(PETSC_SUCCESS);
240: }
242: /* 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).
244: 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
245: 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
246: 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
247: 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
248: 0,1,2 is 1,3,6 respectively. root is 10.
250: We use a simpler implementation. From the same initial state, we copy leafdata to leafupdate
251: rank-0 rank-1 rank-2
252: Root 1
253: Leaf 2 3 4
254: Leafupdate 2 3 4
256: Do MPI_Exscan on leafupdate,
257: rank-0 rank-1 rank-2
258: Root 1
259: Leaf 2 3 4
260: Leafupdate 2 2 5
262: BcastAndOp from root to leafupdate,
263: rank-0 rank-1 rank-2
264: Root 1
265: Leaf 2 3 4
266: Leafupdate 3 3 6
268: Copy root to leafupdate on rank-0
269: rank-0 rank-1 rank-2
270: Root 1
271: Leaf 2 3 4
272: Leafupdate 1 3 6
274: Reduce from leaf to root,
275: rank-0 rank-1 rank-2
276: Root 10
277: Leaf 2 3 4
278: Leafupdate 1 3 6
279: */
280: PETSC_INTERN PetscErrorCode PetscSFFetchAndOpBegin_Allgatherv(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, void *rootdata, PetscMemType leafmtype, const void *leafdata, void *leafupdate, MPI_Op op)
281: {
282: PetscSFLink link;
283: MPI_Comm comm;
284: PetscMPIInt count;
286: PetscFunctionBegin;
287: PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
288: PetscCheck(!PetscMemTypeDevice(rootmtype) && !PetscMemTypeDevice(leafmtype), comm, PETSC_ERR_SUP, "Do FetchAndOp on device");
289: /* Copy leafdata to leafupdate */
290: PetscCall(PetscSFLinkCreate(sf, unit, rootmtype, rootdata, leafmtype, leafdata, op, PETSCSF_FETCH, &link));
291: PetscCall(PetscSFLinkPackLeafData(sf, link, PETSCSF_REMOTE, leafdata)); /* Sync the device */
292: PetscCall((*link->Memcpy)(link, leafmtype, leafupdate, leafmtype, leafdata, sf->nleaves * link->unitbytes));
293: PetscCall(PetscSFLinkGetInUse(sf, unit, rootdata, leafdata, PETSC_OWN_POINTER, &link));
295: /* Exscan on leafupdate and then BcastAndOp rootdata to leafupdate */
296: if (op == MPI_REPLACE) {
297: PetscMPIInt size, rank, prev, next;
298: PetscCallMPI(MPI_Comm_rank(comm, &rank));
299: PetscCallMPI(MPI_Comm_size(comm, &size));
300: prev = rank ? rank - 1 : MPI_PROC_NULL;
301: next = (rank < size - 1) ? rank + 1 : MPI_PROC_NULL;
302: PetscCall(PetscMPIIntCast(sf->nleaves, &count));
303: PetscCallMPI(MPI_Sendrecv_replace(leafupdate, count, unit, next, link->tag, prev, link->tag, comm, MPI_STATUSES_IGNORE));
304: } else {
305: PetscCall(PetscMPIIntCast(sf->nleaves * link->bs, &count));
306: PetscCallMPI(MPI_Exscan(MPI_IN_PLACE, leafupdate, count, link->basicunit, op, comm));
307: }
308: PetscCall(PetscSFLinkReclaim(sf, &link));
309: PetscCall(PetscSFBcastBegin(sf, unit, rootdata, leafupdate, op));
310: PetscCall(PetscSFBcastEnd(sf, unit, rootdata, leafupdate, op));
312: /* Bcast roots to rank 0's leafupdate */
313: PetscCall(PetscSFBcastToZero_Private(sf, unit, rootdata, leafupdate)); /* Using this line makes Allgather SFs able to inherit this routine */
315: /* Reduce leafdata to rootdata */
316: PetscCall(PetscSFReduceBegin(sf, unit, leafdata, rootdata, op));
317: PetscFunctionReturn(PETSC_SUCCESS);
318: }
320: PETSC_INTERN PetscErrorCode PetscSFFetchAndOpEnd_Allgatherv(PetscSF sf, MPI_Datatype unit, void *rootdata, const void *leafdata, void *leafupdate, MPI_Op op)
321: {
322: PetscFunctionBegin;
323: PetscCall(PetscSFReduceEnd(sf, unit, leafdata, rootdata, op));
324: PetscFunctionReturn(PETSC_SUCCESS);
325: }
327: /* Get root ranks accessing my leaves */
328: PETSC_INTERN PetscErrorCode PetscSFGetRootRanks_Allgatherv(PetscSF sf, PetscMPIInt *nranks, const PetscMPIInt **ranks, const PetscInt **roffset, const PetscInt **rmine, const PetscInt **rremote)
329: {
330: PetscInt j, k, size;
331: const PetscInt *range;
333: PetscFunctionBegin;
334: /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */
335: if (sf->nranks && !sf->ranks) { /* On rank!=0, sf->nranks=0. The sf->nranks test makes this routine also works for sfgatherv */
336: size = sf->nranks;
337: PetscCall(PetscLayoutGetRanges(sf->map, &range));
338: PetscCall(PetscMalloc4(size, &sf->ranks, size + 1, &sf->roffset, sf->nleaves, &sf->rmine, sf->nleaves, &sf->rremote));
339: for (PetscMPIInt i = 0; i < size; i++) sf->ranks[i] = i;
340: PetscCall(PetscArraycpy(sf->roffset, range, size + 1));
341: for (PetscInt i = 0; i < sf->nleaves; i++) sf->rmine[i] = i; /*rmine are never NULL even for contiguous leaves */
342: for (PetscMPIInt i = 0; i < size; i++) {
343: for (j = range[i], k = 0; j < range[i + 1]; j++, k++) sf->rremote[j] = k;
344: }
345: }
347: if (nranks) *nranks = sf->nranks;
348: if (ranks) *ranks = sf->ranks;
349: if (roffset) *roffset = sf->roffset;
350: if (rmine) *rmine = sf->rmine;
351: if (rremote) *rremote = sf->rremote;
352: PetscFunctionReturn(PETSC_SUCCESS);
353: }
355: /* Get leaf ranks accessing my roots */
356: PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Allgatherv(PetscSF sf, PetscMPIInt *niranks, const PetscMPIInt **iranks, const PetscInt **ioffset, const PetscInt **irootloc)
357: {
358: PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
359: MPI_Comm comm;
360: PetscMPIInt size, rank;
362: PetscFunctionBegin;
363: /* Lazily construct these large arrays if users really need them for this type of SF. Very likely, they do not */
364: PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
365: PetscCallMPI(MPI_Comm_size(comm, &size));
366: PetscCallMPI(MPI_Comm_rank(comm, &rank));
367: if (niranks) *niranks = size;
369: /* PetscSF_Basic has distinguished incoming ranks. Here we do not need that. But we must put self as the first and
370: sort other ranks. See comments in PetscSFSetUp_Basic about MatGetBrowsOfAoCols_MPIAIJ on why.
371: */
372: if (iranks) {
373: if (!dat->iranks) {
374: PetscCall(PetscMalloc1(size, &dat->iranks));
375: dat->iranks[0] = rank;
376: for (PetscMPIInt i = 0, j = 1; i < size; i++) {
377: if (i == rank) continue;
378: dat->iranks[j++] = i;
379: }
380: }
381: *iranks = dat->iranks; /* dat->iranks was init'ed to NULL by PetscNew */
382: }
384: if (ioffset) {
385: if (!dat->ioffset) {
386: PetscCall(PetscMalloc1(size + 1, &dat->ioffset));
387: for (PetscMPIInt i = 0; i <= size; i++) dat->ioffset[i] = i * sf->nroots;
388: }
389: *ioffset = dat->ioffset;
390: }
392: if (irootloc) {
393: if (!dat->irootloc) {
394: PetscCall(PetscMalloc1(sf->nleaves, &dat->irootloc));
395: for (PetscMPIInt i = 0; i < size; i++) {
396: for (PetscInt j = 0; j < sf->nroots; j++) dat->irootloc[i * sf->nroots + j] = j;
397: }
398: }
399: *irootloc = dat->irootloc;
400: }
401: PetscFunctionReturn(PETSC_SUCCESS);
402: }
404: PETSC_INTERN PetscErrorCode PetscSFCreateLocalSF_Allgatherv(PetscSF sf, PetscSF *out)
405: {
406: PetscInt i, nroots, nleaves, rstart, *ilocal;
407: PetscSFNode *iremote;
408: PetscSF lsf;
410: PetscFunctionBegin;
411: nleaves = sf->nleaves ? sf->nroots : 0; /* sf->nleaves can be zero with SFGather(v) */
412: nroots = nleaves;
413: PetscCall(PetscMalloc1(nleaves, &ilocal));
414: PetscCall(PetscMalloc1(nleaves, &iremote));
415: PetscCall(PetscLayoutGetRange(sf->map, &rstart, NULL));
417: for (i = 0; i < nleaves; i++) {
418: ilocal[i] = rstart + i; /* lsf does not change leave indices */
419: iremote[i].rank = 0; /* rank in PETSC_COMM_SELF */
420: iremote[i].index = i; /* root index */
421: }
423: PetscCall(PetscSFCreate(PETSC_COMM_SELF, &lsf));
424: PetscCall(PetscSFSetGraph(lsf, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
425: PetscCall(PetscSFSetUp(lsf));
426: *out = lsf;
427: PetscFunctionReturn(PETSC_SUCCESS);
428: }
430: PETSC_INTERN PetscErrorCode PetscSFCreate_Allgatherv(PetscSF sf)
431: {
432: PetscSF_Allgatherv *dat = (PetscSF_Allgatherv *)sf->data;
434: PetscFunctionBegin;
435: sf->ops->BcastEnd = PetscSFBcastEnd_Basic;
436: sf->ops->ReduceEnd = PetscSFReduceEnd_Allgatherv;
438: sf->ops->SetUp = PetscSFSetUp_Allgatherv;
439: sf->ops->Reset = PetscSFReset_Allgatherv;
440: sf->ops->Destroy = PetscSFDestroy_Allgatherv;
441: sf->ops->GetRootRanks = PetscSFGetRootRanks_Allgatherv;
442: sf->ops->GetLeafRanks = PetscSFGetLeafRanks_Allgatherv;
443: sf->ops->GetGraph = PetscSFGetGraph_Allgatherv;
444: sf->ops->BcastBegin = PetscSFBcastBegin_Allgatherv;
445: sf->ops->ReduceBegin = PetscSFReduceBegin_Allgatherv;
446: sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Allgatherv;
447: sf->ops->FetchAndOpEnd = PetscSFFetchAndOpEnd_Allgatherv;
448: sf->ops->CreateLocalSF = PetscSFCreateLocalSF_Allgatherv;
449: sf->ops->BcastToZero = PetscSFBcastToZero_Allgatherv;
451: sf->collective = PETSC_TRUE;
453: PetscCall(PetscNew(&dat));
454: dat->bcast_root = -1;
455: sf->data = (void *)dat;
456: PetscFunctionReturn(PETSC_SUCCESS);
457: }