xref: /petsc/src/vec/is/sf/impls/basic/sfmpi.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
171438e86SJunchao Zhang /* Mainly for MPI_Isend in SFBASIC. Once SFNEIGHBOR, SFALLGHATERV etc have a persistent version,
271438e86SJunchao Zhang    we can also do abstractions like Prepare/StartCommunication.
371438e86SJunchao Zhang */
471438e86SJunchao Zhang 
571438e86SJunchao Zhang #include <../src/vec/is/sf/impls/basic/sfpack.h>
671438e86SJunchao Zhang 
771438e86SJunchao Zhang /* Start MPI requests. If use non-GPU aware MPI, we might need to copy data from device buf to host buf */
8d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFLinkStartRequests_MPI(PetscSF sf, PetscSFLink link, PetscSFDirection direction)
9d71ae5a4SJacob Faibussowitsch {
1071438e86SJunchao Zhang   PetscMPIInt    nreqs;
1171438e86SJunchao Zhang   MPI_Request   *reqs = NULL;
1271438e86SJunchao Zhang   PetscSF_Basic *bas  = (PetscSF_Basic *)sf->data;
1371438e86SJunchao Zhang   PetscInt       buflen;
1471438e86SJunchao Zhang 
1571438e86SJunchao Zhang   PetscFunctionBegin;
1671438e86SJunchao Zhang   buflen = (direction == PETSCSF_ROOT2LEAF) ? sf->leafbuflen[PETSCSF_REMOTE] : bas->rootbuflen[PETSCSF_REMOTE];
1771438e86SJunchao Zhang   if (buflen) {
1871438e86SJunchao Zhang     if (direction == PETSCSF_ROOT2LEAF) {
1971438e86SJunchao Zhang       nreqs = sf->nleafreqs;
209566063dSJacob Faibussowitsch       PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, direction, NULL, NULL, NULL, &reqs));
2171438e86SJunchao Zhang     } else { /* leaf to root */
2271438e86SJunchao Zhang       nreqs = bas->nrootreqs;
239566063dSJacob Faibussowitsch       PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, direction, NULL, NULL, &reqs, NULL));
2471438e86SJunchao Zhang     }
259566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Startall_irecv(buflen, link->unit, nreqs, reqs));
2671438e86SJunchao Zhang   }
2771438e86SJunchao Zhang 
2871438e86SJunchao Zhang   buflen = (direction == PETSCSF_ROOT2LEAF) ? bas->rootbuflen[PETSCSF_REMOTE] : sf->leafbuflen[PETSCSF_REMOTE];
2971438e86SJunchao Zhang   if (buflen) {
3071438e86SJunchao Zhang     if (direction == PETSCSF_ROOT2LEAF) {
3171438e86SJunchao Zhang       nreqs = bas->nrootreqs;
329566063dSJacob Faibussowitsch       PetscCall(PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE /*device2host before sending */));
339566063dSJacob Faibussowitsch       PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, direction, NULL, NULL, &reqs, NULL));
3471438e86SJunchao Zhang     } else { /* leaf to root */
3571438e86SJunchao Zhang       nreqs = sf->nleafreqs;
369566063dSJacob Faibussowitsch       PetscCall(PetscSFLinkCopyLeafBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_TRUE));
379566063dSJacob Faibussowitsch       PetscCall(PetscSFLinkGetMPIBuffersAndRequests(sf, link, direction, NULL, NULL, NULL, &reqs));
3871438e86SJunchao Zhang     }
399566063dSJacob Faibussowitsch     PetscCall(PetscSFLinkSyncStreamBeforeCallMPI(sf, link, direction));
409566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Startall_isend(buflen, link->unit, nreqs, reqs));
4171438e86SJunchao Zhang   }
42*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4371438e86SJunchao Zhang }
4471438e86SJunchao Zhang 
45d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFLinkWaitRequests_MPI(PetscSF sf, PetscSFLink link, PetscSFDirection direction)
46d71ae5a4SJacob Faibussowitsch {
4771438e86SJunchao Zhang   PetscSF_Basic     *bas           = (PetscSF_Basic *)sf->data;
4871438e86SJunchao Zhang   const PetscMemType rootmtype_mpi = link->rootmtype_mpi, leafmtype_mpi = link->leafmtype_mpi;
4971438e86SJunchao Zhang   const PetscInt     rootdirect_mpi = link->rootdirect_mpi, leafdirect_mpi = link->leafdirect_mpi;
5071438e86SJunchao Zhang 
5171438e86SJunchao Zhang   PetscFunctionBegin;
529566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(bas->nrootreqs, link->rootreqs[direction][rootmtype_mpi][rootdirect_mpi], MPI_STATUSES_IGNORE));
539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Waitall(sf->nleafreqs, link->leafreqs[direction][leafmtype_mpi][leafdirect_mpi], MPI_STATUSES_IGNORE));
5471438e86SJunchao Zhang   if (direction == PETSCSF_ROOT2LEAF) {
559566063dSJacob Faibussowitsch     PetscCall(PetscSFLinkCopyLeafBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_FALSE /* host2device after recving */));
5671438e86SJunchao Zhang   } else {
579566063dSJacob Faibussowitsch     PetscCall(PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_FALSE));
5871438e86SJunchao Zhang   }
59*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6071438e86SJunchao Zhang }
6171438e86SJunchao Zhang 
6271438e86SJunchao Zhang /*
6371438e86SJunchao Zhang    The routine Creates a communication link for the given operation. It first looks up its link cache. If
6471438e86SJunchao Zhang    there is a free & suitable one, it uses it. Otherwise it creates a new one.
6571438e86SJunchao Zhang 
6671438e86SJunchao Zhang    A link contains buffers and MPI requests for send/recv. It also contains pack/unpack routines to pack/unpack
6771438e86SJunchao Zhang    root/leafdata to/from these buffers. Buffers are allocated at our discretion. When we find root/leafata
6871438e86SJunchao Zhang    can be directly passed to MPI, we won't allocate them. Even we allocate buffers, we only allocate
6971438e86SJunchao Zhang    those that are needed by the given `sfop` and `op`, in other words, we do lazy memory-allocation.
7071438e86SJunchao Zhang 
7171438e86SJunchao Zhang    The routine also allocates buffers on CPU when one does not use gpu-aware MPI but data is on GPU.
7271438e86SJunchao Zhang 
7371438e86SJunchao Zhang    In SFBasic, MPI requests are persistent. They are init'ed until we try to get requests from a link.
7471438e86SJunchao Zhang 
7571438e86SJunchao Zhang    The routine is shared by SFBasic and SFNeighbor based on the fact they all deal with sparse graphs and
7671438e86SJunchao Zhang    need pack/unpack data.
7771438e86SJunchao Zhang */
78d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFLinkCreate_MPI(PetscSF sf, MPI_Datatype unit, PetscMemType xrootmtype, const void *rootdata, PetscMemType xleafmtype, const void *leafdata, MPI_Op op, PetscSFOperation sfop, PetscSFLink *mylink)
79d71ae5a4SJacob Faibussowitsch {
8071438e86SJunchao Zhang   PetscSF_Basic   *bas = (PetscSF_Basic *)sf->data;
8171438e86SJunchao Zhang   PetscInt         i, j, k, nrootreqs, nleafreqs, nreqs;
8271438e86SJunchao Zhang   PetscSFLink     *p, link;
8371438e86SJunchao Zhang   PetscSFDirection direction;
8471438e86SJunchao Zhang   MPI_Request     *reqs = NULL;
8571438e86SJunchao Zhang   PetscBool        match, rootdirect[2], leafdirect[2];
8671438e86SJunchao Zhang   PetscMemType     rootmtype = PetscMemTypeHost(xrootmtype) ? PETSC_MEMTYPE_HOST : PETSC_MEMTYPE_DEVICE; /* Convert to 0/1 as we will use it in subscript */
8771438e86SJunchao Zhang   PetscMemType     leafmtype = PetscMemTypeHost(xleafmtype) ? PETSC_MEMTYPE_HOST : PETSC_MEMTYPE_DEVICE;
8871438e86SJunchao Zhang   PetscMemType     rootmtype_mpi, leafmtype_mpi;   /* mtypes seen by MPI */
8971438e86SJunchao Zhang   PetscInt         rootdirect_mpi, leafdirect_mpi; /* root/leafdirect seen by MPI*/
9071438e86SJunchao Zhang 
9171438e86SJunchao Zhang   PetscFunctionBegin;
9271438e86SJunchao Zhang 
9371438e86SJunchao Zhang   /* Can we directly use root/leafdirect with the given sf, sfop and op? */
9471438e86SJunchao Zhang   for (i = PETSCSF_LOCAL; i <= PETSCSF_REMOTE; i++) {
9571438e86SJunchao Zhang     if (sfop == PETSCSF_BCAST) {
9671438e86SJunchao Zhang       rootdirect[i] = bas->rootcontig[i];                                                  /* Pack roots */
9771438e86SJunchao Zhang       leafdirect[i] = (sf->leafcontig[i] && op == MPI_REPLACE) ? PETSC_TRUE : PETSC_FALSE; /* Unpack leaves */
9871438e86SJunchao Zhang     } else if (sfop == PETSCSF_REDUCE) {
9971438e86SJunchao Zhang       leafdirect[i] = sf->leafcontig[i];                                                    /* Pack leaves */
10071438e86SJunchao Zhang       rootdirect[i] = (bas->rootcontig[i] && op == MPI_REPLACE) ? PETSC_TRUE : PETSC_FALSE; /* Unpack roots */
10171438e86SJunchao Zhang     } else {                                                                                /* PETSCSF_FETCH */
10271438e86SJunchao Zhang       rootdirect[i] = PETSC_FALSE;                                                          /* FETCH always need a separate rootbuf */
10371438e86SJunchao Zhang       leafdirect[i] = PETSC_FALSE;                                                          /* We also force allocating a separate leafbuf so that leafdata and leafupdate can share mpi requests */
10471438e86SJunchao Zhang     }
10571438e86SJunchao Zhang   }
10671438e86SJunchao Zhang 
10771438e86SJunchao Zhang   if (sf->use_gpu_aware_mpi) {
10871438e86SJunchao Zhang     rootmtype_mpi = rootmtype;
10971438e86SJunchao Zhang     leafmtype_mpi = leafmtype;
11071438e86SJunchao Zhang   } else {
11171438e86SJunchao Zhang     rootmtype_mpi = leafmtype_mpi = PETSC_MEMTYPE_HOST;
11271438e86SJunchao Zhang   }
11371438e86SJunchao Zhang   /* Will root/leafdata be directly accessed by MPI?  Without use_gpu_aware_mpi, device data is bufferred on host and then passed to MPI */
11471438e86SJunchao Zhang   rootdirect_mpi = rootdirect[PETSCSF_REMOTE] && (rootmtype_mpi == rootmtype) ? 1 : 0;
11571438e86SJunchao Zhang   leafdirect_mpi = leafdirect[PETSCSF_REMOTE] && (leafmtype_mpi == leafmtype) ? 1 : 0;
11671438e86SJunchao Zhang 
11771438e86SJunchao Zhang   direction = (sfop == PETSCSF_BCAST) ? PETSCSF_ROOT2LEAF : PETSCSF_LEAF2ROOT;
11871438e86SJunchao Zhang   nrootreqs = bas->nrootreqs;
11971438e86SJunchao Zhang   nleafreqs = sf->nleafreqs;
12071438e86SJunchao Zhang 
12171438e86SJunchao Zhang   /* Look for free links in cache */
12271438e86SJunchao Zhang   for (p = &bas->avail; (link = *p); p = &link->next) {
12371438e86SJunchao Zhang     if (!link->use_nvshmem) { /* Only check with MPI links */
1249566063dSJacob Faibussowitsch       PetscCall(MPIPetsc_Type_compare(unit, link->unit, &match));
12571438e86SJunchao Zhang       if (match) {
12671438e86SJunchao Zhang         /* If root/leafdata will be directly passed to MPI, test if the data used to initialized the MPI requests matches with the current.
12771438e86SJunchao Zhang            If not, free old requests. New requests will be lazily init'ed until one calls PetscSFLinkGetMPIBuffersAndRequests().
12871438e86SJunchao Zhang         */
12971438e86SJunchao Zhang         if (rootdirect_mpi && sf->persistent && link->rootreqsinited[direction][rootmtype][1] && link->rootdatadirect[direction][rootmtype] != rootdata) {
13071438e86SJunchao Zhang           reqs = link->rootreqs[direction][rootmtype][1]; /* Here, rootmtype = rootmtype_mpi */
1319371c9d4SSatish Balay           for (i = 0; i < nrootreqs; i++) {
1329371c9d4SSatish Balay             if (reqs[i] != MPI_REQUEST_NULL) PetscCallMPI(MPI_Request_free(&reqs[i]));
1339371c9d4SSatish Balay           }
13471438e86SJunchao Zhang           link->rootreqsinited[direction][rootmtype][1] = PETSC_FALSE;
13571438e86SJunchao Zhang         }
13671438e86SJunchao Zhang         if (leafdirect_mpi && sf->persistent && link->leafreqsinited[direction][leafmtype][1] && link->leafdatadirect[direction][leafmtype] != leafdata) {
13771438e86SJunchao Zhang           reqs = link->leafreqs[direction][leafmtype][1];
1389371c9d4SSatish Balay           for (i = 0; i < nleafreqs; i++) {
1399371c9d4SSatish Balay             if (reqs[i] != MPI_REQUEST_NULL) PetscCallMPI(MPI_Request_free(&reqs[i]));
1409371c9d4SSatish Balay           }
14171438e86SJunchao Zhang           link->leafreqsinited[direction][leafmtype][1] = PETSC_FALSE;
14271438e86SJunchao Zhang         }
14371438e86SJunchao Zhang         *p = link->next; /* Remove from available list */
14471438e86SJunchao Zhang         goto found;
14571438e86SJunchao Zhang       }
14671438e86SJunchao Zhang     }
14771438e86SJunchao Zhang   }
14871438e86SJunchao Zhang 
1499566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
1509566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkSetUp_Host(sf, link, unit));
1519566063dSJacob Faibussowitsch   PetscCall(PetscCommGetNewTag(PetscObjectComm((PetscObject)sf), &link->tag)); /* One tag per link */
15271438e86SJunchao Zhang 
15371438e86SJunchao Zhang   nreqs = (nrootreqs + nleafreqs) * 8;
1549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nreqs, &link->reqs));
15571438e86SJunchao Zhang   for (i = 0; i < nreqs; i++) link->reqs[i] = MPI_REQUEST_NULL; /* Initialized to NULL so that we know which need to be freed in Destroy */
15671438e86SJunchao Zhang 
15771438e86SJunchao Zhang   for (i = 0; i < 2; i++) {     /* Two communication directions */
15871438e86SJunchao Zhang     for (j = 0; j < 2; j++) {   /* Two memory types */
15971438e86SJunchao Zhang       for (k = 0; k < 2; k++) { /* root/leafdirect 0 or 1 */
16071438e86SJunchao Zhang         link->rootreqs[i][j][k] = link->reqs + nrootreqs * (4 * i + 2 * j + k);
16171438e86SJunchao Zhang         link->leafreqs[i][j][k] = link->reqs + nrootreqs * 8 + nleafreqs * (4 * i + 2 * j + k);
16271438e86SJunchao Zhang       }
16371438e86SJunchao Zhang     }
16471438e86SJunchao Zhang   }
16571438e86SJunchao Zhang   link->StartCommunication  = PetscSFLinkStartRequests_MPI;
16671438e86SJunchao Zhang   link->FinishCommunication = PetscSFLinkWaitRequests_MPI;
16771438e86SJunchao Zhang 
16871438e86SJunchao Zhang found:
16971438e86SJunchao Zhang 
17071438e86SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
17171438e86SJunchao Zhang   if ((PetscMemTypeDevice(xrootmtype) || PetscMemTypeDevice(xleafmtype)) && !link->deviceinited) {
17271438e86SJunchao Zhang   #if defined(PETSC_HAVE_CUDA)
1739566063dSJacob Faibussowitsch     if (sf->backend == PETSCSF_BACKEND_CUDA) PetscCall(PetscSFLinkSetUp_CUDA(sf, link, unit)); /* Setup streams etc */
17471438e86SJunchao Zhang   #endif
17571438e86SJunchao Zhang   #if defined(PETSC_HAVE_HIP)
1769566063dSJacob Faibussowitsch     if (sf->backend == PETSCSF_BACKEND_HIP) PetscCall(PetscSFLinkSetUp_HIP(sf, link, unit)); /* Setup streams etc */
17771438e86SJunchao Zhang   #endif
17871438e86SJunchao Zhang   #if defined(PETSC_HAVE_KOKKOS)
1799566063dSJacob Faibussowitsch     if (sf->backend == PETSCSF_BACKEND_KOKKOS) PetscCall(PetscSFLinkSetUp_Kokkos(sf, link, unit));
18071438e86SJunchao Zhang   #endif
18171438e86SJunchao Zhang   }
18271438e86SJunchao Zhang #endif
18371438e86SJunchao Zhang 
18471438e86SJunchao Zhang   /* Allocate buffers along root/leafdata */
18571438e86SJunchao Zhang   for (i = PETSCSF_LOCAL; i <= PETSCSF_REMOTE; i++) {
18671438e86SJunchao Zhang     /* For local communication, buffers are only needed when roots and leaves have different mtypes */
18771438e86SJunchao Zhang     if (i == PETSCSF_LOCAL && rootmtype == leafmtype) continue;
18871438e86SJunchao Zhang     if (bas->rootbuflen[i]) {
18971438e86SJunchao Zhang       if (rootdirect[i]) { /* Aha, we disguise rootdata as rootbuf */
19071438e86SJunchao Zhang         link->rootbuf[i][rootmtype] = (char *)rootdata + bas->rootstart[i] * link->unitbytes;
19171438e86SJunchao Zhang       } else { /* Have to have a separate rootbuf */
19248a46eb9SPierre Jolivet         if (!link->rootbuf_alloc[i][rootmtype]) PetscCall(PetscSFMalloc(sf, rootmtype, bas->rootbuflen[i] * link->unitbytes, (void **)&link->rootbuf_alloc[i][rootmtype]));
19371438e86SJunchao Zhang         link->rootbuf[i][rootmtype] = link->rootbuf_alloc[i][rootmtype];
19471438e86SJunchao Zhang       }
19571438e86SJunchao Zhang     }
19671438e86SJunchao Zhang 
19771438e86SJunchao Zhang     if (sf->leafbuflen[i]) {
19871438e86SJunchao Zhang       if (leafdirect[i]) {
19971438e86SJunchao Zhang         link->leafbuf[i][leafmtype] = (char *)leafdata + sf->leafstart[i] * link->unitbytes;
20071438e86SJunchao Zhang       } else {
20148a46eb9SPierre Jolivet         if (!link->leafbuf_alloc[i][leafmtype]) PetscCall(PetscSFMalloc(sf, leafmtype, sf->leafbuflen[i] * link->unitbytes, (void **)&link->leafbuf_alloc[i][leafmtype]));
20271438e86SJunchao Zhang         link->leafbuf[i][leafmtype] = link->leafbuf_alloc[i][leafmtype];
20371438e86SJunchao Zhang       }
20471438e86SJunchao Zhang     }
20571438e86SJunchao Zhang   }
20671438e86SJunchao Zhang 
20771438e86SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
20871438e86SJunchao Zhang   /* Allocate buffers on host for buffering data on device in cast not use_gpu_aware_mpi */
20971438e86SJunchao Zhang   if (PetscMemTypeDevice(rootmtype) && PetscMemTypeHost(rootmtype_mpi)) {
21048a46eb9SPierre Jolivet     if (!link->rootbuf_alloc[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST]) PetscCall(PetscMalloc(bas->rootbuflen[PETSCSF_REMOTE] * link->unitbytes, &link->rootbuf_alloc[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST]));
21171438e86SJunchao Zhang     link->rootbuf[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST] = link->rootbuf_alloc[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST];
21271438e86SJunchao Zhang   }
21371438e86SJunchao Zhang   if (PetscMemTypeDevice(leafmtype) && PetscMemTypeHost(leafmtype_mpi)) {
21448a46eb9SPierre Jolivet     if (!link->leafbuf_alloc[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST]) PetscCall(PetscMalloc(sf->leafbuflen[PETSCSF_REMOTE] * link->unitbytes, &link->leafbuf_alloc[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST]));
21571438e86SJunchao Zhang     link->leafbuf[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST] = link->leafbuf_alloc[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST];
21671438e86SJunchao Zhang   }
21771438e86SJunchao Zhang #endif
21871438e86SJunchao Zhang 
21971438e86SJunchao Zhang   /* Set `current` state of the link. They may change between different SF invocations with the same link */
22071438e86SJunchao Zhang   if (sf->persistent) { /* If data is directly passed to MPI and inits MPI requests, record the data for comparison on future invocations */
22171438e86SJunchao Zhang     if (rootdirect_mpi) link->rootdatadirect[direction][rootmtype] = rootdata;
22271438e86SJunchao Zhang     if (leafdirect_mpi) link->leafdatadirect[direction][leafmtype] = leafdata;
22371438e86SJunchao Zhang   }
22471438e86SJunchao Zhang 
22571438e86SJunchao Zhang   link->rootdata = rootdata; /* root/leafdata are keys to look up links in PetscSFXxxEnd */
22671438e86SJunchao Zhang   link->leafdata = leafdata;
22771438e86SJunchao Zhang   for (i = PETSCSF_LOCAL; i <= PETSCSF_REMOTE; i++) {
22871438e86SJunchao Zhang     link->rootdirect[i] = rootdirect[i];
22971438e86SJunchao Zhang     link->leafdirect[i] = leafdirect[i];
23071438e86SJunchao Zhang   }
23171438e86SJunchao Zhang   link->rootdirect_mpi = rootdirect_mpi;
23271438e86SJunchao Zhang   link->leafdirect_mpi = leafdirect_mpi;
23371438e86SJunchao Zhang   link->rootmtype      = rootmtype;
23471438e86SJunchao Zhang   link->leafmtype      = leafmtype;
23571438e86SJunchao Zhang   link->rootmtype_mpi  = rootmtype_mpi;
23671438e86SJunchao Zhang   link->leafmtype_mpi  = leafmtype_mpi;
23771438e86SJunchao Zhang 
23871438e86SJunchao Zhang   link->next = bas->inuse;
23971438e86SJunchao Zhang   bas->inuse = link;
24071438e86SJunchao Zhang   *mylink    = link;
241*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24271438e86SJunchao Zhang }
243