xref: /petsc/src/vec/is/sf/impls/basic/sfmpi.c (revision 44be734963a6af33561331b2d89af9cfd1caf25e)
171438e86SJunchao Zhang #include <../src/vec/is/sf/impls/basic/sfpack.h>
271438e86SJunchao Zhang 
3d8b4a066SPierre Jolivet // Though there is no default mechanism to start a communication, we have a
4f5d27ee7SJunchao Zhang // default to finish communication, which is just waiting on the requests.
5f5d27ee7SJunchao Zhang // It should work for both non-blocking or persistent send/recvs or collectivwes.
PetscSFLinkFinishCommunication_Default(PetscSF sf,PetscSFLink link,PetscSFDirection direction)6f5d27ee7SJunchao Zhang static PetscErrorCode PetscSFLinkFinishCommunication_Default(PetscSF sf, PetscSFLink link, PetscSFDirection direction)
7d71ae5a4SJacob Faibussowitsch {
871438e86SJunchao Zhang   PetscSF_Basic     *bas           = (PetscSF_Basic *)sf->data;
971438e86SJunchao Zhang   const PetscMemType rootmtype_mpi = link->rootmtype_mpi, leafmtype_mpi = link->leafmtype_mpi;
1071438e86SJunchao Zhang   const PetscInt     rootdirect_mpi = link->rootdirect_mpi, leafdirect_mpi = link->leafdirect_mpi;
1171438e86SJunchao Zhang 
1271438e86SJunchao Zhang   PetscFunctionBegin;
13*f9334340SJunchao Zhang   if (sf->monitor) {
14*f9334340SJunchao Zhang     PetscMPIInt rank;
15*f9334340SJunchao Zhang     const char *rootaction = (direction == PETSCSF_ROOT2LEAF) ? "sending to  " : "recving from";
16*f9334340SJunchao Zhang     const char *leafaction = (direction == PETSCSF_ROOT2LEAF) ? "recving from" : "sending to  ";
17*f9334340SJunchao Zhang     const char *sfaction   = (direction == PETSCSF_ROOT2LEAF) ? "PetscSFBcast" : "PetscSFReduce";
18*f9334340SJunchao Zhang 
19*f9334340SJunchao Zhang     PetscCall(PetscPrintf(PETSC_COMM_SELF, "---------------- Begin %s communication -------------------\n", sfaction));
20*f9334340SJunchao Zhang     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
21*f9334340SJunchao Zhang 
22*f9334340SJunchao Zhang     for (PetscMPIInt i = 0; i < bas->nrootreqs; i++) {
23*f9334340SJunchao Zhang       size_t size = (bas->ioffset[i + bas->ndiranks + 1] - bas->ioffset[i + bas->ndiranks]) * link->unitbytes;
24*f9334340SJunchao Zhang       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Rank %6d %s Rank %6d (%16zu bytes) with MPI tag %10d ... ", rank, rootaction, bas->iranks[i + bas->ndiranks], size, link->tag));
25*f9334340SJunchao Zhang       PetscCallMPI(MPI_Wait(link->rootreqs[direction][rootmtype_mpi][rootdirect_mpi] + i, MPI_STATUS_IGNORE));
26*f9334340SJunchao Zhang       PetscCall(PetscPrintf(PETSC_COMM_SELF, "DONE\n"));
27*f9334340SJunchao Zhang     }
28*f9334340SJunchao Zhang     for (PetscMPIInt i = 0; i < sf->nleafreqs; i++) {
29*f9334340SJunchao Zhang       size_t size = (sf->roffset[i + sf->ndranks + 1] - sf->roffset[i + sf->ndranks]) * link->unitbytes;
30*f9334340SJunchao Zhang       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Rank %6d %s Rank %6d (%16zu bytes) with MPI tag %10d ... ", rank, leafaction, sf->ranks[i + sf->ndranks], size, link->tag));
31*f9334340SJunchao Zhang       PetscCallMPI(MPI_Wait(link->leafreqs[direction][leafmtype_mpi][leafdirect_mpi] + i, MPI_STATUS_IGNORE));
32*f9334340SJunchao Zhang       PetscCall(PetscPrintf(PETSC_COMM_SELF, "DONE\n"));
33*f9334340SJunchao Zhang     }
34*f9334340SJunchao Zhang     PetscCall(PetscPrintf(PETSC_COMM_SELF, "---------------- End   %s communication -------------------\n\n", sfaction));
35*f9334340SJunchao Zhang   } else {
36f5d27ee7SJunchao Zhang     if (bas->nrootreqs) PetscCallMPI(MPI_Waitall(bas->nrootreqs, link->rootreqs[direction][rootmtype_mpi][rootdirect_mpi], MPI_STATUSES_IGNORE));
37f5d27ee7SJunchao Zhang     if (sf->nleafreqs) PetscCallMPI(MPI_Waitall(sf->nleafreqs, link->leafreqs[direction][leafmtype_mpi][leafdirect_mpi], MPI_STATUSES_IGNORE));
38*f9334340SJunchao Zhang   }
39*f9334340SJunchao Zhang 
4071438e86SJunchao Zhang   if (direction == PETSCSF_ROOT2LEAF) {
419566063dSJacob Faibussowitsch     PetscCall(PetscSFLinkCopyLeafBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_FALSE /* host2device after recving */));
4271438e86SJunchao Zhang   } else {
439566063dSJacob Faibussowitsch     PetscCall(PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI(sf, link, PETSC_FALSE));
4471438e86SJunchao Zhang   }
453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4671438e86SJunchao Zhang }
4771438e86SJunchao Zhang 
4871438e86SJunchao Zhang /*
4971438e86SJunchao Zhang    The routine Creates a communication link for the given operation. It first looks up its link cache. If
5071438e86SJunchao Zhang    there is a free & suitable one, it uses it. Otherwise it creates a new one.
5171438e86SJunchao Zhang 
5271438e86SJunchao Zhang    A link contains buffers and MPI requests for send/recv. It also contains pack/unpack routines to pack/unpack
5371438e86SJunchao Zhang    root/leafdata to/from these buffers. Buffers are allocated at our discretion. When we find root/leafata
5471438e86SJunchao Zhang    can be directly passed to MPI, we won't allocate them. Even we allocate buffers, we only allocate
5571438e86SJunchao Zhang    those that are needed by the given `sfop` and `op`, in other words, we do lazy memory-allocation.
5671438e86SJunchao Zhang 
5771438e86SJunchao Zhang    The routine also allocates buffers on CPU when one does not use gpu-aware MPI but data is on GPU.
5871438e86SJunchao Zhang 
5971438e86SJunchao Zhang    In SFBasic, MPI requests are persistent. They are init'ed until we try to get requests from a link.
6071438e86SJunchao Zhang 
6171438e86SJunchao Zhang    The routine is shared by SFBasic and SFNeighbor based on the fact they all deal with sparse graphs and
6271438e86SJunchao Zhang    need pack/unpack data.
6371438e86SJunchao Zhang */
PetscSFLinkCreate_MPI(PetscSF sf,MPI_Datatype unit,PetscMemType xrootmtype,const void * rootdata,PetscMemType xleafmtype,const void * leafdata,MPI_Op op,PetscSFOperation sfop,PetscSFLink * mylink)64d71ae5a4SJacob 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)
65d71ae5a4SJacob Faibussowitsch {
6671438e86SJunchao Zhang   PetscSF_Basic   *bas = (PetscSF_Basic *)sf->data;
6771438e86SJunchao Zhang   PetscInt         i, j, k, nrootreqs, nleafreqs, nreqs;
6871438e86SJunchao Zhang   PetscSFLink     *p, link;
6971438e86SJunchao Zhang   PetscSFDirection direction;
7071438e86SJunchao Zhang   MPI_Request     *reqs = NULL;
7171438e86SJunchao Zhang   PetscBool        match, rootdirect[2], leafdirect[2];
7271438e86SJunchao Zhang   PetscMemType     rootmtype = PetscMemTypeHost(xrootmtype) ? PETSC_MEMTYPE_HOST : PETSC_MEMTYPE_DEVICE; /* Convert to 0/1 as we will use it in subscript */
7371438e86SJunchao Zhang   PetscMemType     leafmtype = PetscMemTypeHost(xleafmtype) ? PETSC_MEMTYPE_HOST : PETSC_MEMTYPE_DEVICE;
7471438e86SJunchao Zhang   PetscMemType     rootmtype_mpi, leafmtype_mpi;   /* mtypes seen by MPI */
7571438e86SJunchao Zhang   PetscInt         rootdirect_mpi, leafdirect_mpi; /* root/leafdirect seen by MPI*/
7671438e86SJunchao Zhang 
7771438e86SJunchao Zhang   PetscFunctionBegin;
7871438e86SJunchao Zhang   /* Can we directly use root/leafdirect with the given sf, sfop and op? */
7971438e86SJunchao Zhang   for (i = PETSCSF_LOCAL; i <= PETSCSF_REMOTE; i++) {
8071438e86SJunchao Zhang     if (sfop == PETSCSF_BCAST) {
8171438e86SJunchao Zhang       rootdirect[i] = bas->rootcontig[i];                                                  /* Pack roots */
8271438e86SJunchao Zhang       leafdirect[i] = (sf->leafcontig[i] && op == MPI_REPLACE) ? PETSC_TRUE : PETSC_FALSE; /* Unpack leaves */
8371438e86SJunchao Zhang     } else if (sfop == PETSCSF_REDUCE) {
8471438e86SJunchao Zhang       leafdirect[i] = sf->leafcontig[i];                                                    /* Pack leaves */
8571438e86SJunchao Zhang       rootdirect[i] = (bas->rootcontig[i] && op == MPI_REPLACE) ? PETSC_TRUE : PETSC_FALSE; /* Unpack roots */
8671438e86SJunchao Zhang     } else {                                                                                /* PETSCSF_FETCH */
8771438e86SJunchao Zhang       rootdirect[i] = PETSC_FALSE;                                                          /* FETCH always need a separate rootbuf */
8871438e86SJunchao Zhang       leafdirect[i] = PETSC_FALSE;                                                          /* We also force allocating a separate leafbuf so that leafdata and leafupdate can share mpi requests */
8971438e86SJunchao Zhang     }
9071438e86SJunchao Zhang   }
9171438e86SJunchao Zhang 
926677b1c1SJunchao Zhang   // NEVER use root/leafdirect[] for persistent collectives. Otherwise, suppose for the first time, all ranks build
936677b1c1SJunchao Zhang   // a persistent MPI request in a collective call. Then in a second call to PetscSFBcast, one rank uses root/leafdirect
946677b1c1SJunchao Zhang   // but with new rootdata/leafdata pointers. Other ranks keep using the same rootdata/leafdata pointers as last time.
956677b1c1SJunchao Zhang   // Only that rank will try to rebuild the request with a collective call, resulting in hanging. We could to call
966677b1c1SJunchao Zhang   // MPI_Allreduce() every time to detect changes in root/leafdata, but that is too expensive for sparse communication.
976677b1c1SJunchao Zhang   // So we always set root/leafdirect[] to false and allocate additional root/leaf buffers for persistent collectives.
986677b1c1SJunchao Zhang   if (sf->persistent && sf->collective) {
996677b1c1SJunchao Zhang     rootdirect[PETSCSF_REMOTE] = PETSC_FALSE;
1006677b1c1SJunchao Zhang     leafdirect[PETSCSF_REMOTE] = PETSC_FALSE;
1016677b1c1SJunchao Zhang   }
1026677b1c1SJunchao Zhang 
10371438e86SJunchao Zhang   if (sf->use_gpu_aware_mpi) {
10471438e86SJunchao Zhang     rootmtype_mpi = rootmtype;
10571438e86SJunchao Zhang     leafmtype_mpi = leafmtype;
10671438e86SJunchao Zhang   } else {
10771438e86SJunchao Zhang     rootmtype_mpi = leafmtype_mpi = PETSC_MEMTYPE_HOST;
10871438e86SJunchao Zhang   }
109da81f932SPierre Jolivet   /* Will root/leafdata be directly accessed by MPI?  Without use_gpu_aware_mpi, device data is buffered on host and then passed to MPI */
11071438e86SJunchao Zhang   rootdirect_mpi = rootdirect[PETSCSF_REMOTE] && (rootmtype_mpi == rootmtype) ? 1 : 0;
11171438e86SJunchao Zhang   leafdirect_mpi = leafdirect[PETSCSF_REMOTE] && (leafmtype_mpi == leafmtype) ? 1 : 0;
11271438e86SJunchao Zhang 
11371438e86SJunchao Zhang   direction = (sfop == PETSCSF_BCAST) ? PETSCSF_ROOT2LEAF : PETSCSF_LEAF2ROOT;
11471438e86SJunchao Zhang   nrootreqs = bas->nrootreqs;
11571438e86SJunchao Zhang   nleafreqs = sf->nleafreqs;
11671438e86SJunchao Zhang 
11771438e86SJunchao Zhang   /* Look for free links in cache */
11871438e86SJunchao Zhang   for (p = &bas->avail; (link = *p); p = &link->next) {
11971438e86SJunchao Zhang     if (!link->use_nvshmem) { /* Only check with MPI links */
1209566063dSJacob Faibussowitsch       PetscCall(MPIPetsc_Type_compare(unit, link->unit, &match));
12171438e86SJunchao Zhang       if (match) {
12271438e86SJunchao Zhang         /* If root/leafdata will be directly passed to MPI, test if the data used to initialized the MPI requests matches with the current.
1236677b1c1SJunchao Zhang            If not, free old requests. New requests will be lazily init'ed until one calls PetscSFLinkGetMPIBuffersAndRequests() with the same tag.
12471438e86SJunchao Zhang         */
12571438e86SJunchao Zhang         if (rootdirect_mpi && sf->persistent && link->rootreqsinited[direction][rootmtype][1] && link->rootdatadirect[direction][rootmtype] != rootdata) {
12671438e86SJunchao Zhang           reqs = link->rootreqs[direction][rootmtype][1]; /* Here, rootmtype = rootmtype_mpi */
1279371c9d4SSatish Balay           for (i = 0; i < nrootreqs; i++) {
1289371c9d4SSatish Balay             if (reqs[i] != MPI_REQUEST_NULL) PetscCallMPI(MPI_Request_free(&reqs[i]));
1299371c9d4SSatish Balay           }
13071438e86SJunchao Zhang           link->rootreqsinited[direction][rootmtype][1] = PETSC_FALSE;
13171438e86SJunchao Zhang         }
13271438e86SJunchao Zhang         if (leafdirect_mpi && sf->persistent && link->leafreqsinited[direction][leafmtype][1] && link->leafdatadirect[direction][leafmtype] != leafdata) {
13371438e86SJunchao Zhang           reqs = link->leafreqs[direction][leafmtype][1];
1349371c9d4SSatish Balay           for (i = 0; i < nleafreqs; i++) {
1359371c9d4SSatish Balay             if (reqs[i] != MPI_REQUEST_NULL) PetscCallMPI(MPI_Request_free(&reqs[i]));
1369371c9d4SSatish Balay           }
13771438e86SJunchao Zhang           link->leafreqsinited[direction][leafmtype][1] = PETSC_FALSE;
13871438e86SJunchao Zhang         }
13971438e86SJunchao Zhang         *p = link->next; /* Remove from available list */
14071438e86SJunchao Zhang         goto found;
14171438e86SJunchao Zhang       }
14271438e86SJunchao Zhang     }
14371438e86SJunchao Zhang   }
14471438e86SJunchao Zhang 
1459566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
1469566063dSJacob Faibussowitsch   PetscCall(PetscSFLinkSetUp_Host(sf, link, unit));
1479566063dSJacob Faibussowitsch   PetscCall(PetscCommGetNewTag(PetscObjectComm((PetscObject)sf), &link->tag)); /* One tag per link */
14871438e86SJunchao Zhang 
14971438e86SJunchao Zhang   nreqs = (nrootreqs + nleafreqs) * 8;
1509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nreqs, &link->reqs));
15171438e86SJunchao 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 */
15271438e86SJunchao Zhang 
1535c0db29aSPierre Jolivet   if (nreqs)
15471438e86SJunchao Zhang     for (i = 0; i < 2; i++) {     /* Two communication directions */
15571438e86SJunchao Zhang       for (j = 0; j < 2; j++) {   /* Two memory types */
15671438e86SJunchao Zhang         for (k = 0; k < 2; k++) { /* root/leafdirect 0 or 1 */
15771438e86SJunchao Zhang           link->rootreqs[i][j][k] = link->reqs + nrootreqs * (4 * i + 2 * j + k);
15871438e86SJunchao Zhang           link->leafreqs[i][j][k] = link->reqs + nrootreqs * 8 + nleafreqs * (4 * i + 2 * j + k);
15971438e86SJunchao Zhang         }
16071438e86SJunchao Zhang       }
16171438e86SJunchao Zhang     }
162f5d27ee7SJunchao Zhang 
163f5d27ee7SJunchao Zhang   link->FinishCommunication = PetscSFLinkFinishCommunication_Default;
164f5d27ee7SJunchao Zhang   // each SF type could customize their communication by setting function pointers in the link.
165f5d27ee7SJunchao Zhang   // Currently only BASIC and NEIGHBOR use this abstraction.
166f5d27ee7SJunchao Zhang   PetscTryTypeMethod(sf, SetCommunicationOps, link);
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;
2413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24271438e86SJunchao Zhang }
243