Actual source code: sf.c

petsc-master 2019-11-16
Report Typos and Errors
  1:  #include <petsc/private/sfimpl.h>
  2:  #include <petsc/private/hashseti.h>
  3:  #include <petscctable.h>

  5: #if defined(PETSC_USE_DEBUG)
  6: #  define PetscSFCheckGraphSet(sf,arg) do {                          \
  7:     if (PetscUnlikely(!(sf)->graphset))                              \
  8:       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() or PetscSFSetGraphWithPattern() on argument %D \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME); \
  9:   } while (0)
 10: #else
 11: #  define PetscSFCheckGraphSet(sf,arg) do {} while (0)
 12: #endif

 14: const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",0};

 16: /*@
 17:    PetscSFCreate - create a star forest communication context

 19:    Collective

 21:    Input Arguments:
 22: .  comm - communicator on which the star forest will operate

 24:    Output Arguments:
 25: .  sf - new star forest context

 27:    Options Database Keys:
 28: +  -sf_type basic     -Use MPI persistent Isend/Irecv for communication (Default)
 29: .  -sf_type window    -Use MPI-3 one-sided window for communication
 30: -  -sf_type neighbor  -Use MPI-3 neighborhood collectives for communication

 32:    Level: intermediate

 34:    Notes:
 35:    When one knows the communication graph is one of the predefined graph, such as MPI_Alltoall, MPI_Allgatherv,
 36:    MPI_Gatherv, one can create a PetscSF and then set its graph with PetscSFSetGraphWithPattern(). These special
 37:    SFs are optimized and they have better performance than general SFs.

 39: .seealso: PetscSFSetGraph(), PetscSFSetGraphWithPattern(), PetscSFDestroy()
 40: @*/
 41: PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf)
 42: {
 44:   PetscSF        b;

 48:   PetscSFInitializePackage();

 50:   PetscHeaderCreate(b,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);

 52:   b->nroots    = -1;
 53:   b->nleaves   = -1;
 54:   b->minleaf   = PETSC_MAX_INT;
 55:   b->maxleaf   = PETSC_MIN_INT;
 56:   b->nranks    = -1;
 57:   b->rankorder = PETSC_TRUE;
 58:   b->ingroup   = MPI_GROUP_NULL;
 59:   b->outgroup  = MPI_GROUP_NULL;
 60:   b->graphset  = PETSC_FALSE;

 62:   *sf = b;
 63:   return(0);
 64: }

 66: /*@
 67:    PetscSFReset - Reset a star forest so that different sizes or neighbors can be used

 69:    Collective

 71:    Input Arguments:
 72: .  sf - star forest

 74:    Level: advanced

 76: .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy()
 77: @*/
 78: PetscErrorCode PetscSFReset(PetscSF sf)
 79: {

 84:   if (sf->ops->Reset) {(*sf->ops->Reset)(sf);}
 85:   sf->nroots   = -1;
 86:   sf->nleaves  = -1;
 87:   sf->minleaf  = PETSC_MAX_INT;
 88:   sf->maxleaf  = PETSC_MIN_INT;
 89:   sf->mine     = NULL;
 90:   sf->remote   = NULL;
 91:   sf->graphset = PETSC_FALSE;
 92:   PetscFree(sf->mine_alloc);
 93:   PetscFree(sf->remote_alloc);
 94:   sf->nranks = -1;
 95:   PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);
 96: #if defined(PETSC_HAVE_CUDA)
 97:   if (sf->rmine_d) {cudaError_t err = cudaFree(sf->rmine_d);CHKERRCUDA(err);sf->rmine_d=NULL;}
 98: #endif
 99:   sf->degreeknown = PETSC_FALSE;
100:   PetscFree(sf->degree);
101:   if (sf->ingroup  != MPI_GROUP_NULL) {MPI_Group_free(&sf->ingroup);}
102:   if (sf->outgroup != MPI_GROUP_NULL) {MPI_Group_free(&sf->outgroup);}
103:   PetscSFDestroy(&sf->multi);
104:   PetscLayoutDestroy(&sf->map);
105:   sf->setupcalled = PETSC_FALSE;
106:   return(0);
107: }

109: /*@C
110:    PetscSFSetType - Set the PetscSF communication implementation

112:    Collective on PetscSF

114:    Input Parameters:
115: +  sf - the PetscSF context
116: -  type - a known method

118:    Options Database Key:
119: .  -sf_type <type> - Sets the method; use -help for a list
120:    of available methods (for instance, window, pt2pt, neighbor)

122:    Notes:
123:    See "include/petscsf.h" for available methods (for instance)
124: +    PETSCSFWINDOW - MPI-2/3 one-sided
125: -    PETSCSFBASIC - basic implementation using MPI-1 two-sided

127:   Level: intermediate

129: .seealso: PetscSFType, PetscSFCreate()
130: @*/
131: PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type)
132: {
133:   PetscErrorCode ierr,(*r)(PetscSF);
134:   PetscBool      match;


140:   PetscObjectTypeCompare((PetscObject)sf,type,&match);
141:   if (match) return(0);

143:   PetscFunctionListFind(PetscSFList,type,&r);
144:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type);
145:   /* Destroy the previous PetscSF implementation context */
146:   if (sf->ops->Destroy) {(*(sf)->ops->Destroy)(sf);}
147:   PetscMemzero(sf->ops,sizeof(*sf->ops));
148:   PetscObjectChangeTypeName((PetscObject)sf,type);
149:   (*r)(sf);
150:   return(0);
151: }

153: /*@C
154:   PetscSFGetType - Get the PetscSF communication implementation

156:   Not Collective

158:   Input Parameter:
159: . sf  - the PetscSF context

161:   Output Parameter:
162: . type - the PetscSF type name

164:   Level: intermediate

166: .seealso: PetscSFSetType(), PetscSFCreate()
167: @*/
168: PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type)
169: {
173:   *type = ((PetscObject)sf)->type_name;
174:   return(0);
175: }

177: /*@
178:    PetscSFDestroy - destroy star forest

180:    Collective

182:    Input Arguments:
183: .  sf - address of star forest

185:    Level: intermediate

187: .seealso: PetscSFCreate(), PetscSFReset()
188: @*/
189: PetscErrorCode PetscSFDestroy(PetscSF *sf)
190: {

194:   if (!*sf) return(0);
196:   if (--((PetscObject)(*sf))->refct > 0) {*sf = NULL; return(0);}
197:   PetscSFReset(*sf);
198:   if ((*sf)->ops->Destroy) {(*(*sf)->ops->Destroy)(*sf);}
199:   PetscHeaderDestroy(sf);
200:   return(0);
201: }

203: static PetscErrorCode PetscSFCheckGraphValid_Private(PetscSF sf)
204: {
205: #if defined(PETSC_USE_DEBUG)
206:   PetscInt           i, nleaves;
207:   PetscMPIInt        size;
208:   const PetscInt    *ilocal;
209:   const PetscSFNode *iremote;
210:   PetscErrorCode     ierr;

213:   if (!sf->graphset) return(0);
214:   PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,&iremote);
215:   MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);
216:   for (i = 0; i < nleaves; i++) {
217:     const PetscInt rank = iremote[i].rank;
218:     const PetscInt remote = iremote[i].index;
219:     const PetscInt leaf = ilocal ? ilocal[i] : i;
220:     if (rank < 0 || rank >= size) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided rank (%D) for remote %D is invalid, should be in [0, %d)",rank,i,size);
221:     if (remote < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided index (%D) for remote %D is invalid, should be >= 0",remote,i);
222:     if (leaf < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided location (%D) for leaf %D is invalid, should be >= 0",leaf,i);
223:   }
224:   return(0);
225: #else
227:   return(0);
228: #endif
229: }

231: /*@
232:    PetscSFSetUp - set up communication structures

234:    Collective

236:    Input Arguments:
237: .  sf - star forest communication object

239:    Level: beginner

241: .seealso: PetscSFSetFromOptions(), PetscSFSetType()
242: @*/
243: PetscErrorCode PetscSFSetUp(PetscSF sf)
244: {

249:   PetscSFCheckGraphSet(sf,1);
250:   if (sf->setupcalled) return(0);
251:   PetscSFCheckGraphValid_Private(sf);
252:   if (!((PetscObject)sf)->type_name) {PetscSFSetType(sf,PETSCSFBASIC);}
253:   PetscLogEventBegin(PETSCSF_SetUp,sf,0,0,0);
254:   if (sf->ops->SetUp) {(*sf->ops->SetUp)(sf);}
255:   PetscLogEventEnd(PETSCSF_SetUp,sf,0,0,0);
256:   sf->setupcalled = PETSC_TRUE;
257:   return(0);
258: }

260: /*@
261:    PetscSFSetFromOptions - set PetscSF options using the options database

263:    Logically Collective

265:    Input Arguments:
266: .  sf - star forest

268:    Options Database Keys:
269: +  -sf_type - implementation type, see PetscSFSetType()
270: -  -sf_rank_order - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise

272:    Level: intermediate

274: .seealso: PetscSFWindowSetSyncType()
275: @*/
276: PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
277: {
278:   PetscSFType    deft;
279:   char           type[256];
281:   PetscBool      flg;

285:   PetscObjectOptionsBegin((PetscObject)sf);
286:   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
287:   PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,sizeof(type),&flg);
288:   PetscSFSetType(sf,flg ? type : deft);
289:   PetscOptionsBool("-sf_rank_order","sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise","PetscSFSetRankOrder",sf->rankorder,&sf->rankorder,NULL);
290:   if (sf->ops->SetFromOptions) {(*sf->ops->SetFromOptions)(PetscOptionsObject,sf);}
291:   PetscOptionsEnd();
292:   return(0);
293: }

