Actual source code: sfbasic.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2:  #include <../src/vec/is/sf/impls/basic/sfbasic.h>
  3:  #include <../src/vec/is/sf/impls/basic/sfpack.h>

  5: /*===================================================================================*/
  6: /*              SF public interface implementations                                  */
  7: /*===================================================================================*/
  8: PETSC_INTERN PetscErrorCode PetscSFSetUp_Basic(PetscSF sf)
  9: {
 11:   PetscSF_Basic  *bas = (PetscSF_Basic*)sf->data;
 12:   PetscInt       *rlengths,*ilengths,i;
 13:   PetscMPIInt    rank,niranks,*iranks,tag;
 14:   MPI_Comm       comm;
 15:   MPI_Group      group;
 16:   MPI_Request    *rootreqs,*leafreqs;

 19:   MPI_Comm_group(PETSC_COMM_SELF,&group);
 20:   PetscSFSetUpRanks(sf,group);
 21:   MPI_Group_free(&group);
 22:   PetscObjectGetComm((PetscObject)sf,&comm);
 23:   PetscObjectGetNewTag((PetscObject)sf,&tag);
 24:   MPI_Comm_rank(comm,&rank);
 25:   /*
 26:    * Inform roots about how many leaves and from which ranks
 27:    */
 28:   PetscMalloc1(sf->nranks,&rlengths);
 29:   /* Determine number, sending ranks and length of incoming */
 30:   for (i=0; i<sf->nranks; i++) {
 31:     rlengths[i] = sf->roffset[i+1] - sf->roffset[i]; /* Number of roots referenced by my leaves; for rank sf->ranks[i] */
 32:   }
 33:   PetscCommBuildTwoSided(comm,1,MPIU_INT,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,rlengths+sf->ndranks,&niranks,&iranks,(void**)&ilengths);

 35:   /* Sort iranks. See use of VecScatterGetRemoteOrdered_Private() in MatGetBrowsOfAoCols_MPIAIJ() on why.
 36:      We could sort ranks there at the price of allocating extra working arrays. Presumably, niranks is
 37:      small and the sorting is cheap.
 38:    */
 39:   PetscSortMPIIntWithIntArray(niranks,iranks,ilengths);

 41:   /* Partition into distinguished and non-distinguished incoming ranks */
 42:   bas->ndiranks = sf->ndranks;
 43:   bas->niranks = bas->ndiranks + niranks;
 44:   PetscMalloc2(bas->niranks,&bas->iranks,bas->niranks+1,&bas->ioffset);
 45:   bas->ioffset[0] = 0;
 46:   for (i=0; i<bas->ndiranks; i++) {
 47:     bas->iranks[i] = sf->ranks[i];
 48:     bas->ioffset[i+1] = bas->ioffset[i] + rlengths[i];
 49:   }
 50:   if (bas->ndiranks > 1 || (bas->ndiranks == 1 && bas->iranks[0] != rank)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Broken setup for shared ranks");
 51:   for ( ; i<bas->niranks; i++) {
 52:     bas->iranks[i] = iranks[i-bas->ndiranks];
 53:     bas->ioffset[i+1] = bas->ioffset[i] + ilengths[i-bas->ndiranks];
 54:   }
 55:   bas->itotal = bas->ioffset[i];
 56:   PetscFree(rlengths);
 57:   PetscFree(iranks);
 58:   PetscFree(ilengths);

 60:   /* Send leaf identities to roots */
 61:   PetscMalloc1(bas->itotal,&bas->irootloc);
 62:   PetscMalloc2(bas->niranks-bas->ndiranks,&rootreqs,sf->nranks-sf->ndranks,&leafreqs);
 63:   for (i=bas->ndiranks; i<bas->niranks; i++) {
 64:     MPI_Irecv(bas->irootloc+bas->ioffset[i],bas->ioffset[i+1]-bas->ioffset[i],MPIU_INT,bas->iranks[i],tag,comm,&rootreqs[i-bas->ndiranks]);
 65:   }
 66:   for (i=0; i<sf->nranks; i++) {
 67:     PetscMPIInt npoints;
 68:     PetscMPIIntCast(sf->roffset[i+1] - sf->roffset[i],&npoints);
 69:     if (i < sf->ndranks) {
 70:       if (sf->ranks[i] != rank) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot interpret distinguished leaf rank");
 71:       if (bas->iranks[0] != rank) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot interpret distinguished root rank");
 72:       if (npoints != bas->ioffset[1]-bas->ioffset[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Distinguished rank exchange has mismatched lengths");
 73:       PetscArraycpy(bas->irootloc+bas->ioffset[0],sf->rremote+sf->roffset[i],npoints);
 74:       continue;
 75:     }
 76:     MPI_Isend(sf->rremote+sf->roffset[i],npoints,MPIU_INT,sf->ranks[i],tag,comm,&leafreqs[i-sf->ndranks]);
 77:   }
 78:   MPI_Waitall(bas->niranks-bas->ndiranks,rootreqs,MPI_STATUSES_IGNORE);
 79:   MPI_Waitall(sf->nranks-sf->ndranks,leafreqs,MPI_STATUSES_IGNORE);
 80:   PetscFree2(rootreqs,leafreqs);

 82:   sf->nleafreqs  = sf->nranks - sf->ndranks;
 83:   bas->nrootreqs = bas->niranks - bas->ndiranks;
 84:   sf->persistent = PETSC_TRUE;

 86:   /* Setup fields related to packing */
 87:   PetscSFSetUpPackFields(sf);
 88:   return(0);
 89: }

 91: PETSC_INTERN PetscErrorCode PetscSFReset_Basic(PetscSF sf)
 92: {
 93:   PetscErrorCode    ierr;
 94:   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;

 97:   if (bas->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed");
 98:   PetscFree2(bas->iranks,bas->ioffset);
 99:   PetscFree(bas->irootloc);
100: #if defined(PETSC_HAVE_CUDA)
101:   {
102:   PetscInt  i;
103:   for (i=0; i<2; i++) {if (bas->irootloc_d[i]) {cudaError_t err = cudaFree(bas->irootloc_d[i]);CHKERRCUDA(err);bas->irootloc_d[i]=NULL;}}
104:   }
105: #endif
106:   PetscSFLinkDestroy(sf,&bas->avail);
107:   PetscSFResetPackFields(sf);
108:   return(0);
109: }

111: PETSC_INTERN PetscErrorCode PetscSFDestroy_Basic(PetscSF sf)
112: {

116:   PetscSFReset_Basic(sf);
117:   PetscFree(sf->data);
118:   return(0);
119: }

121: PETSC_INTERN PetscErrorCode PetscSFView_Basic(PetscSF sf,PetscViewer viewer)
122: {
124:   PetscBool      iascii;

127:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
128:   if (iascii) {PetscViewerASCIIPrintf(viewer,"  sort=%s\n",sf->rankorder ? "rank-order" : "unordered");}
129:   return(0);
130: }

132: static PetscErrorCode PetscSFBcastAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,const void *rootdata,PetscMemType leafmtype,void *leafdata,MPI_Op op)
133: {
134:   PetscErrorCode    ierr;
135:   PetscSFLink       link = NULL;
136:   MPI_Request       *rootreqs = NULL,*leafreqs = NULL;
137:   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;

140:   /* Create a communication link, which provides buffers & MPI requests etc */
141:   PetscSFLinkCreate(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op,PETSCSF_BCAST,&link);
142:   /* Get MPI requests from the link. We do not need buffers explicitly since we use persistent MPI */
143:   PetscSFLinkGetMPIBuffersAndRequests(sf,link,PETSCSF_../../../../../..2LEAF,NULL,NULL,&rootreqs,&leafreqs);
144:   /* Post Irecv for remote */
145:   MPI_Startall_irecv(sf->leafbuflen[PETSCSF_REMOTE],unit,sf->nleafreqs,leafreqs);
146:   /* Pack rootdata and do Isend for remote */
147:   PetscSFLinkPackRootData(sf,link,PETSCSF_REMOTE,rootdata);
148:   MPI_Startall_isend(bas->rootbuflen[PETSCSF_REMOTE],unit,bas->nrootreqs,rootreqs);
149:   /* Do local BcastAndOp, which overlaps with the irecv/isend above */
150:   PetscSFLinkBcastAndOpLocal(sf,link,rootdata,leafdata,op);
151:   return(0);
152: }

154: PETSC_INTERN PetscErrorCode PetscSFBcastAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
155: {
156:   PetscErrorCode    ierr;
157:   PetscSFLink       link = NULL;

160:   /* Retrieve the link used in XxxBegin() with root/leafdata as key */
161:   PetscSFLinkGetInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,&link);
162:   /* Wait for the completion of mpi */
163:   PetscSFLinkMPIWaitall(sf,link,PETSCSF_../../../../../..2LEAF);
164:   /* Unpack leafdata and reclaim the link */
165:   PetscSFLinkUnpackLeafData(sf,link,PETSCSF_REMOTE,leafdata,op);
166:   PetscSFLinkReclaim(sf,&link);
167:   return(0);
168: }

170: /* Shared by ReduceBegin and FetchAndOpBegin */
171: PETSC_STATIC_INLINE PetscErrorCode PetscSFLeafToRootBegin_Basic(PetscSF sf,MPI_Datatype unit,PetscMemType leafmtype,const void *leafdata,PetscMemType rootmtype,void *rootdata,MPI_Op op,PetscSFOperation sfop,PetscSFLink *out)
172: {
173:   PetscErrorCode    ierr;
174:   PetscSFLink       link;
175:   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;
176:   MPI_Request       *rootreqs = NULL,*leafreqs = NULL;

179:   PetscSFLinkCreate(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op,sfop,&link);
180:   PetscSFLinkGetMPIBuffersAndRequests(sf,link,PETSCSF_LEAF2../../../../../..,NULL,NULL,&rootreqs,&leafreqs);
181:   MPI_Startall_irecv(bas->rootbuflen[PETSCSF_REMOTE],unit,bas->nrootreqs,rootreqs);
182:   PetscSFLinkPackLeafData(sf,link,PETSCSF_REMOTE,leafdata);
183:   MPI_Startall_isend(sf->leafbuflen[PETSCSF_REMOTE],unit,sf->nleafreqs,leafreqs);
184:   *out = link;
185:   return(0);
186: }

188: /* leaf -> root with reduction */
189: static PetscErrorCode PetscSFReduceBegin_Basic(PetscSF sf,MPI_Datatype unit,PetscMemType leafmtype,const void *leafdata,PetscMemType rootmtype,void *rootdata,MPI_Op op)
190: {
191:   PetscErrorCode    ierr;
192:   PetscSFLink       link = NULL;

195:   PetscSFLeafToRootBegin_Basic(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op,PETSCSF_REDUCE,&link);
196:   PetscSFLinkReduceLocal(sf,link,leafdata,rootdata,op);
197:   return(0);
198: }

200: PETSC_INTERN PetscErrorCode PetscSFReduceEnd_Basic(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
201: {
202:   PetscErrorCode    ierr;
203:   PetscSFLink       link = NULL;

206:   PetscSFLinkGetInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,&link);
207:   PetscSFLinkMPIWaitall(sf,link,PETSCSF_LEAF2../../../../../..);
208:   PetscSFLinkUnpackRootData(sf,link,PETSCSF_REMOTE,rootdata,op);
209:   PetscSFLinkReclaim(sf,&link);
210:   return(0);
211: }

213: PETSC_INTERN PetscErrorCode PetscSFFetchAndOpBegin_Basic(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,void *rootdata,PetscMemType leafmtype,const void *leafdata,void *leafupdate,MPI_Op op)
214: {
215:   PetscErrorCode    ierr;
216:   PetscSFLink       link = NULL;

219:   PetscSFLeafToRootBegin_Basic(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op,PETSCSF_FETCH,&link);
220:   PetscSFLinkFetchAndOpLocal(sf,link,rootdata,leafdata,leafupdate,op);
221:   return(0);
222: }

224: static PetscErrorCode PetscSFFetchAndOpEnd_Basic(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
225: {
226:   PetscErrorCode    ierr;
227:   PetscSFLink       link = NULL;
228:   MPI_Request       *rootreqs = NULL,*leafreqs = NULL;
229:   PetscSF_Basic     *bas = (PetscSF_Basic*)sf->data;

232:   PetscSFLinkGetInUse(sf,unit,rootdata,leafdata,PETSC_OWN_POINTER,&link);
233:   /* This implementation could be changed to unpack as receives arrive, at the cost of non-determinism */
234:   PetscSFLinkMPIWaitall(sf,link,PETSCSF_LEAF2../../../../../..);
235:   /* Do fetch-and-op, the (remote) update results are in rootbuf */
236:   PetscSFLinkFetchRootData(sf,link,PETSCSF_REMOTE,rootdata,op);

238:   /* Bcast rootbuf to leafupdate */
239:   PetscSFLinkGetMPIBuffersAndRequests(sf,link,PETSCSF_../../../../../..2LEAF,NULL,NULL,&rootreqs,&leafreqs);
240:   /* Post leaf receives and root sends */
241:   MPI_Startall_irecv(sf->leafbuflen[PETSCSF_REMOTE],unit,sf->nleafreqs,leafreqs);
242:   MPI_Startall_isend(bas->rootbuflen[PETSCSF_REMOTE],unit,bas->nrootreqs,rootreqs);
243:   /* Unpack and insert fetched data into leaves */
244:   PetscSFLinkMPIWaitall(sf,link,PETSCSF_../../../../../..2LEAF);
245:   PetscSFLinkUnpackLeafData(sf,link,PETSCSF_REMOTE,leafupdate,MPIU_REPLACE);
246:   PetscSFLinkReclaim(sf,&link);
247:   return(0);
248: }

250: PETSC_INTERN PetscErrorCode PetscSFGetLeafRanks_Basic(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
251: {
252:   PetscSF_Basic *bas = (PetscSF_Basic*)sf->data;

255:   if (niranks)  *niranks  = bas->niranks;
256:   if (iranks)   *iranks   = bas->iranks;
257:   if (ioffset)  *ioffset  = bas->ioffset;
258:   if (irootloc) *irootloc = bas->irootloc;
259:   return(0);
260: }

262: /* An optimized PetscSFCreateEmbeddedSF. We aggresively make use of the established communication on sf.
263:    We need one bcast on sf, and no communication anymore to build the embedded sf. Note that selected[]
264:    was sorted before calling the routine.
265:  */
266: PETSC_INTERN PetscErrorCode PetscSFCreateEmbeddedSF_Basic(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
267: {
268:   PetscSF           esf;
269:   PetscInt          esf_nranks,esf_ndranks,*esf_roffset,*esf_rmine,*esf_rremote;
270:   PetscInt          i,j,p,q,nroots,esf_nleaves,*new_ilocal,nranks,ndranks,niranks,ndiranks,minleaf,maxleaf,maxlocal;
271:   char              *rootdata,*leafdata,*leafmem; /* Only stores 0 or 1, so we can save memory with char */
272:   PetscMPIInt       *esf_ranks;
273:   const PetscMPIInt *ranks,*iranks;
274:   const PetscInt    *roffset,*rmine,*rremote,*ioffset,*irootloc;
275:   PetscBool         connected;
276:   PetscSFNode       *new_iremote;
277:   PetscSF_Basic     *bas;
278:   PetscErrorCode    ierr;

281:   PetscSFCreate(PetscObjectComm((PetscObject)sf),&esf);
282:   PetscSFSetType(esf,PETSCSFBASIC); /* This optimized routine can only create a basic sf */

284:   /* Find out which leaves are still connected to roots in the embedded sf by doing a Bcast */
285:   PetscSFGetGraph(sf,&nroots,NULL,NULL,NULL);
286:   PetscSFGetLeafRange(sf,&minleaf,&maxleaf);
287:   maxlocal = maxleaf - minleaf + 1;
288:   PetscCalloc2(nroots,&rootdata,maxlocal,&leafmem);
289:   leafdata = leafmem - minleaf;
290:   /* Tag selected roots */
291:   for (i=0; i<nselected; ++i) rootdata[selected[i]] = 1;

293:   PetscSFBcastBegin(sf,MPI_CHAR,rootdata,leafdata);
294:   PetscSFBcastEnd(sf,MPI_CHAR,rootdata,leafdata);
295:   PetscSFGetLeafInfo_Basic(sf,&nranks,&ndranks,&ranks,&roffset,&rmine,&rremote); /* Get send info */
296:   esf_nranks = esf_ndranks = esf_nleaves = 0;
297:   for (i=0; i<nranks; i++) {
298:     connected = PETSC_FALSE; /* Is this process still connected to this remote root rank? */
299:     for (j=roffset[i]; j<roffset[i+1]; j++) {if (leafdata[rmine[j]]) {esf_nleaves++; connected = PETSC_TRUE;}}
300:     if (connected) {esf_nranks++; if (i < ndranks) esf_ndranks++;}
301:   }

303:   /* Set graph of esf and also set up its outgoing communication (i.e., send info), which is usually done by PetscSFSetUpRanks */
304:   PetscMalloc1(esf_nleaves,&new_ilocal);
305:   PetscMalloc1(esf_nleaves,&new_iremote);
306:   PetscMalloc4(esf_nranks,&esf_ranks,esf_nranks+1,&esf_roffset,esf_nleaves,&esf_rmine,esf_nleaves,&esf_rremote);
307:   p    = 0; /* Counter for connected root ranks */
308:   q    = 0; /* Counter for connected leaves */
309:   esf_roffset[0] = 0;
310:   for (i=0; i<nranks; i++) { /* Scan leaf data again to fill esf arrays */
311:     connected = PETSC_FALSE;
312:     for (j=roffset[i]; j<roffset[i+1]; j++) {
313:       if (leafdata[rmine[j]]) {
314:         esf_rmine[q]         = new_ilocal[q] = rmine[j];
315:         esf_rremote[q]       = rremote[j];
316:         new_iremote[q].index = rremote[j];
317:         new_iremote[q].rank  = ranks[i];
318:         connected            = PETSC_TRUE;
319:         q++;
320:       }
321:     }
322:     if (connected) {
323:       esf_ranks[p]     = ranks[i];
324:       esf_roffset[p+1] = q;
325:       p++;
326:     }
327:   }

329:   /* SetGraph internally resets the SF, so we only set its fields after the call */
330:   PetscSFSetGraph(esf,nroots,esf_nleaves,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);
331:   esf->nranks    = esf_nranks;
332:   esf->ndranks   = esf_ndranks;
333:   esf->ranks     = esf_ranks;
334:   esf->roffset   = esf_roffset;
335:   esf->rmine     = esf_rmine;
336:   esf->rremote   = esf_rremote;
337:   esf->nleafreqs = esf_nranks - esf_ndranks;

339:   /* Set up the incoming communication (i.e., recv info) stored in esf->data, which is usually done by PetscSFSetUp_Basic */
340:   bas  = (PetscSF_Basic*)esf->data;
341:   PetscSFGetRootInfo_Basic(sf,&niranks,&ndiranks,&iranks,&ioffset,&irootloc); /* Get recv info */
342:   /* Embedded sf always has simpler communication than the original one. We might allocate longer arrays than needed here. But we
343:      we do not care since these arrays are usually short. The benefit is we can fill these arrays by just parsing irootloc once.
344:    */
345:   PetscMalloc2(niranks,&bas->iranks,niranks+1,&bas->ioffset);
346:   PetscMalloc1(ioffset[niranks],&bas->irootloc);
347:   bas->niranks = bas->ndiranks = bas->ioffset[0] = 0;
348:   p = 0; /* Counter for connected leaf ranks */
349:   q = 0; /* Counter for connected roots */
350:   for (i=0; i<niranks; i++) {
351:     connected = PETSC_FALSE; /* Is the current process still connected to this remote leaf rank? */
352:     for (j=ioffset[i]; j<ioffset[i+1]; j++) {
353:       if (rootdata[irootloc[j]]) {
354:         bas->irootloc[q++] = irootloc[j];
355:         connected = PETSC_TRUE;
356:       }
357:     }
358:     if (connected) {
359:       bas->niranks++;
360:       if (i<ndiranks) bas->ndiranks++; /* Note that order of ranks (including distinguished ranks) is kept */
361:       bas->iranks[p]    = iranks[i];
362:       bas->ioffset[p+1] = q;
363:       p++;
364:     }
365:   }
366:   bas->itotal     = q;
367:   bas->nrootreqs  = bas->niranks - bas->ndiranks;
368:   esf->persistent = PETSC_TRUE;
369:   /* Setup packing related fields */
370:   PetscSFSetUpPackFields(esf);

372:   esf->setupcalled = PETSC_TRUE; /* We have done setup ourselves! */
373:   PetscFree2(rootdata,leafmem);
374:   *newsf = esf;
375:   return(0);
376: }

378: PETSC_EXTERN PetscErrorCode PetscSFCreate_Basic(PetscSF sf)
379: {
380:   PetscSF_Basic  *dat;

384:   sf->ops->SetUp                = PetscSFSetUp_Basic;
385:   sf->ops->Reset                = PetscSFReset_Basic;
386:   sf->ops->Destroy              = PetscSFDestroy_Basic;
387:   sf->ops->View                 = PetscSFView_Basic;
388:   sf->ops->BcastAndOpBegin      = PetscSFBcastAndOpBegin_Basic;
389:   sf->ops->BcastAndOpEnd        = PetscSFBcastAndOpEnd_Basic;
390:   sf->ops->ReduceBegin          = PetscSFReduceBegin_Basic;
391:   sf->ops->ReduceEnd            = PetscSFReduceEnd_Basic;
392:   sf->ops->FetchAndOpBegin      = PetscSFFetchAndOpBegin_Basic;
393:   sf->ops->FetchAndOpEnd        = PetscSFFetchAndOpEnd_Basic;
394:   sf->ops->GetLeafRanks         = PetscSFGetLeafRanks_Basic;
395:   sf->ops->CreateEmbeddedSF     = PetscSFCreateEmbeddedSF_Basic;

397:   PetscNewLog(sf,&dat);
398:   sf->data = (void*)dat;
399:   return(0);
400: }