xref: /petsc/src/vec/is/sf/impls/basic/sfmpi.c (revision 71438e86cf8dea9f708c08733b37fef8eb68dc06)
1*71438e86SJunchao Zhang /* Mainly for MPI_Isend in SFBASIC. Once SFNEIGHBOR, SFALLGHATERV etc have a persistent version,
2*71438e86SJunchao Zhang    we can also do abstractions like Prepare/StartCommunication.
3*71438e86SJunchao Zhang */
4*71438e86SJunchao Zhang 
5*71438e86SJunchao Zhang #include <../src/vec/is/sf/impls/basic/sfpack.h>
6*71438e86SJunchao Zhang 
7*71438e86SJunchao Zhang /* Start MPI requests. If use non-GPU aware MPI, we might need to copy data from device buf to host buf */
8*71438e86SJunchao Zhang static PetscErrorCode PetscSFLinkStartRequests_MPI(PetscSF sf,PetscSFLink link,PetscSFDirection direction)
9*71438e86SJunchao Zhang {
10*71438e86SJunchao Zhang   PetscErrorCode    ierr;
11*71438e86SJunchao Zhang   PetscMPIInt       nreqs;
12*71438e86SJunchao Zhang   MPI_Request       *reqs = NULL;
13*71438e86SJunchao Zhang   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
14*71438e86SJunchao Zhang   PetscInt          buflen;
15*71438e86SJunchao Zhang 
16*71438e86SJunchao Zhang   PetscFunctionBegin;
17*71438e86SJunchao Zhang   buflen = (direction == PETSCSF_ROOT2LEAF) ? sf->leafbuflen[PETSCSF_REMOTE] : bas->rootbuflen[PETSCSF_REMOTE];
18*71438e86SJunchao Zhang   if (buflen) {
19*71438e86SJunchao Zhang     if (direction == PETSCSF_ROOT2LEAF) {
20*71438e86SJunchao Zhang       nreqs = sf->nleafreqs;
21*71438e86SJunchao Zhang       ierr = PetscSFLinkGetMPIBuffersAndRequests(sf,link,direction,NULL,NULL,NULL,&reqs);CHKERRQ(ierr);
22*71438e86SJunchao Zhang     } else { /* leaf to root */
23*71438e86SJunchao Zhang       nreqs = bas->nrootreqs;
24*71438e86SJunchao Zhang       ierr = PetscSFLinkGetMPIBuffersAndRequests(sf,link,direction,NULL,NULL,&reqs,NULL);CHKERRQ(ierr);
25*71438e86SJunchao Zhang     }
26*71438e86SJunchao Zhang     ierr = MPI_Startall_irecv(buflen,link->unit,nreqs,reqs);CHKERRMPI(ierr);
27*71438e86SJunchao Zhang   }
28*71438e86SJunchao Zhang 
29*71438e86SJunchao Zhang   buflen = (direction == PETSCSF_ROOT2LEAF) ? bas->rootbuflen[PETSCSF_REMOTE] : sf->leafbuflen[PETSCSF_REMOTE];
30*71438e86SJunchao Zhang   if (buflen) {
31*71438e86SJunchao Zhang     if (direction == PETSCSF_ROOT2LEAF) {
32*71438e86SJunchao Zhang       nreqs  = bas->nrootreqs;
33*71438e86SJunchao Zhang       ierr   = PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI(sf,link,PETSC_TRUE/*device2host before sending */);CHKERRQ(ierr);
34*71438e86SJunchao Zhang       ierr   = PetscSFLinkGetMPIBuffersAndRequests(sf,link,direction,NULL,NULL,&reqs,NULL);CHKERRQ(ierr);
35*71438e86SJunchao Zhang     } else { /* leaf to root */
36*71438e86SJunchao Zhang       nreqs  = sf->nleafreqs;
37*71438e86SJunchao Zhang       ierr   = PetscSFLinkCopyLeafBufferInCaseNotUseGpuAwareMPI(sf,link,PETSC_TRUE);CHKERRQ(ierr);
38*71438e86SJunchao Zhang       ierr   = PetscSFLinkGetMPIBuffersAndRequests(sf,link,direction,NULL,NULL,NULL,&reqs);CHKERRQ(ierr);
39*71438e86SJunchao Zhang     }
40*71438e86SJunchao Zhang     ierr = PetscSFLinkSyncStreamBeforeCallMPI(sf,link,direction);CHKERRQ(ierr);
41*71438e86SJunchao Zhang     ierr = MPI_Startall_isend(buflen,link->unit,nreqs,reqs);CHKERRMPI(ierr);
42*71438e86SJunchao Zhang   }
43*71438e86SJunchao Zhang   PetscFunctionReturn(0);
44*71438e86SJunchao Zhang }
45*71438e86SJunchao Zhang 
46*71438e86SJunchao Zhang static PetscErrorCode PetscSFLinkWaitRequests_MPI(PetscSF sf,PetscSFLink link,PetscSFDirection direction)
47*71438e86SJunchao Zhang {
48*71438e86SJunchao Zhang   PetscErrorCode       ierr;
49*71438e86SJunchao Zhang   PetscSF_Basic        *bas = (PetscSF_Basic*)sf->data;
50*71438e86SJunchao Zhang   const PetscMemType   rootmtype_mpi = link->rootmtype_mpi,leafmtype_mpi = link->leafmtype_mpi;
51*71438e86SJunchao Zhang   const PetscInt       rootdirect_mpi = link->rootdirect_mpi,leafdirect_mpi = link->leafdirect_mpi;
52*71438e86SJunchao Zhang 
53*71438e86SJunchao Zhang   PetscFunctionBegin;
54*71438e86SJunchao Zhang   ierr = MPI_Waitall(bas->nrootreqs,link->rootreqs[direction][rootmtype_mpi][rootdirect_mpi],MPI_STATUSES_IGNORE);CHKERRMPI(ierr);
55*71438e86SJunchao Zhang   ierr = MPI_Waitall(sf->nleafreqs, link->leafreqs[direction][leafmtype_mpi][leafdirect_mpi],MPI_STATUSES_IGNORE);CHKERRMPI(ierr);
56*71438e86SJunchao Zhang   if (direction == PETSCSF_ROOT2LEAF) {
57*71438e86SJunchao Zhang     ierr = PetscSFLinkCopyLeafBufferInCaseNotUseGpuAwareMPI(sf,link,PETSC_FALSE/* host2device after recving */);CHKERRQ(ierr);
58*71438e86SJunchao Zhang   } else {
59*71438e86SJunchao Zhang     ierr = PetscSFLinkCopyRootBufferInCaseNotUseGpuAwareMPI(sf,link,PETSC_FALSE);CHKERRQ(ierr);
60*71438e86SJunchao Zhang   }
61*71438e86SJunchao Zhang   PetscFunctionReturn(0);
62*71438e86SJunchao Zhang }
63*71438e86SJunchao Zhang 
64*71438e86SJunchao Zhang /*
65*71438e86SJunchao Zhang    The routine Creates a communication link for the given operation. It first looks up its link cache. If
66*71438e86SJunchao Zhang    there is a free & suitable one, it uses it. Otherwise it creates a new one.
67*71438e86SJunchao Zhang 
68*71438e86SJunchao Zhang    A link contains buffers and MPI requests for send/recv. It also contains pack/unpack routines to pack/unpack
69*71438e86SJunchao Zhang    root/leafdata to/from these buffers. Buffers are allocated at our discretion. When we find root/leafata
70*71438e86SJunchao Zhang    can be directly passed to MPI, we won't allocate them. Even we allocate buffers, we only allocate
71*71438e86SJunchao Zhang    those that are needed by the given `sfop` and `op`, in other words, we do lazy memory-allocation.
72*71438e86SJunchao Zhang 
73*71438e86SJunchao Zhang    The routine also allocates buffers on CPU when one does not use gpu-aware MPI but data is on GPU.
74*71438e86SJunchao Zhang 
75*71438e86SJunchao Zhang    In SFBasic, MPI requests are persistent. They are init'ed until we try to get requests from a link.
76*71438e86SJunchao Zhang 
77*71438e86SJunchao Zhang    The routine is shared by SFBasic and SFNeighbor based on the fact they all deal with sparse graphs and
78*71438e86SJunchao Zhang    need pack/unpack data.
79*71438e86SJunchao Zhang */
80*71438e86SJunchao Zhang 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)
81*71438e86SJunchao Zhang {
82*71438e86SJunchao Zhang   PetscErrorCode    ierr;
83*71438e86SJunchao Zhang   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
84*71438e86SJunchao Zhang   PetscInt          i,j,k,nrootreqs,nleafreqs,nreqs;
85*71438e86SJunchao Zhang   PetscSFLink       *p,link;
86*71438e86SJunchao Zhang   PetscSFDirection  direction;
87*71438e86SJunchao Zhang   MPI_Request       *reqs = NULL;
88*71438e86SJunchao Zhang   PetscBool         match,rootdirect[2],leafdirect[2];
89*71438e86SJunchao Zhang   PetscMemType      rootmtype = PetscMemTypeHost(xrootmtype) ? PETSC_MEMTYPE_HOST : PETSC_MEMTYPE_DEVICE; /* Convert to 0/1 as we will use it in subscript */
90*71438e86SJunchao Zhang   PetscMemType      leafmtype = PetscMemTypeHost(xleafmtype) ? PETSC_MEMTYPE_HOST : PETSC_MEMTYPE_DEVICE;
91*71438e86SJunchao Zhang   PetscMemType      rootmtype_mpi,leafmtype_mpi;   /* mtypes seen by MPI */
92*71438e86SJunchao Zhang   PetscInt          rootdirect_mpi,leafdirect_mpi; /* root/leafdirect seen by MPI*/
93*71438e86SJunchao Zhang 
94*71438e86SJunchao Zhang   PetscFunctionBegin;
95*71438e86SJunchao Zhang 
96*71438e86SJunchao Zhang   /* Can we directly use root/leafdirect with the given sf, sfop and op? */
97*71438e86SJunchao Zhang   for (i=PETSCSF_LOCAL; i<=PETSCSF_REMOTE; i++) {
98*71438e86SJunchao Zhang     if (sfop == PETSCSF_BCAST) {
99*71438e86SJunchao Zhang       rootdirect[i] = bas->rootcontig[i]; /* Pack roots */
100*71438e86SJunchao Zhang       leafdirect[i] = (sf->leafcontig[i] && op == MPI_REPLACE) ? PETSC_TRUE : PETSC_FALSE;  /* Unpack leaves */
101*71438e86SJunchao Zhang     } else if (sfop == PETSCSF_REDUCE) {
102*71438e86SJunchao Zhang       leafdirect[i] = sf->leafcontig[i];  /* Pack leaves */
103*71438e86SJunchao Zhang       rootdirect[i] = (bas->rootcontig[i] && op == MPI_REPLACE) ? PETSC_TRUE : PETSC_FALSE; /* Unpack roots */
104*71438e86SJunchao Zhang     } else { /* PETSCSF_FETCH */
105*71438e86SJunchao Zhang       rootdirect[i] = PETSC_FALSE; /* FETCH always need a separate rootbuf */
106*71438e86SJunchao Zhang       leafdirect[i] = PETSC_FALSE; /* We also force allocating a separate leafbuf so that leafdata and leafupdate can share mpi requests */
107*71438e86SJunchao Zhang     }
108*71438e86SJunchao Zhang   }
109*71438e86SJunchao Zhang 
110*71438e86SJunchao Zhang   if (sf->use_gpu_aware_mpi) {
111*71438e86SJunchao Zhang     rootmtype_mpi = rootmtype;
112*71438e86SJunchao Zhang     leafmtype_mpi = leafmtype;
113*71438e86SJunchao Zhang   } else {
114*71438e86SJunchao Zhang     rootmtype_mpi = leafmtype_mpi = PETSC_MEMTYPE_HOST;
115*71438e86SJunchao Zhang   }
116*71438e86SJunchao 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 */
117*71438e86SJunchao Zhang   rootdirect_mpi = rootdirect[PETSCSF_REMOTE] && (rootmtype_mpi == rootmtype)? 1 : 0;
118*71438e86SJunchao Zhang   leafdirect_mpi = leafdirect[PETSCSF_REMOTE] && (leafmtype_mpi == leafmtype)? 1 : 0;
119*71438e86SJunchao Zhang 
120*71438e86SJunchao Zhang   direction = (sfop == PETSCSF_BCAST)? PETSCSF_ROOT2LEAF : PETSCSF_LEAF2ROOT;
121*71438e86SJunchao Zhang   nrootreqs = bas->nrootreqs;
122*71438e86SJunchao Zhang   nleafreqs = sf->nleafreqs;
123*71438e86SJunchao Zhang 
124*71438e86SJunchao Zhang   /* Look for free links in cache */
125*71438e86SJunchao Zhang   for (p=&bas->avail; (link=*p); p=&link->next) {
126*71438e86SJunchao Zhang     if (!link->use_nvshmem) { /* Only check with MPI links */
127*71438e86SJunchao Zhang       ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr);
128*71438e86SJunchao Zhang       if (match) {
129*71438e86SJunchao Zhang         /* If root/leafdata will be directly passed to MPI, test if the data used to initialized the MPI requests matches with the current.
130*71438e86SJunchao Zhang            If not, free old requests. New requests will be lazily init'ed until one calls PetscSFLinkGetMPIBuffersAndRequests().
131*71438e86SJunchao Zhang         */
132*71438e86SJunchao Zhang         if (rootdirect_mpi && sf->persistent && link->rootreqsinited[direction][rootmtype][1] && link->rootdatadirect[direction][rootmtype] != rootdata) {
133*71438e86SJunchao Zhang           reqs = link->rootreqs[direction][rootmtype][1]; /* Here, rootmtype = rootmtype_mpi */
134*71438e86SJunchao Zhang           for (i=0; i<nrootreqs; i++) {if (reqs[i] != MPI_REQUEST_NULL) {ierr = MPI_Request_free(&reqs[i]);CHKERRMPI(ierr);}}
135*71438e86SJunchao Zhang           link->rootreqsinited[direction][rootmtype][1] = PETSC_FALSE;
136*71438e86SJunchao Zhang         }
137*71438e86SJunchao Zhang         if (leafdirect_mpi && sf->persistent && link->leafreqsinited[direction][leafmtype][1] && link->leafdatadirect[direction][leafmtype] != leafdata) {
138*71438e86SJunchao Zhang           reqs = link->leafreqs[direction][leafmtype][1];
139*71438e86SJunchao Zhang           for (i=0; i<nleafreqs; i++) {if (reqs[i] != MPI_REQUEST_NULL) {ierr = MPI_Request_free(&reqs[i]);CHKERRMPI(ierr);}}
140*71438e86SJunchao Zhang           link->leafreqsinited[direction][leafmtype][1] = PETSC_FALSE;
141*71438e86SJunchao Zhang         }
142*71438e86SJunchao Zhang         *p = link->next; /* Remove from available list */
143*71438e86SJunchao Zhang         goto found;
144*71438e86SJunchao Zhang       }
145*71438e86SJunchao Zhang     }
146*71438e86SJunchao Zhang   }
147*71438e86SJunchao Zhang 
148*71438e86SJunchao Zhang   ierr = PetscNew(&link);CHKERRQ(ierr);
149*71438e86SJunchao Zhang   ierr = PetscSFLinkSetUp_Host(sf,link,unit);CHKERRQ(ierr);
150*71438e86SJunchao Zhang   ierr = PetscCommGetNewTag(PetscObjectComm((PetscObject)sf),&link->tag);CHKERRQ(ierr); /* One tag per link */
151*71438e86SJunchao Zhang 
152*71438e86SJunchao Zhang   nreqs = (nrootreqs+nleafreqs)*8;
153*71438e86SJunchao Zhang   ierr  = PetscMalloc1(nreqs,&link->reqs);CHKERRQ(ierr);
154*71438e86SJunchao 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 */
155*71438e86SJunchao Zhang 
156*71438e86SJunchao Zhang   for (i=0; i<2; i++) { /* Two communication directions */
157*71438e86SJunchao Zhang     for (j=0; j<2; j++) { /* Two memory types */
158*71438e86SJunchao Zhang       for (k=0; k<2; k++) { /* root/leafdirect 0 or 1 */
159*71438e86SJunchao Zhang         link->rootreqs[i][j][k] = link->reqs + nrootreqs*(4*i+2*j+k);
160*71438e86SJunchao Zhang         link->leafreqs[i][j][k] = link->reqs + nrootreqs*8 + nleafreqs*(4*i+2*j+k);
161*71438e86SJunchao Zhang       }
162*71438e86SJunchao Zhang     }
163*71438e86SJunchao Zhang   }
164*71438e86SJunchao Zhang   link->StartCommunication    = PetscSFLinkStartRequests_MPI;
165*71438e86SJunchao Zhang   link->FinishCommunication   = PetscSFLinkWaitRequests_MPI;
166*71438e86SJunchao Zhang 
167*71438e86SJunchao Zhang found:
168*71438e86SJunchao Zhang 
169*71438e86SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
170*71438e86SJunchao Zhang   if ((PetscMemTypeDevice(xrootmtype) || PetscMemTypeDevice(xleafmtype)) && !link->deviceinited) {
171*71438e86SJunchao Zhang     #if defined(PETSC_HAVE_CUDA)
172*71438e86SJunchao Zhang       if (sf->backend == PETSCSF_BACKEND_CUDA)   {ierr = PetscSFLinkSetUp_CUDA(sf,link,unit);CHKERRQ(ierr);} /* Setup streams etc */
173*71438e86SJunchao Zhang     #endif
174*71438e86SJunchao Zhang     #if defined(PETSC_HAVE_HIP)
175*71438e86SJunchao Zhang       if (sf->backend == PETSCSF_BACKEND_HIP)    {ierr = PetscSFLinkSetUp_HIP(sf,link,unit);CHKERRQ(ierr);} /* Setup streams etc */
176*71438e86SJunchao Zhang     #endif
177*71438e86SJunchao Zhang     #if defined(PETSC_HAVE_KOKKOS)
178*71438e86SJunchao Zhang       if (sf->backend == PETSCSF_BACKEND_KOKKOS) {ierr = PetscSFLinkSetUp_Kokkos(sf,link,unit);CHKERRQ(ierr);}
179*71438e86SJunchao Zhang     #endif
180*71438e86SJunchao Zhang   }
181*71438e86SJunchao Zhang #endif
182*71438e86SJunchao Zhang 
183*71438e86SJunchao Zhang   /* Allocate buffers along root/leafdata */
184*71438e86SJunchao Zhang   for (i=PETSCSF_LOCAL; i<=PETSCSF_REMOTE; i++) {
185*71438e86SJunchao Zhang     /* For local communication, buffers are only needed when roots and leaves have different mtypes */
186*71438e86SJunchao Zhang     if (i == PETSCSF_LOCAL && rootmtype == leafmtype) continue;
187*71438e86SJunchao Zhang     if (bas->rootbuflen[i]) {
188*71438e86SJunchao Zhang       if (rootdirect[i]) { /* Aha, we disguise rootdata as rootbuf */
189*71438e86SJunchao Zhang         link->rootbuf[i][rootmtype] = (char*)rootdata + bas->rootstart[i]*link->unitbytes;
190*71438e86SJunchao Zhang       } else { /* Have to have a separate rootbuf */
191*71438e86SJunchao Zhang         if (!link->rootbuf_alloc[i][rootmtype]) {
192*71438e86SJunchao Zhang           ierr = PetscSFMalloc(sf,rootmtype,bas->rootbuflen[i]*link->unitbytes,(void**)&link->rootbuf_alloc[i][rootmtype]);CHKERRQ(ierr);
193*71438e86SJunchao Zhang         }
194*71438e86SJunchao Zhang         link->rootbuf[i][rootmtype] = link->rootbuf_alloc[i][rootmtype];
195*71438e86SJunchao Zhang       }
196*71438e86SJunchao Zhang     }
197*71438e86SJunchao Zhang 
198*71438e86SJunchao Zhang     if (sf->leafbuflen[i]) {
199*71438e86SJunchao Zhang       if (leafdirect[i]) {
200*71438e86SJunchao Zhang         link->leafbuf[i][leafmtype] = (char*)leafdata + sf->leafstart[i]*link->unitbytes;
201*71438e86SJunchao Zhang       } else {
202*71438e86SJunchao Zhang         if (!link->leafbuf_alloc[i][leafmtype]) {
203*71438e86SJunchao Zhang           ierr = PetscSFMalloc(sf,leafmtype,sf->leafbuflen[i]*link->unitbytes,(void**)&link->leafbuf_alloc[i][leafmtype]);CHKERRQ(ierr);
204*71438e86SJunchao Zhang         }
205*71438e86SJunchao Zhang         link->leafbuf[i][leafmtype] = link->leafbuf_alloc[i][leafmtype];
206*71438e86SJunchao Zhang       }
207*71438e86SJunchao Zhang     }
208*71438e86SJunchao Zhang   }
209*71438e86SJunchao Zhang 
210*71438e86SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
211*71438e86SJunchao Zhang   /* Allocate buffers on host for buffering data on device in cast not use_gpu_aware_mpi */
212*71438e86SJunchao Zhang   if (PetscMemTypeDevice(rootmtype) && PetscMemTypeHost(rootmtype_mpi)) {
213*71438e86SJunchao Zhang     if (!link->rootbuf_alloc[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST]) {
214*71438e86SJunchao Zhang       ierr = PetscMalloc(bas->rootbuflen[PETSCSF_REMOTE]*link->unitbytes,&link->rootbuf_alloc[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST]);CHKERRQ(ierr);
215*71438e86SJunchao Zhang     }
216*71438e86SJunchao Zhang     link->rootbuf[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST] = link->rootbuf_alloc[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST];
217*71438e86SJunchao Zhang   }
218*71438e86SJunchao Zhang   if (PetscMemTypeDevice(leafmtype) && PetscMemTypeHost(leafmtype_mpi)) {
219*71438e86SJunchao Zhang     if (!link->leafbuf_alloc[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST]) {
220*71438e86SJunchao Zhang       ierr = PetscMalloc(sf->leafbuflen[PETSCSF_REMOTE]*link->unitbytes,&link->leafbuf_alloc[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST]);CHKERRQ(ierr);
221*71438e86SJunchao Zhang     }
222*71438e86SJunchao Zhang     link->leafbuf[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST] = link->leafbuf_alloc[PETSCSF_REMOTE][PETSC_MEMTYPE_HOST];
223*71438e86SJunchao Zhang   }
224*71438e86SJunchao Zhang #endif
225*71438e86SJunchao Zhang 
226*71438e86SJunchao Zhang   /* Set `current` state of the link. They may change between different SF invocations with the same link */
227*71438e86SJunchao Zhang   if (sf->persistent) { /* If data is directly passed to MPI and inits MPI requests, record the data for comparison on future invocations */
228*71438e86SJunchao Zhang     if (rootdirect_mpi) link->rootdatadirect[direction][rootmtype] = rootdata;
229*71438e86SJunchao Zhang     if (leafdirect_mpi) link->leafdatadirect[direction][leafmtype] = leafdata;
230*71438e86SJunchao Zhang   }
231*71438e86SJunchao Zhang 
232*71438e86SJunchao Zhang   link->rootdata = rootdata; /* root/leafdata are keys to look up links in PetscSFXxxEnd */
233*71438e86SJunchao Zhang   link->leafdata = leafdata;
234*71438e86SJunchao Zhang   for (i=PETSCSF_LOCAL; i<=PETSCSF_REMOTE; i++) {
235*71438e86SJunchao Zhang     link->rootdirect[i] = rootdirect[i];
236*71438e86SJunchao Zhang     link->leafdirect[i] = leafdirect[i];
237*71438e86SJunchao Zhang   }
238*71438e86SJunchao Zhang   link->rootdirect_mpi  = rootdirect_mpi;
239*71438e86SJunchao Zhang   link->leafdirect_mpi  = leafdirect_mpi;
240*71438e86SJunchao Zhang   link->rootmtype       = rootmtype;
241*71438e86SJunchao Zhang   link->leafmtype       = leafmtype;
242*71438e86SJunchao Zhang   link->rootmtype_mpi   = rootmtype_mpi;
243*71438e86SJunchao Zhang   link->leafmtype_mpi   = leafmtype_mpi;
244*71438e86SJunchao Zhang 
245*71438e86SJunchao Zhang   link->next            = bas->inuse;
246*71438e86SJunchao Zhang   bas->inuse            = link;
247*71438e86SJunchao Zhang   *mylink               = link;
248*71438e86SJunchao Zhang   PetscFunctionReturn(0);
249*71438e86SJunchao Zhang }
250*71438e86SJunchao Zhang 
251