295: /*@
296:    PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order

298:    Logically Collective

300:    Input Arguments:
301: +  sf - star forest
302: -  flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic)

304:    Level: advanced

306: .seealso: PetscSFGatherBegin(), PetscSFScatterBegin()
307: @*/
308: PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg)
309: {

314:   if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
315:   sf->rankorder = flg;
316:   return(0);
317: }

319: /*@
320:    PetscSFSetGraph - Set a parallel star forest

322:    Collective

324:    Input Arguments:
325: +  sf - star forest
326: .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
327: .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
328: .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage (locations must be >= 0, enforced
329: during setup in debug mode)
330: .  localmode - copy mode for ilocal
331: .  iremote - remote locations of root vertices for each leaf on the current process (locations must be >= 0, enforced
332: during setup in debug mode)
333: -  remotemode - copy mode for iremote

335:    Level: intermediate

337:    Notes:
338:     In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode

340:    Developers Note: Local indices which are the identity permutation in the range [0,nleaves) are discarded as they
341:    encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not
342:    needed

344:    Developers Note: This object does not necessarily encode a true star forest in the graph theoretic sense, since leaf
345:    indices are not required to be unique. Some functions, however, rely on unique leaf indices (checked in debug mode).

347: .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph()
348: @*/
349: PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode)
350: {

357:   if (nroots  < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nroots %D, cannot be negative",nroots);
358:   if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves);

360:   PetscSFReset(sf);
361:   PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);

363:   sf->nroots  = nroots;
364:   sf->nleaves = nleaves;

366:   if (nleaves && ilocal) {
367:     PetscInt i;
368:     PetscInt minleaf = PETSC_MAX_INT;
369:     PetscInt maxleaf = PETSC_MIN_INT;
370:     int      contiguous = 1;
371:     for (i=0; i<nleaves; i++) {
372:       minleaf = PetscMin(minleaf,ilocal[i]);
373:       maxleaf = PetscMax(maxleaf,ilocal[i]);
374:       contiguous &= (ilocal[i] == i);
375:     }
376:     sf->minleaf = minleaf;
377:     sf->maxleaf = maxleaf;
378:     if (contiguous) {
379:       if (localmode == PETSC_OWN_POINTER) {
380:         PetscFree(ilocal);
381:       }
382:       ilocal = NULL;
383:     }
384:   } else {
385:     sf->minleaf = 0;
386:     sf->maxleaf = nleaves - 1;
387:   }

389:   if (ilocal) {
390:     switch (localmode) {
391:     case PETSC_COPY_VALUES:
392:       PetscMalloc1(nleaves,&sf->mine_alloc);
393:       PetscArraycpy(sf->mine_alloc,ilocal,nleaves);
394:       sf->mine = sf->mine_alloc;
395:       break;
396:     case PETSC_OWN_POINTER:
397:       sf->mine_alloc = (PetscInt*)ilocal;
398:       sf->mine       = sf->mine_alloc;
399:       break;
400:     case PETSC_USE_POINTER:
401:       sf->mine_alloc = NULL;
402:       sf->mine       = (PetscInt*)ilocal;
403:       break;
404:     default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode");
405:     }
406:   }

408:   switch (remotemode) {
409:   case PETSC_COPY_VALUES:
410:     PetscMalloc1(nleaves,&sf->remote_alloc);
411:     PetscArraycpy(sf->remote_alloc,iremote,nleaves);
412:     sf->remote = sf->remote_alloc;
413:     break;
414:   case PETSC_OWN_POINTER:
415:     sf->remote_alloc = (PetscSFNode*)iremote;
416:     sf->remote       = sf->remote_alloc;
417:     break;
418:   case PETSC_USE_POINTER:
419:     sf->remote_alloc = NULL;
420:     sf->remote       = (PetscSFNode*)iremote;
421:     break;
422:   default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode");
423:   }

425:   PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);
426:   sf->graphset = PETSC_TRUE;
427:   return(0);
428: }

430: /*@
431:   PetscSFSetGraphWithPattern - Sets the graph of an SF with a specific pattern

433:   Collective

435:   Input Parameters:
436: + sf      - The PetscSF
437: . map     - Layout of roots over all processes (insignificant when pattern is PETSCSF_PATTERN_ALLTOALL)
438: - pattern - One of PETSCSF_PATTERN_ALLGATHER, PETSCSF_PATTERN_GATHER, PETSCSF_PATTERN_ALLTOALL

440:   Notes:
441:   It is easier to explain PetscSFPattern using vectors. Suppose we have an MPI vector x and its layout is map.
442:   n and N are local and global sizes of x respectively.

444:   With PETSCSF_PATTERN_ALLGATHER, the routine creates a graph that if one does Bcast on it, it will copy x to
445:   sequential vectors y on all ranks.

447:   With PETSCSF_PATTERN_GATHER, the routine creates a graph that if one does Bcast on it, it will copy x to a
448:   sequential vector y on rank 0.

450:   In above cases, entries of x are roots and entries of y are leaves.

452:   With PETSCSF_PATTERN_ALLTOALL, map is insignificant. Suppose NP is size of sf's communicator. The routine
453:   creates a graph that every rank has NP leaves and NP roots. On rank i, its leaf j is connected to root i
454:   of rank j. Here 0 <=i,j<NP. It is a kind of MPI_Alltoall with sendcount/recvcount being 1. Note that it does
455:   not mean one can not send multiple items. One just needs to create a new MPI datatype for the mulptiple data
456:   items with MPI_Type_contiguous() and use that as the <unit> argument in SF routines.

458:   In this case, roots and leaves are symmetric.

460:   Level: intermediate
461:  @*/
462: PetscErrorCode PetscSFSetGraphWithPattern(PetscSF sf,PetscLayout map,PetscSFPattern pattern)
463: {
464:   MPI_Comm       comm;
465:   PetscInt       n,N,res[2];
466:   PetscMPIInt    rank,size;
467:   PetscSFType    type;

471:   PetscObjectGetComm((PetscObject)sf, &comm);
472:   if (pattern < PETSCSF_PATTERN_ALLGATHER || pattern > PETSCSF_PATTERN_ALLTOALL) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Unsupported PetscSFPattern %D\n",pattern);
473:   MPI_Comm_rank(comm,&rank);
474:   MPI_Comm_size(comm,&size);

476:   if (pattern == PETSCSF_PATTERN_ALLTOALL) {
477:     type = PETSCSFALLTOALL;
478:     PetscLayoutCreate(comm,&sf->map);
479:     PetscLayoutSetLocalSize(sf->map,size);
480:     PetscLayoutSetSize(sf->map,((PetscInt)size)*size);
481:     PetscLayoutSetUp(sf->map);
482:   } else {
483:     PetscLayoutGetLocalSize(map,&n);
484:     PetscLayoutGetSize(map,&N);
485:     res[0] = n;
486:     res[1] = -n;
487:     /* Check if n are same over all ranks so that we can optimize it */
488:     MPIU_Allreduce(MPI_IN_PLACE,res,2,MPIU_INT,MPI_MAX,comm);
489:     if (res[0] == -res[1]) { /* same n */
490:       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHER  : PETSCSFGATHER;
491:     } else {
492:       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHERV : PETSCSFGATHERV;
493:     }
494:     PetscLayoutReference(map,&sf->map);
495:   }
496:   PetscSFSetType(sf,type);

498:   sf->pattern = pattern;
499:   sf->mine    = NULL; /* Contiguous */

501:   /* Set nleaves, nroots here in case user calls PetscSFGetGraph, which is legal to call even before PetscSFSetUp is called.
502:      Also set other easy stuff.
503:    */
504:   if (pattern == PETSCSF_PATTERN_ALLGATHER) {
505:     sf->nleaves      = N;
506:     sf->nroots       = n;
507:     sf->nranks       = size;
508:     sf->minleaf      = 0;
509:     sf->maxleaf      = N - 1;
510:   } else if (pattern == PETSCSF_PATTERN_GATHER) {
511:     sf->nleaves      = rank ? 0 : N;
512:     sf->nroots       = n;
513:     sf->nranks       = rank ? 0 : size;
514:     sf->minleaf      = 0;
515:     sf->maxleaf      = rank ? -1 : N - 1;
516:   } else if (pattern == PETSCSF_PATTERN_ALLTOALL) {
517:     sf->nleaves      = size;
518:     sf->nroots       = size;
519:     sf->nranks       = size;
520:     sf->minleaf      = 0;
521:     sf->maxleaf      = size - 1;
522:   }
523:   sf->ndranks  = 0; /* We do not need to separate out distinguished ranks for patterned graphs to improve communication performance */
524:   sf->graphset = PETSC_TRUE;
525:   return(0);
526: }

528: /*@
529:    PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map

531:    Collective

533:    Input Arguments:

535: .  sf - star forest to invert

537:    Output Arguments:
538: .  isf - inverse of sf
539:    Level: advanced

541:    Notes:
542:    All roots must have degree 1.

544:    The local space may be a permutation, but cannot be sparse.

546: .seealso: PetscSFSetGraph()
547: @*/
548: PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
549: {
551:   PetscMPIInt    rank;
552:   PetscInt       i,nroots,nleaves,maxlocal,count,*newilocal;
553:   const PetscInt *ilocal;
554:   PetscSFNode    *roots,*leaves;

558:   PetscSFCheckGraphSet(sf,1);

561:   PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);
562:   maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */

564:   MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);
565:   PetscMalloc2(nroots,&roots,maxlocal,&leaves);
566:   for (i=0; i<maxlocal; i++) {
567:     leaves[i].rank  = rank;
568:     leaves[i].index = i;
569:   }
570:   for (i=0; i <nroots; i++) {
571:     roots[i].rank  = -1;
572:     roots[i].index = -1;
573:   }
574:   PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);
575:   PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);

577:   /* Check whether our leaves are sparse */
578:   for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
579:   if (count == nroots) newilocal = NULL;
580:   else {                        /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
581:     PetscMalloc1(count,&newilocal);
582:     for (i=0,count=0; i<nroots; i++) {
583:       if (roots[i].rank >= 0) {
584:         newilocal[count]   = i;
585:         roots[count].rank  = roots[i].rank;
586:         roots[count].index = roots[i].index;
587:         count++;
588:       }
589:     }
590:   }

592:   PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);
593:   PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);
594:   PetscFree2(roots,leaves);
595:   return(0);
596: }

598: /*@
599:    PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph

601:    Collective

603:    Input Arguments:
604: +  sf - communication object to duplicate
605: -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)

607:    Output Arguments:
608: .  newsf - new communication object

610:    Level: beginner

612: .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph()
613: @*/
614: PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
615: {
616:   PetscSFType    type;

623:   PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);
624:   PetscSFGetType(sf,&type);
625:   if (type) {PetscSFSetType(*newsf,type);}
626:   if (opt == PETSCSF_DUPLICATE_GRAPH) {
627:     PetscSFCheckGraphSet(sf,1);
628:     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
629:       PetscInt          nroots,nleaves;
630:       const PetscInt    *ilocal;
631:       const PetscSFNode *iremote;
632:       PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);
633:       PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);
634:     } else {
635:       PetscSFSetGraphWithPattern(*newsf,sf->map,sf->pattern);
636:     }
637:   }
638:   if (sf->ops->Duplicate) {(*sf->ops->Duplicate)(sf,opt,*newsf);}
639:   return(0);
640: }

642: /*@C
643:    PetscSFGetGraph - Get the graph specifying a parallel star forest

645:    Not Collective

647:    Input Arguments:
648: .  sf - star forest

650:    Output Arguments:
651: +  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
652: .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
653: .  ilocal - locations of leaves in leafdata buffers
654: -  iremote - remote locations of root vertices for each leaf on the current process

656:    Notes:
657:    We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set yet

659:    When called from Fortran, the returned iremote array is a copy and must deallocated after use. Consequently, if you
660:    want to update the graph, you must call PetscSFSetGraph after modifying the iremote array.

662:    Level: intermediate

664: .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph()
665: @*/
666: PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
667: {

672:   if (sf->ops->GetGraph) {
673:     (sf->ops->GetGraph)(sf,nroots,nleaves,ilocal,iremote);
674:   } else {
675:     if (nroots) *nroots = sf->nroots;
676:     if (nleaves) *nleaves = sf->nleaves;
677:     if (ilocal) *ilocal = sf->mine;
678:     if (iremote) *iremote = sf->remote;
679:   }
680:   return(0);
681: }

683: /*@
684:    PetscSFGetLeafRange - Get the active leaf ranges

686:    Not Collective

688:    Input Arguments:
689: .  sf - star forest

691:    Output Arguments:
692: +  minleaf - minimum active leaf on this process. Return 0 if there are no leaves.
693: -  maxleaf - maximum active leaf on this process. Return -1 if there are no leaves.

695:    Level: developer

697: .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
698: @*/
699: PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
700: {

704:   PetscSFCheckGraphSet(sf,1);
705:   if (minleaf) *minleaf = sf->minleaf;
706:   if (maxleaf) *maxleaf = sf->maxleaf;
707:   return(0);
708: }

710: /*@C
711:    PetscSFViewFromOptions - View from Options

713:    Collective on PetscSF

715:    Input Parameters:
716: +  A - the star forest
717: .  obj - Optional object
718: -  name - command line option

720:    Level: intermediate
721: .seealso:  PetscSF, PetscSFView, PetscObjectViewFromOptions(), PetscSFCreate()
722: @*/
723: PetscErrorCode  PetscSFViewFromOptions(PetscSF A,PetscObject obj,const char name[])
724: {

729:   PetscObjectViewFromOptions((PetscObject)A,obj,name);
730:   return(0);
731: }

733: /*@C
734:    PetscSFView - view a star forest

736:    Collective

738:    Input Arguments:
739: +  sf - star forest
740: -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD

742:    Level: beginner

744: .seealso: PetscSFCreate(), PetscSFSetGraph()
745: @*/
746: PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
747: {
748:   PetscErrorCode    ierr;
749:   PetscBool         iascii;
750:   PetscViewerFormat format;

754:   if (!viewer) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);}
757:   if (sf->graphset) {PetscSFSetUp(sf);}
758:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
759:   if (iascii) {
760:     PetscMPIInt rank;
761:     PetscInt    ii,i,j;

763:     PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);
764:     PetscViewerASCIIPushTab(viewer);
765:     if (sf->ops->View) {(*sf->ops->View)(sf,viewer);}
766:     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
767:       if (!sf->graphset) {
768:         PetscViewerASCIIPrintf(viewer,"PetscSFSetGraph() has not been called yet\n");
769:         PetscViewerASCIIPopTab(viewer);
770:         return(0);
771:       }
772:       MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);
773:       PetscViewerASCIIPushSynchronized(viewer);
774:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);
775:       for (i=0; i<sf->nleaves; i++) {
776:         PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- (%D,%D)\n",rank,sf->mine ? sf->mine[i] : i,sf->remote[i].rank,sf->remote[i].index);
777:       }
778:       PetscViewerFlush(viewer);
779:       PetscViewerGetFormat(viewer,&format);
780:       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
781:         PetscMPIInt *tmpranks,*perm;
782:         PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm);
783:         PetscArraycpy(tmpranks,sf->ranks,sf->nranks);
784:         for (i=0; i<sf->nranks; i++) perm[i] = i;
785:         PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm);
786:         PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);
787:         for (ii=0; ii<sf->nranks; ii++) {
788:           i = perm[ii];
789:           PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);
790:           for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
791:             PetscViewerASCIISynchronizedPrintf(viewer,"[%d]    %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);
792:           }
793:         }
794:         PetscFree2(tmpranks,perm);
795:       }
796:       PetscViewerFlush(viewer);
797:       PetscViewerASCIIPopSynchronized(viewer);
798:     }
799:     PetscViewerASCIIPopTab(viewer);
800:   }
801:   return(0);
802: }

804: /*@C
805:    PetscSFGetRootRanks - Get root ranks and number of vertices referenced by leaves on this process

807:    Not Collective

809:    Input Arguments:
810: .  sf - star forest

812:    Output Arguments:
813: +  nranks - number of ranks referenced by local part
814: .  ranks - array of ranks
815: .  roffset - offset in rmine/rremote for each rank (length nranks+1)
816: .  rmine - concatenated array holding local indices referencing each remote rank
817: -  rremote - concatenated array holding remote indices referenced for each remote rank

819:    Level: developer

821: .seealso: PetscSFGetLeafRanks()
822: @*/
823: PetscErrorCode PetscSFGetRootRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
824: {

829:   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
830:   if (sf->ops->GetRootRanks) {
831:     (sf->ops->GetRootRanks)(sf,nranks,ranks,roffset,rmine,rremote);
832:   } else {
833:     /* The generic implementation */
834:     if (nranks)  *nranks  = sf->nranks;
835:     if (ranks)   *ranks   = sf->ranks;
836:     if (roffset) *roffset = sf->roffset;
837:     if (rmine)   *rmine   = sf->rmine;
838:     if (rremote) *rremote = sf->rremote;
839:   }
840:   return(0);
841: }

843: /*@C
844:    PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process

846:    Not Collective

848:    Input Arguments:
849: .  sf - star forest

851:    Output Arguments:
852: +  niranks - number of leaf ranks referencing roots on this process
853: .  iranks - array of ranks
854: .  ioffset - offset in irootloc for each rank (length niranks+1)
855: -  irootloc - concatenated array holding local indices of roots referenced by each leaf rank

857:    Level: developer

859: .seealso: PetscSFGetRootRanks()
860: @*/
861: PetscErrorCode PetscSFGetLeafRanks(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
862: {

867:   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
868:   if (sf->ops->GetLeafRanks) {
869:     (sf->ops->GetLeafRanks)(sf,niranks,iranks,ioffset,irootloc);
870:   } else {
871:     PetscSFType type;
872:     PetscSFGetType(sf,&type);
873:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type);
874:   }
875:   return(0);
876: }

878: static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) {
879:   PetscInt i;
880:   for (i=0; i<n; i++) {
881:     if (needle == list[i]) return PETSC_TRUE;
882:   }
883:   return PETSC_FALSE;
884: }

886: /*@C
887:    PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations.

889:    Collective

891:    Input Arguments:
892: +  sf - PetscSF to set up; PetscSFSetGraph() must have been called
893: -  dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange)

895:    Level: developer

897: .seealso: PetscSFGetRootRanks()
898: @*/
899: PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup)
900: {
901:   PetscErrorCode     ierr;
902:   PetscTable         table;
903:   PetscTablePosition pos;
904:   PetscMPIInt        size,groupsize,*groupranks;
905:   PetscInt           *rcount,*ranks;
906:   PetscInt           i, irank = -1,orank = -1;

910:   PetscSFCheckGraphSet(sf,1);
911:   MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);
912:   PetscTableCreate(10,size,&table);
913:   for (i=0; i<sf->nleaves; i++) {
914:     /* Log 1-based rank */
915:     PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);
916:   }
917:   PetscTableGetCount(table,&sf->nranks);
918:   PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);
919:   PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);
920:   PetscTableGetHeadPosition(table,&pos);
921:   for (i=0; i<sf->nranks; i++) {
922:     PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);
923:     ranks[i]--;             /* Convert back to 0-based */
924:   }
925:   PetscTableDestroy(&table);

927:   /* We expect that dgroup is reliably "small" while nranks could be large */
928:   {
929:     MPI_Group group = MPI_GROUP_NULL;
930:     PetscMPIInt *dgroupranks;
931:     MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);
932:     MPI_Group_size(dgroup,&groupsize);
933:     PetscMalloc1(groupsize,&dgroupranks);
934:     PetscMalloc1(groupsize,&groupranks);
935:     for (i=0; i<groupsize; i++) dgroupranks[i] = i;
936:     if (groupsize) {MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);}
937:     MPI_Group_free(&group);
938:     PetscFree(dgroupranks);
939:   }

941:   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
942:   for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i; ) {
943:     for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */
944:       if (InList(ranks[i],groupsize,groupranks)) break;
945:     }
946:     for ( ; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
947:       if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break;
948:     }
949:     if (sf->ndranks < i) {                         /* Swap ranks[sf->ndranks] with ranks[i] */
950:       PetscInt    tmprank,tmpcount;

952:       tmprank             = ranks[i];
953:       tmpcount            = rcount[i];
954:       ranks[i]            = ranks[sf->ndranks];
955:       rcount[i]           = rcount[sf->ndranks];
956:       ranks[sf->ndranks]  = tmprank;
957:       rcount[sf->ndranks] = tmpcount;
958:       sf->ndranks++;
959:     }
960:   }
961:   PetscFree(groupranks);
962:   PetscSortIntWithArray(sf->ndranks,ranks,rcount);
963:   PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);
964:   sf->roffset[0] = 0;
965:   for (i=0; i<sf->nranks; i++) {
966:     PetscMPIIntCast(ranks[i],sf->ranks+i);
967:     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
968:     rcount[i]        = 0;
969:   }
970:   for (i=0, irank = -1, orank = -1; i<sf->nleaves; i++) {
971:     /* short circuit */
972:     if (orank != sf->remote[i].rank) {
973:       /* Search for index of iremote[i].rank in sf->ranks */
974:       PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);
975:       if (irank < 0) {
976:         PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);
977:         if (irank >= 0) irank += sf->ndranks;
978:       }
979:       orank = sf->remote[i].rank;
980:     }
981:     if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank);
982:     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
983:     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
984:     rcount[irank]++;
985:   }
986:   PetscFree2(rcount,ranks);
987:   return(0);
988: }

990: /*@C
991:    PetscSFGetGroups - gets incoming and outgoing process groups

993:    Collective

995:    Input Argument:
996: .  sf - star forest

998:    Output Arguments:
999: +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
1000: -  outgoing - group of destination processes for outgoing edges (roots that I reference)

1002:    Level: developer

1004: .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
1005: @*/
1006: PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
1007: {
1009:   MPI_Group      group = MPI_GROUP_NULL;

1012:   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining groups");
1013:   if (sf->ingroup == MPI_GROUP_NULL) {
1014:     PetscInt       i;
1015:     const PetscInt *indegree;
1016:     PetscMPIInt    rank,*outranks,*inranks;
1017:     PetscSFNode    *remote;
1018:     PetscSF        bgcount;

1020:     /* Compute the number of incoming ranks */
1021:     PetscMalloc1(sf->nranks,&remote);
1022:     for (i=0; i<sf->nranks; i++) {
1023:       remote[i].rank  = sf->ranks[i];
1024:       remote[i].index = 0;
1025:     }
1026:     PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);
1027:     PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);
1028:     PetscSFComputeDegreeBegin(bgcount,&indegree);
1029:     PetscSFComputeDegreeEnd(bgcount,&indegree);

1031:     /* Enumerate the incoming ranks */
1032:     PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);
1033:     MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);
1034:     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
1035:     PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);
1036:     PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);
1037:     MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);
1038:     MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);
1039:     MPI_Group_free(&group);
1040:     PetscFree2(inranks,outranks);
1041:     PetscSFDestroy(&bgcount);
1042:   }
1043:   *incoming = sf->ingroup;

1045:   if (sf->outgroup == MPI_GROUP_NULL) {
1046:     MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);
1047:     MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);
1048:     MPI_Group_free(&group);
1049:   }
1050:   *outgoing = sf->outgroup;
1051:   return(0);
1052: }

1054: /*@
1055:    PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters

1057:    Collective

1059:    Input Argument:
1060: .  sf - star forest that may contain roots with 0 or with more than 1 vertex

1062:    Output Arguments:
1063: .  multi - star forest with split roots, such that each root has degree exactly 1

1065:    Level: developer

1067:    Notes:

1069:    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
1070:    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
1071:    edge, it is a candidate for future optimization that might involve its removal.

1073: .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin(), PetscSFComputeMultiRootOriginalNumbering()
1074: @*/
1075: PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
1076: {

1082:   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
1083:     PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);
1084:     *multi = sf->multi;
1085:     return(0);
1086:   }
1087:   if (!sf->multi) {
1088:     const PetscInt *indegree;
1089:     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
1090:     PetscSFNode    *remote;
1091:     maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
1092:     PetscSFComputeDegreeBegin(sf,&indegree);
1093:     PetscSFComputeDegreeEnd(sf,&indegree);
1094:     PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);
1095:     inoffset[0] = 0;
1096:     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
1097:     for (i=0; i<maxlocal; i++) outones[i] = 1;
1098:     PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);
1099:     PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);
1100:     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
1101: #if 0
1102: #if defined(PETSC_USE_DEBUG)                                 /* Check that the expected number of increments occurred */
1103:     for (i=0; i<sf->nroots; i++) {
1104:       if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
1105:     }
1106: #endif
1107: #endif
1108:     PetscMalloc1(sf->nleaves,&remote);
1109:     for (i=0; i<sf->nleaves; i++) {
1110:       remote[i].rank  = sf->remote[i].rank;
1111:       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
1112:     }
1113:     PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);
1114:     PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);
1115:     if (sf->rankorder) {        /* Sort the ranks */
1116:       PetscMPIInt rank;
1117:       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
1118:       PetscSFNode *newremote;
1119:       MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);
1120:       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
1121:       PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);
1122:       for (i=0; i<maxlocal; i++) outranks[i] = rank;
1123:       PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);
1124:       PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);
1125:       /* Sort the incoming ranks at each vertex, build the inverse map */
1126:       for (i=0; i<sf->nroots; i++) {
1127:         PetscInt j;
1128:         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
1129:         PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);
1130:         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
1131:       }
1132:       PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);
1133:       PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);
1134:       PetscMalloc1(sf->nleaves,&newremote);
1135:       for (i=0; i<sf->nleaves; i++) {
1136:         newremote[i].rank  = sf->remote[i].rank;
1137:         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
1138:       }
1139:       PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);
1140:       PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);
1141:     }
1142:     PetscFree3(inoffset,outones,outoffset);
1143:   }
1144:   *multi = sf->multi;
1145:   return(0);
1146: }

1148: /*@C
1149:    PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices

1151:    Collective

1153:    Input Arguments:
1154: +  sf - original star forest
1155: .  nselected  - number of selected roots on this process
1156: -  selected   - indices of the selected roots on this process

1158:    Output Arguments:
1159: .  newsf - new star forest

1161:    Level: advanced

1163:    Note:
1164:    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
1165:    be done by calling PetscSFGetGraph().

1167: .seealso: PetscSFSetGraph(), PetscSFGetGraph()
1168: @*/
1169: PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
1170: {
1171:   PetscInt          i,n,*roots,*rootdata,*leafdata,nroots,nleaves,connected_leaves,*new_ilocal;
1172:   const PetscSFNode *iremote;
1173:   PetscSFNode       *new_iremote;
1174:   PetscSF           tmpsf;
1175:   MPI_Comm          comm;
1176:   PetscErrorCode    ierr;

1180:   PetscSFCheckGraphSet(sf,1);

1184:   PetscSFSetUp(sf);

1186:   /* Uniq selected[] and put results in roots[] */
1187:   PetscObjectGetComm((PetscObject)sf,&comm);
1188:   PetscMalloc1(nselected,&roots);
1189:   PetscArraycpy(roots,selected,nselected);
1190:   PetscSortedRemoveDupsInt(&nselected,roots);
1191:   if (nselected && (roots[0] < 0 || roots[nselected-1] >= sf->nroots)) SETERRQ3(comm,PETSC_ERR_ARG_OUTOFRANGE,"Min/Max root indices %D/%D are not in [0,%D)",roots[0],roots[nselected-1],sf->nroots);

1193:   if (sf->ops->CreateEmbeddedSF) {
1194:     (*sf->ops->CreateEmbeddedSF)(sf,nselected,roots,newsf);
1195:   } else {
1196:     /* A generic version of creating embedded sf. Note that we called PetscSFSetGraph() twice, which is certainly expensive */
1197:     /* Find out which leaves (not leaf data items) are still connected to roots in the embedded sf */
1198:     PetscSFGetGraph(sf,&nroots,&nleaves,NULL,&iremote);
1199:     PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&tmpsf);
1200:     PetscSFSetGraph(tmpsf,nroots,nleaves,NULL/*contiguous*/,PETSC_USE_POINTER,iremote,PETSC_USE_POINTER);
1201:     PetscCalloc2(nroots,&rootdata,nleaves,&leafdata);
1202:     for (i=0; i<nselected; ++i) rootdata[roots[i]] = 1;
1203:     PetscSFBcastBegin(tmpsf,MPIU_INT,rootdata,leafdata);
1204:     PetscSFBcastEnd(tmpsf,MPIU_INT,rootdata,leafdata);
1205:     PetscSFDestroy(&tmpsf);

1207:     /* Build newsf with leaves that are still connected */
1208:     connected_leaves = 0;
1209:     for (i=0; i<nleaves; ++i) connected_leaves += leafdata[i];
1210:     PetscMalloc1(connected_leaves,&new_ilocal);
1211:     PetscMalloc1(connected_leaves,&new_iremote);
1212:     for (i=0, n=0; i<nleaves; ++i) {
1213:       if (leafdata[i]) {
1214:         new_ilocal[n]        = sf->mine ? sf->mine[i] : i;
1215:         new_iremote[n].rank  = sf->remote[i].rank;
1216:         new_iremote[n].index = sf->remote[i].index;
1217:         ++n;
1218:       }
1219:     }

1221:     if (n != connected_leaves) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"There is a size mismatch in the SF embedding, %d != %d",n,connected_leaves);
1222:     PetscSFCreate(comm,newsf);
1223:     PetscSFSetFromOptions(*newsf);
1224:     PetscSFSetGraph(*newsf,nroots,connected_leaves,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);
1225:     PetscFree2(rootdata,leafdata);
1226:   }
1227:   PetscFree(roots);
1228:   return(0);
1229: }

1231: /*@C
1232:   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices

1234:   Collective

1236:   Input Arguments:
1237: + sf - original star forest
1238: . nselected  - number of selected leaves on this process
1239: - selected   - indices of the selected leaves on this process

1241:   Output Arguments:
1242: .  newsf - new star forest

1244:   Level: advanced

1246: .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph()
1247: @*/
1248: PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
1249: {
1250:   const PetscSFNode *iremote;
1251:   PetscSFNode       *new_iremote;
1252:   const PetscInt    *ilocal;
1253:   PetscInt          i,nroots,*leaves,*new_ilocal;
1254:   MPI_Comm          comm;
1255:   PetscErrorCode    ierr;

1259:   PetscSFCheckGraphSet(sf,1);

1263:   /* Uniq selected[] and put results in leaves[] */
1264:   PetscObjectGetComm((PetscObject)sf,&comm);
1265:   PetscMalloc1(nselected,&leaves);
1266:   PetscArraycpy(leaves,selected,nselected);
1267:   PetscSortedRemoveDupsInt(&nselected,leaves);
1268:   if (nselected && (leaves[0] < 0 || leaves[nselected-1] >= sf->nleaves)) SETERRQ3(comm,PETSC_ERR_ARG_OUTOFRANGE,"Min/Max leaf indices %D/%D are not in [0,%D)",leaves[0],leaves[nselected-1],sf->nleaves);

1270:   /* Optimize the routine only when sf is setup and hence we can reuse sf's communication pattern */
1271:   if (sf->setupcalled && sf->ops->CreateEmbeddedLeafSF) {
1272:     (*sf->ops->CreateEmbeddedLeafSF)(sf,nselected,leaves,newsf);
1273:   } else {
1274:     PetscSFGetGraph(sf,&nroots,NULL,&ilocal,&iremote);
1275:     PetscMalloc1(nselected,&new_ilocal);
1276:     PetscMalloc1(nselected,&new_iremote);
1277:     for (i=0; i<nselected; ++i) {
1278:       const PetscInt l     = leaves[i];
1279:       new_ilocal[i]        = ilocal ? ilocal[l] : l;
1280:       new_iremote[i].rank  = iremote[l].rank;
1281:       new_iremote[i].index = iremote[l].index;
1282:     }
1283:     PetscSFCreate(comm,newsf);
1284:     PetscSFSetFromOptions(*newsf);
1285:     PetscSFSetGraph(*newsf,nroots,nselected,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);
1286:   }
1287:   PetscFree(leaves);
1288:   return(0);
1289: }

1291: /*@C
1292:    PetscSFBcastAndOpBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastAndOpEnd()

1294:    Collective on PetscSF

1296:    Input Arguments:
1297: +  sf - star forest on which to communicate
1298: .  unit - data type associated with each node
1299: .  rootdata - buffer to broadcast
1300: -  op - operation to use for reduction

1302:    Output Arguments:
1303: .  leafdata - buffer to be reduced with values from each leaf's respective root

1305:    Level: intermediate

1307: .seealso: PetscSFBcastAndOpEnd(), PetscSFBcastBegin(), PetscSFBcastEnd()
1308: @*/
1309: PetscErrorCode PetscSFBcastAndOpBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
1310: {
1312:   PetscMemType   rootmtype,leafmtype;

1316:   PetscSFSetUp(sf);
1317:   PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);
1318:   PetscGetMemType(rootdata,&rootmtype);
1319:   PetscGetMemType(leafdata,&leafmtype);
1320: #if defined(PETSC_HAVE_CUDA)
1321:   /*  Shall we assume rootdata, leafdata are ready to use (instead of being computed by some asynchronous kernels)?
1322:     To be similar to MPI, I'd like to have this assumption, since MPI does not have a concept of stream.
1323:     But currently this assumption is not enforecd in Petsc. To be safe, I do synchronization here. Otherwise, if
1324:     we do not sync now and call the Pack kernel directly on the default NULL stream (assume petsc objects are also
1325:     computed on it), we have to sync the NULL stream before calling MPI routines. So, it looks a cudaDeviceSynchronize
1326:     is inevitable. We do it now and put pack/unpack kernels to non-NULL streams.
1327:    */
1328:   if (rootmtype == PETSC_MEMTYPE_DEVICE || leafmtype == PETSC_MEMTYPE_DEVICE) {cudaError_t err = cudaDeviceSynchronize();CHKERRCUDA(err);}
1329: #endif
1330:   (*sf->ops->BcastAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op);
1331:   PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);
1332:   return(0);
1333: }

1335: /*@C
1336:    PetscSFBcastAndOpEnd - end a broadcast & reduce operation started with PetscSFBcastAndOpBegin()

1338:    Collective

1340:    Input Arguments:
1341: +  sf - star forest
1342: .  unit - data type
1343: .  rootdata - buffer to broadcast
1344: -  op - operation to use for reduction

1346:    Output Arguments:
1347: .  leafdata - buffer to be reduced with values from each leaf's respective root

1349:    Level: intermediate

1351: .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
1352: @*/
1353: PetscErrorCode PetscSFBcastAndOpEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
1354: {
1356:   PetscMemType   rootmtype,leafmtype;

1360:   PetscSFSetUp(sf);
1361:   PetscLogEventBegin(PETSCSF_BcastAndOpEnd,sf,0,0,0);
1362:   PetscGetMemType(rootdata,&rootmtype);
1363:   PetscGetMemType(leafdata,&leafmtype);
1364:   (*sf->ops->BcastAndOpEnd)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op);
1365:   PetscLogEventEnd(PETSCSF_BcastAndOpEnd,sf,0,0,0);
1366:   return(0);
1367: }

1369: /*@C
1370:    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()

1372:    Collective

1374:    Input Arguments:
1375: +  sf - star forest
1376: .  unit - data type
1377: .  leafdata - values to reduce
1378: -  op - reduction operation

1380:    Output Arguments:
1381: .  rootdata - result of reduction of values from all leaves of each root

1383:    Level: intermediate

1385: .seealso: PetscSFBcastBegin()
1386: @*/
1387: PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
1388: {
1390:   PetscMemType   rootmtype,leafmtype;

1394:   PetscSFSetUp(sf);
1395:   PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);
1396:   PetscGetMemType(rootdata,&rootmtype);
1397:   PetscGetMemType(leafdata,&leafmtype);
1398: #if defined(PETSC_HAVE_CUDA)
1399:   if (rootmtype == PETSC_MEMTYPE_DEVICE || leafmtype == PETSC_MEMTYPE_DEVICE) {cudaError_t err = cudaDeviceSynchronize();CHKERRCUDA(err);}
1400: #endif
1401:   (sf->ops->ReduceBegin)(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op);
1402:   PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);
1403:   return(0);
1404: }

1406: /*@C
1407:    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()

1409:    Collective

1411:    Input Arguments:
1412: +  sf - star forest
1413: .  unit - data type
1414: .  leafdata - values to reduce
1415: -  op - reduction operation

1417:    Output Arguments:
1418: .  rootdata - result of reduction of values from all leaves of each root

1420:    Level: intermediate

1422: .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
1423: @*/
1424: PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
1425: {
1427:   PetscMemType   rootmtype,leafmtype;

1431:   PetscSFSetUp(sf);
1432:   PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);
1433:   PetscGetMemType(rootdata,&rootmtype);
1434:   PetscGetMemType(leafdata,&leafmtype);
1435:   (*sf->ops->ReduceEnd)(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op);
1436:   PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);
1437:   return(0);
1438: }

1440: /*@C
1441:    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()

1443:    Collective

1445:    Input Arguments:
1446: +  sf - star forest
1447: .  unit - data type
1448: .  leafdata - leaf values to use in reduction
1449: -  op - operation to use for reduction

1451:    Output Arguments:
1452: +  rootdata - root values to be updated, input state is seen by first process to perform an update
1453: -  leafupdate - state at each leaf's respective root immediately prior to my atomic update

1455:    Level: advanced

1457:    Note:
1458:    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
1459:    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
1460:    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
1461:    integers.

1463: .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
1464: @*/
1465: PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1466: {
1468:   PetscMemType   rootmtype,leafmtype,leafupdatemtype;

1472:   PetscSFSetUp(sf);
1473:   PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);
1474:   PetscGetMemType(rootdata,&rootmtype);
1475:   PetscGetMemType(leafdata,&leafmtype);
1476:   PetscGetMemType(leafupdate,&leafupdatemtype);
1477:   if (leafmtype != leafupdatemtype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for leafdata and leafupdate in different memory types");
1478: #if defined(PETSC_HAVE_CUDA)
1479:   if (rootmtype == PETSC_MEMTYPE_DEVICE || leafmtype == PETSC_MEMTYPE_DEVICE) {cudaError_t err = cudaDeviceSynchronize();CHKERRCUDA(err);}
1480: #endif
1481:   (*sf->ops->FetchAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,leafupdate,op);
1482:   PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);
1483:   return(0);
1484: }

1486: /*@C
1487:    PetscSFFetchAndOpEnd - end operation started in matching call to PetscSFFetchAndOpBegin() to fetch values from roots and update atomically by applying operation using my leaf value

1489:    Collective

1491:    Input Arguments:
1492: +  sf - star forest
1493: .  unit - data type
1494: .  leafdata - leaf values to use in reduction
1495: -  op - operation to use for reduction

1497:    Output Arguments:
1498: +  rootdata - root values to be updated, input state is seen by first process to perform an update
1499: -  leafupdate - state at each leaf's respective root immediately prior to my atomic update

1501:    Level: advanced

1503: .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
1504: @*/
1505: PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1506: {
1508:   PetscMemType   rootmtype,leafmtype,leafupdatemtype;

1512:   PetscSFSetUp(sf);
1513:   PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);
1514:   PetscGetMemType(rootdata,&rootmtype);
1515:   PetscGetMemType(leafdata,&leafmtype);
1516:   PetscGetMemType(leafupdate,&leafupdatemtype);
1517:   if (leafmtype != leafupdatemtype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for leafdata and leafupdate in different memory types");
1518:   (*sf->ops->FetchAndOpEnd)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,leafupdate,op);
1519:   PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);
1520:   return(0);
1521: }

1523: /*@C
1524:    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()

1526:    Collective

1528:    Input Arguments:
1529: .  sf - star forest

1531:    Output Arguments:
1532: .  degree - degree of each root vertex

1534:    Level: advanced

1536:    Notes:
1537:    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.

1539: .seealso: PetscSFGatherBegin()
1540: @*/
1541: PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
1542: {

1547:   PetscSFCheckGraphSet(sf,1);
1549:   if (!sf->degreeknown) {
1550:     PetscInt i, nroots = sf->nroots, maxlocal = sf->maxleaf+1;  /* TODO: We should use PetscSFGetLeafRange() */
1551:     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
1552:     PetscMalloc1(nroots,&sf->degree);
1553:     PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
1554:     for (i=0; i<nroots; i++) sf->degree[i] = 0;
1555:     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
1556:     PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);
1557:   }
1558:   *degree = NULL;
1559:   return(0);
1560: }

1562: /*@C
1563:    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()

1565:    Collective

1567:    Input Arguments:
1568: .  sf - star forest

1570:    Output Arguments:
1571: .  degree - degree of each root vertex

1573:    Level: developer

1575:    Notes:
1576:    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.

1578: .seealso:
1579: @*/
1580: PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
1581: {

1586:   PetscSFCheckGraphSet(sf,1);
1588:   if (!sf->degreeknown) {
1589:     if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
1590:     PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);
1591:     PetscFree(sf->degreetmp);
1592:     sf->degreeknown = PETSC_TRUE;
1593:   }
1594:   *degree = sf->degree;
1595:   return(0);
1596: }


1599: /*@C
1600:    PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()).
1601:    Each multi-root is assigned index of the corresponding original root.

1603:    Collective

1605:    Input Arguments:
1606: +  sf - star forest
1607: -  degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd()

1609:    Output Arguments:
1610: +  nMultiRoots - (optional) number of multi-roots (roots of multi-SF)
1611: -  multiRootsOrigNumbering - original indices of multi-roots; length of this array is nMultiRoots

1613:    Level: developer

1615:    Notes:
1616:    The returned array multiRootsOrigNumbering is newly allocated and should be destroyed with PetscFree() when no longer needed.

1618: .seealso: PetscSFComputeDegreeBegin(), PetscSFComputeDegreeEnd(), PetscSFGetMultiSF()
1619: @*/
1620: PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[])
1621: {
1622:   PetscSF             msf;
1623:   PetscInt            i, j, k, nroots, nmroots;
1624:   PetscErrorCode      ierr;

1628:   PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);
1632:   PetscSFGetMultiSF(sf,&msf);
1633:   PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL);
1634:   PetscMalloc1(nmroots, multiRootsOrigNumbering);
1635:   for (i=0,j=0,k=0; i<nroots; i++) {
1636:     if (!degree[i]) continue;
1637:     for (j=0; j<degree[i]; j++,k++) {
1638:       (*multiRootsOrigNumbering)[k] = i;
1639:     }
1640:   }
1641: #if defined(PETSC_USE_DEBUG)
1642:   if (PetscUnlikely(k != nmroots)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail");
1643: #endif
1644:   if (nMultiRoots) *nMultiRoots = nmroots;
1645:   return(0);
1646: }

1648: /*@C
1649:    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()

1651:    Collective

1653:    Input Arguments:
1654: +  sf - star forest
1655: .  unit - data type
1656: -  leafdata - leaf data to gather to roots

1658:    Output Argument:
1659: .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree

1661:    Level: intermediate

1663: .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
1664: @*/
1665: PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
1666: {
1668:   PetscSF        multi;

1672:   PetscSFSetUp(sf);
1673:   PetscSFGetMultiSF(sf,&multi);
1674:   PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);
1675:   return(0);
1676: }

1678: /*@C
1679:    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()

1681:    Collective

1683:    Input Arguments:
1684: +  sf - star forest
1685: .  unit - data type
1686: -  leafdata - leaf data to gather to roots

1688:    Output Argument:
1689: .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree

1691:    Level: intermediate

1693: .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
1694: @*/
1695: PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
1696: {
1698:   PetscSF        multi;

1702:   PetscSFSetUp(sf);
1703:   PetscSFGetMultiSF(sf,&multi);
1704:   PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);
1705:   return(0);
1706: }

1708: /*@C
1709:    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()

1711:    Collective

1713:    Input Arguments:
1714: +  sf - star forest
1715: .  unit - data type
1716: -  multirootdata - root buffer to send to each leaf, one unit of data per leaf

1718:    Output Argument:
1719: .  leafdata - leaf data to be update with personal data from each respective root

1721:    Level: intermediate

1723: .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
1724: @*/
1725: PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
1726: {
1728:   PetscSF        multi;

1732:   PetscSFSetUp(sf);
1733:   PetscSFGetMultiSF(sf,&multi);
1734:   PetscSFBcastBegin(multi,unit,multirootdata,leafdata);
1735:   return(0);
1736: }

1738: /*@C
1739:    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()

1741:    Collective

1743:    Input Arguments:
1744: +  sf - star forest
1745: .  unit - data type
1746: -  multirootdata - root buffer to send to each leaf, one unit of data per leaf

1748:    Output Argument:
1749: .  leafdata - leaf data to be update with personal data from each respective root

1751:    Level: intermediate

1753: .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
1754: @*/
1755: PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
1756: {
1758:   PetscSF        multi;

1762:   PetscSFSetUp(sf);
1763:   PetscSFGetMultiSF(sf,&multi);
1764:   PetscSFBcastEnd(multi,unit,multirootdata,leafdata);
1765:   return(0);
1766: }

1768: static PetscErrorCode PetscSFCheckLeavesUnique_Private(PetscSF sf)
1769: {
1770: #if defined(PETSC_USE_DEBUG)
1771:   PetscInt        i, n, nleaves;
1772:   const PetscInt *ilocal = NULL;
1773:   PetscHSetI      seen;
1774:   PetscErrorCode  ierr;

1777:   PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,NULL);
1778:   PetscHSetICreate(&seen);
1779:   for (i = 0; i < nleaves; i++) {
1780:     const PetscInt leaf = ilocal ? ilocal[i] : i;
1781:     PetscHSetIAdd(seen,leaf);
1782:   }
1783:   PetscHSetIGetSize(seen,&n);
1784:   if (n != nleaves) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided leaves have repeated values: all leaves must be unique");
1785:   PetscHSetIDestroy(&seen);
1786:   return(0);
1787: #else
1789:   return(0);
1790: #endif
1791: }

1793: /*@
1794:   PetscSFCompose - Compose a new PetscSF by putting the second SF under the first one in a top (roots) down (leaves) view

1796:   Input Parameters:
1797: + sfA - The first PetscSF
1798: - sfB - The second PetscSF

1800:   Output Parameters:
1801: . sfBA - The composite SF

1803:   Level: developer

1805:   Notes:
1806:   Currently, the two SFs must be defined on congruent communicators and they must be true star
1807:   forests, i.e. the same leaf is not connected with different roots.

1809:   sfA's leaf space and sfB's root space might be partially overlapped. The composition builds
1810:   a graph with sfA's roots and sfB's leaves only when there is a path between them. Unconnected
1811:   nodes (roots or leaves) are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a
1812:   Bcast on sfA, then a Bcast on sfB, on connected nodes.

1814: .seealso: PetscSF, PetscSFComposeInverse(), PetscSFGetGraph(), PetscSFSetGraph()
1815: @*/
1816: PetscErrorCode PetscSFCompose(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
1817: {
1818:   PetscErrorCode    ierr;
1819:   const PetscSFNode *remotePointsA,*remotePointsB;
1820:   PetscSFNode       *remotePointsBA=NULL,*reorderedRemotePointsA = NULL,*leafdataB;
1821:   const PetscInt    *localPointsA,*localPointsB;
1822:   PetscInt          *localPointsBA;
1823:   PetscInt          i,numRootsA,numLeavesA,numRootsB,numLeavesB,minleaf,maxleaf,numLeavesBA;
1824:   PetscBool         denseB;

1828:   PetscSFCheckGraphSet(sfA,1);
1830:   PetscSFCheckGraphSet(sfB,2);
1833:   PetscSFCheckLeavesUnique_Private(sfA);
1834:   PetscSFCheckLeavesUnique_Private(sfB);

1836:   PetscSFGetGraph(sfA,&numRootsA,&numLeavesA,&localPointsA,&remotePointsA);
1837:   PetscSFGetGraph(sfB,&numRootsB,&numLeavesB,&localPointsB,&remotePointsB);
1838:   if (localPointsA) {
1839:     PetscMalloc1(numRootsB,&reorderedRemotePointsA);
1840:     for (i=0; i<numRootsB; i++) {
1841:       reorderedRemotePointsA[i].rank = -1;
1842:       reorderedRemotePointsA[i].index = -1;
1843:     }
1844:     for (i=0; i<numLeavesA; i++) {
1845:       if (localPointsA[i] >= numRootsB) continue;
1846:       reorderedRemotePointsA[localPointsA[i]] = remotePointsA[i];
1847:     }
1848:     remotePointsA = reorderedRemotePointsA;
1849:   }
1850:   PetscSFGetLeafRange(sfB,&minleaf,&maxleaf);
1851:   PetscMalloc1(maxleaf-minleaf+1,&leafdataB);
1852:   PetscSFBcastBegin(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf);
1853:   PetscSFBcastEnd(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf);
1854:   PetscFree(reorderedRemotePointsA);

1856:   denseB = (PetscBool)!localPointsB;
1857:   for (i=0,numLeavesBA=0; i<numLeavesB; i++) {
1858:     if (leafdataB[localPointsB ? localPointsB[i]-minleaf : i].rank == -1) denseB = PETSC_FALSE;
1859:     else numLeavesBA++;
1860:   }
1861:   if (denseB) {
1862:     localPointsBA  = NULL;
1863:     remotePointsBA = leafdataB;
1864:   } else {
1865:     PetscMalloc1(numLeavesBA,&localPointsBA);
1866:     PetscMalloc1(numLeavesBA,&remotePointsBA);
1867:     for (i=0,numLeavesBA=0; i<numLeavesB; i++) {
1868:       const PetscInt l = localPointsB ? localPointsB[i] : i;

1870:       if (leafdataB[l-minleaf].rank == -1) continue;
1871:       remotePointsBA[numLeavesBA] = leafdataB[l-minleaf];
1872:       localPointsBA[numLeavesBA] = l;
1873:       numLeavesBA++;
1874:     }
1875:     PetscFree(leafdataB);
1876:   }
1877:   PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA);
1878:   PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER);
1879:   return(0);
1880: }

1882: /*@
1883:   PetscSFComposeInverse - Compose a new PetscSF by putting the inverse of the second SF under the first one

1885:   Input Parameters:
1886: + sfA - The first PetscSF
1887: - sfB - The second PetscSF

1889:   Output Parameters:
1890: . sfBA - The composite SF.

1892:   Level: developer

1894:   Notes:
1895:   Currently, the two SFs must be defined on congruent communicators and they must be true star
1896:   forests, i.e. the same leaf is not connected with different roots. Even more, all roots of the
1897:   second SF must have a degree of 1, i.e., no roots have more than one leaf connected.

1899:   sfA's leaf space and sfB's leaf space might be partially overlapped. The composition builds
1900:   a graph with sfA's roots and sfB's roots only when there is a path between them. Unconnected
1901:   roots are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a Bcast on sfA, then
1902:   a Reduce on sfB, on connected roots.

1904: .seealso: PetscSF, PetscSFCompose(), PetscSFGetGraph(), PetscSFSetGraph(), PetscSFCreateInverseSF()
1905: @*/
1906: PetscErrorCode PetscSFComposeInverse(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
1907: {
1908:   PetscErrorCode    ierr;
1909:   const PetscSFNode *remotePointsA,*remotePointsB;
1910:   PetscSFNode       *remotePointsBA;
1911:   const PetscInt    *localPointsA,*localPointsB;
1912:   PetscSFNode       *reorderedRemotePointsA = NULL;
1913:   PetscInt          i,numRootsA,numLeavesA,numLeavesBA,numRootsB,numLeavesB,minleaf,maxleaf,*localPointsBA;

1917:   PetscSFCheckGraphSet(sfA,1);
1919:   PetscSFCheckGraphSet(sfB,2);
1922:   PetscSFCheckLeavesUnique_Private(sfA);
1923:   PetscSFCheckLeavesUnique_Private(sfB);
1924:   /* TODO: Check roots of sfB have degree of 1 */

1926:   PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);
1927:   PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);
1928:   PetscSFGetLeafRange(sfB, &minleaf, &maxleaf);
1929:   PetscMalloc1(maxleaf - minleaf + 1,&reorderedRemotePointsA);
1930:   for (i=0; i<maxleaf - minleaf + 1; i++) {
1931:     reorderedRemotePointsA[i].rank = -1;
1932:     reorderedRemotePointsA[i].index = -1;
1933:   }
1934:   if (localPointsA) {
1935:     for (i=0; i<numLeavesA; i++) {
1936:       if (localPointsA[i] > maxleaf || localPointsA[i] < minleaf) continue;
1937:       reorderedRemotePointsA[localPointsA[i] - minleaf] = remotePointsA[i];
1938:     }
1939:   } else {
1940:     for (i=0; i<numLeavesA; i++) {
1941:       if (i > maxleaf || i < minleaf) continue;
1942:       reorderedRemotePointsA[i - minleaf] = remotePointsA[i];
1943:     }
1944:   }

1946:   PetscMalloc1(numRootsB,&localPointsBA);
1947:   PetscMalloc1(numRootsB,&remotePointsBA);
1948:   for (i=0; i<numRootsB; i++) {
1949:     remotePointsBA[i].rank = -1;
1950:     remotePointsBA[i].index = -1;
1951:   }

1953:   /* Once we implement the TODO above (check all roots of sfB have degree of 1), we can replace the MPI_MAXLOC
1954:      with MPIU_REPLACE. In that case, MPI_MAXLOC and MPIU_REPLACE have the same effect.
1955:    */
1956:   PetscSFReduceBegin(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,MPI_MAXLOC);
1957:   PetscSFReduceEnd(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,MPI_MAXLOC);
1958:   PetscFree(reorderedRemotePointsA);
1959:   for (i=0,numLeavesBA=0; i<numRootsB; i++) {
1960:     if (remotePointsBA[i].rank == -1) continue;
1961:     remotePointsBA[numLeavesBA].rank = remotePointsBA[i].rank;
1962:     remotePointsBA[numLeavesBA].index = remotePointsBA[i].index;
1963:     localPointsBA[numLeavesBA] = i;
1964:     numLeavesBA++;
1965:   }
1966:   PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA);
1967:   PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER);
1968:   return(0);
1969: }

1971: /*
1972:   PetscSFCreateLocalSF_Private - Creates a local PetscSF that only has intra-process edges of the global PetscSF

1974:   Input Parameters:
1975: . sf - The global PetscSF

1977:   Output Parameters:
1978: . out - The local PetscSF
1979:  */
1980: PetscErrorCode PetscSFCreateLocalSF_Private(PetscSF sf,PetscSF *out)
1981: {
1982:   MPI_Comm           comm;
1983:   PetscMPIInt        myrank;
1984:   const PetscInt     *ilocal;
1985:   const PetscSFNode  *iremote;
1986:   PetscInt           i,j,nroots,nleaves,lnleaves,*lilocal;
1987:   PetscSFNode        *liremote;
1988:   PetscSF            lsf;
1989:   PetscErrorCode     ierr;

1993:   if (sf->ops->CreateLocalSF) {
1994:     (*sf->ops->CreateLocalSF)(sf,out);
1995:   } else {
1996:     /* Could use PetscSFCreateEmbeddedLeafSF, but since we know the comm is PETSC_COMM_SELF, we can make it fast */
1997:     PetscObjectGetComm((PetscObject)sf,&comm);
1998:     MPI_Comm_rank(comm,&myrank);

2000:     /* Find out local edges and build a local SF */
2001:     PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);
2002:     for (i=lnleaves=0; i<nleaves; i++) {if (iremote[i].rank == (PetscInt)myrank) lnleaves++;}
2003:     PetscMalloc1(lnleaves,&lilocal);
2004:     PetscMalloc1(lnleaves,&liremote);

2006:     for (i=j=0; i<nleaves; i++) {
2007:       if (iremote[i].rank == (PetscInt)myrank) {
2008:         lilocal[j]        = ilocal? ilocal[i] : i; /* ilocal=NULL for contiguous storage */
2009:         liremote[j].rank  = 0; /* rank in PETSC_COMM_SELF */
2010:         liremote[j].index = iremote[i].index;
2011:         j++;
2012:       }
2013:     }
2014:     PetscSFCreate(PETSC_COMM_SELF,&lsf);
2015:     PetscSFSetGraph(lsf,nroots,lnleaves,lilocal,PETSC_OWN_POINTER,liremote,PETSC_OWN_POINTER);
2016:     PetscSFSetUp(lsf);
2017:     *out = lsf;
2018:   }
2019:   return(0);
2020: }

2022: /* Similar to PetscSFBcast, but only Bcast to leaves on rank 0 */
2023: PetscErrorCode PetscSFBcastToZero_Private(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
2024: {
2025:   PetscErrorCode     ierr;
2026:   PetscMemType       rootmtype,leafmtype;

2030:   PetscSFSetUp(sf);
2031:   PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);
2032:   PetscGetMemType(rootdata,&rootmtype);
2033:   PetscGetMemType(leafdata,&leafmtype);
2034:   if (sf->ops->BcastToZero) {
2035:     (*sf->ops->BcastToZero)(sf,unit,rootmtype,rootdata,leafmtype,leafdata);
2036:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFBcastToZero_Private is not supported on this SF type");
2037:   PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);
2038:   return(0);
2039: }