Actual source code: sf.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1:  #include <petsc/private/sfimpl.h>
  2:  #include <petscctable.h>

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

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

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

 18:    Collective on MPI_Comm

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

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

 26:    Level: intermediate

 28: .seealso: PetscSFSetGraph(), PetscSFDestroy()
 29: @*/
 30: PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf)
 31: {
 33:   PetscSF        b;

 37:   PetscSFInitializePackage();

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

 41:   b->nroots    = -1;
 42:   b->nleaves   = -1;
 43:   b->minleaf   = PETSC_MAX_INT;
 44:   b->maxleaf   = PETSC_MIN_INT;
 45:   b->nranks    = -1;
 46:   b->rankorder = PETSC_TRUE;
 47:   b->ingroup   = MPI_GROUP_NULL;
 48:   b->outgroup  = MPI_GROUP_NULL;
 49:   b->graphset  = PETSC_FALSE;

 51:   *sf = b;
 52:   return(0);
 53: }

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

 58:    Collective

 60:    Input Arguments:
 61: .  sf - star forest

 63:    Level: advanced

 65: .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy()
 66: @*/
 67: PetscErrorCode PetscSFReset(PetscSF sf)
 68: {

 73:   if (sf->ops->Reset) {(*sf->ops->Reset)(sf);}
 74:   sf->nroots   = -1;
 75:   sf->nleaves  = -1;
 76:   sf->minleaf  = PETSC_MAX_INT;
 77:   sf->maxleaf  = PETSC_MIN_INT;
 78:   sf->mine     = NULL;
 79:   sf->remote   = NULL;
 80:   sf->graphset = PETSC_FALSE;
 81:   PetscFree(sf->mine_alloc);
 82:   PetscFree(sf->remote_alloc);
 83:   sf->nranks = -1;
 84:   PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);
 85:   sf->degreeknown = PETSC_FALSE;
 86:   PetscFree(sf->degree);
 87:   if (sf->ingroup  != MPI_GROUP_NULL) {MPI_Group_free(&sf->ingroup);}
 88:   if (sf->outgroup != MPI_GROUP_NULL) {MPI_Group_free(&sf->outgroup);}
 89:   PetscSFDestroy(&sf->multi);
 90:   sf->setupcalled = PETSC_FALSE;
 91:   return(0);
 92: }

 94: /*@C
 95:    PetscSFSetType - Set the PetscSF communication implementation

 97:    Collective on PetscSF

 99:    Input Parameters:
100: +  sf - the PetscSF context
101: -  type - a known method

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

107:    Notes:
108:    See "include/petscsf.h" for available methods (for instance)
109: +    PETSCSFWINDOW - MPI-2/3 one-sided
110: -    PETSCSFBASIC - basic implementation using MPI-1 two-sided

112:   Level: intermediate

114: .keywords: PetscSF, set, type

116: .seealso: PetscSFType, PetscSFCreate()
117: @*/
118: PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type)
119: {
120:   PetscErrorCode ierr,(*r)(PetscSF);
121:   PetscBool      match;


127:   PetscObjectTypeCompare((PetscObject)sf,type,&match);
128:   if (match) return(0);

130:   PetscFunctionListFind(PetscSFList,type,&r);
131:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type);
132:   /* Destroy the previous PetscSF implementation context */
133:   if (sf->ops->Destroy) {(*(sf)->ops->Destroy)(sf);}
134:   PetscMemzero(sf->ops,sizeof(*sf->ops));
135:   PetscObjectChangeTypeName((PetscObject)sf,type);
136:   (*r)(sf);
137:   return(0);
138: }

140: /*@C
141:   PetscSFGetType - Get the PetscSF communication implementation

143:   Not Collective

145:   Input Parameter:
146: . sf  - the PetscSF context

148:   Output Parameter:
149: . type - the PetscSF type name

151:   Level: intermediate

153: .keywords: PetscSF, get, type
154: .seealso: PetscSFSetType(), PetscSFCreate()
155: @*/
156: PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type)
157: {
161:   *type = ((PetscObject)sf)->type_name;
162:   return(0);
163: }

165: /*@
166:    PetscSFDestroy - destroy star forest

168:    Collective

170:    Input Arguments:
171: .  sf - address of star forest

173:    Level: intermediate

175: .seealso: PetscSFCreate(), PetscSFReset()
176: @*/
177: PetscErrorCode PetscSFDestroy(PetscSF *sf)
178: {

182:   if (!*sf) return(0);
184:   if (--((PetscObject)(*sf))->refct > 0) {*sf = NULL; return(0);}
185:   PetscSFReset(*sf);
186:   if ((*sf)->ops->Destroy) {(*(*sf)->ops->Destroy)(*sf);}
187:   PetscHeaderDestroy(sf);
188:   return(0);
189: }

191: /*@
192:    PetscSFSetUp - set up communication structures

194:    Collective

196:    Input Arguments:
197: .  sf - star forest communication object

199:    Level: beginner

201: .seealso: PetscSFSetFromOptions(), PetscSFSetType()
202: @*/
203: PetscErrorCode PetscSFSetUp(PetscSF sf)
204: {

209:   PetscSFCheckGraphSet(sf,1);
210:   if (sf->setupcalled) return(0);
211:   if (!((PetscObject)sf)->type_name) {PetscSFSetType(sf,PETSCSFBASIC);}
212:   PetscLogEventBegin(PETSCSF_SetUp,sf,0,0,0);
213:   if (sf->ops->SetUp) {(*sf->ops->SetUp)(sf);}
214:   PetscLogEventEnd(PETSCSF_SetUp,sf,0,0,0);
215:   sf->setupcalled = PETSC_TRUE;
216:   return(0);
217: }

219: /*@
220:    PetscSFSetFromOptions - set PetscSF options using the options database

222:    Logically Collective

224:    Input Arguments:
225: .  sf - star forest

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

231:    Level: intermediate

233: .keywords: KSP, set, from, options, database

235: .seealso: PetscSFWindowSetSyncType()
236: @*/
237: PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
238: {
239:   PetscSFType    deft;
240:   char           type[256];
242:   PetscBool      flg;

246:   PetscObjectOptionsBegin((PetscObject)sf);
247:   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
248:   PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,sizeof(type),&flg);
249:   PetscSFSetType(sf,flg ? type : deft);
250:   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);
251:   if (sf->ops->SetFromOptions) {(*sf->ops->SetFromOptions)(PetscOptionsObject,sf);}
252:   PetscOptionsEnd();
253:   return(0);
254: }

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

259:    Logically Collective

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

265:    Level: advanced

267: .seealso: PetscSFGatherBegin(), PetscSFScatterBegin()
268: @*/
269: PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg)
270: {

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

280: /*@
281:    PetscSFSetGraph - Set a parallel star forest

283:    Collective

285:    Input Arguments:
286: +  sf - star forest
287: .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
288: .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
289: .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage
290: .  localmode - copy mode for ilocal
291: .  iremote - remote locations of root vertices for each leaf on the current process
292: -  remotemode - copy mode for iremote

294:    Level: intermediate

296:    Notes:
297:     In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode

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

303: .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph()
304: @*/
305: PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode)
306: {

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

316:   PetscSFReset(sf);
317:   PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);

319:   sf->nroots  = nroots;
320:   sf->nleaves = nleaves;

322:   if (nleaves && ilocal) {
323:     PetscInt i;
324:     PetscInt minleaf = PETSC_MAX_INT;
325:     PetscInt maxleaf = PETSC_MIN_INT;
326:     int      contiguous = 1;
327:     for (i=0; i<nleaves; i++) {
328:       minleaf = PetscMin(minleaf,ilocal[i]);
329:       maxleaf = PetscMax(maxleaf,ilocal[i]);
330:       contiguous &= (ilocal[i] == i);
331:     }
332:     sf->minleaf = minleaf;
333:     sf->maxleaf = maxleaf;
334:     if (contiguous) {
335:       if (localmode == PETSC_OWN_POINTER) {
336:         PetscFree(ilocal);
337:       }
338:       ilocal = NULL;
339:     }
340:   } else {
341:     sf->minleaf = 0;
342:     sf->maxleaf = nleaves - 1;
343:   }

345:   if (ilocal) {
346:     switch (localmode) {
347:     case PETSC_COPY_VALUES:
348:       PetscMalloc1(nleaves,&sf->mine_alloc);
349:       PetscMemcpy(sf->mine_alloc,ilocal,nleaves*sizeof(*ilocal));
350:       sf->mine = sf->mine_alloc;
351:       break;
352:     case PETSC_OWN_POINTER:
353:       sf->mine_alloc = (PetscInt*)ilocal;
354:       sf->mine       = sf->mine_alloc;
355:       break;
356:     case PETSC_USE_POINTER:
357:       sf->mine_alloc = NULL;
358:       sf->mine       = (PetscInt*)ilocal;
359:       break;
360:     default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode");
361:     }
362:   }

364:   switch (remotemode) {
365:   case PETSC_COPY_VALUES:
366:     PetscMalloc1(nleaves,&sf->remote_alloc);
367:     PetscMemcpy(sf->remote_alloc,iremote,nleaves*sizeof(*iremote));
368:     sf->remote = sf->remote_alloc;
369:     break;
370:   case PETSC_OWN_POINTER:
371:     sf->remote_alloc = (PetscSFNode*)iremote;
372:     sf->remote       = sf->remote_alloc;
373:     break;
374:   case PETSC_USE_POINTER:
375:     sf->remote_alloc = NULL;
376:     sf->remote       = (PetscSFNode*)iremote;
377:     break;
378:   default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode");
379:   }

381:   PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);
382:   sf->graphset = PETSC_TRUE;
383:   return(0);
384: }

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

389:    Collective

391:    Input Arguments:
392: .  sf - star forest to invert

394:    Output Arguments:
395: .  isf - inverse of sf

397:    Level: advanced

399:    Notes:
400:    All roots must have degree 1.

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

404: .seealso: PetscSFSetGraph()
405: @*/
406: PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
407: {
409:   PetscMPIInt    rank;
410:   PetscInt       i,nroots,nleaves,maxlocal,count,*newilocal;
411:   const PetscInt *ilocal;
412:   PetscSFNode    *roots,*leaves;

416:   PetscSFCheckGraphSet(sf,1);

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

422:   MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);
423:   PetscMalloc2(nroots,&roots,maxlocal,&leaves);
424:   for (i=0; i<maxlocal; i++) {
425:     leaves[i].rank  = rank;
426:     leaves[i].index = i;
427:   }
428:   for (i=0; i <nroots; i++) {
429:     roots[i].rank  = -1;
430:     roots[i].index = -1;
431:   }
432:   PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);
433:   PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);

435:   /* Check whether our leaves are sparse */
436:   for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
437:   if (count == nroots) newilocal = NULL;
438:   else {                        /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
439:     PetscMalloc1(count,&newilocal);
440:     for (i=0,count=0; i<nroots; i++) {
441:       if (roots[i].rank >= 0) {
442:         newilocal[count]   = i;
443:         roots[count].rank  = roots[i].rank;
444:         roots[count].index = roots[i].index;
445:         count++;
446:       }
447:     }
448:   }

450:   PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);
451:   PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);
452:   PetscFree2(roots,leaves);
453:   return(0);
454: }

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

459:    Collective

461:    Input Arguments:
462: +  sf - communication object to duplicate
463: -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)

465:    Output Arguments:
466: .  newsf - new communication object

468:    Level: beginner

470: .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph()
471: @*/
472: PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
473: {
474:   PetscSFType    type;

481:   PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);
482:   PetscSFGetType(sf,&type);
483:   if (type) {PetscSFSetType(*newsf,type);}
484:   if (opt == PETSCSF_DUPLICATE_GRAPH) {
485:     PetscInt          nroots,nleaves;
486:     const PetscInt    *ilocal;
487:     const PetscSFNode *iremote;
488:     PetscSFCheckGraphSet(sf,1);
489:     PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);
490:     PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);
491:   }
492:   if (sf->ops->Duplicate) {(*sf->ops->Duplicate)(sf,opt,*newsf);}
493:   return(0);
494: }

496: /*@C
497:    PetscSFGetGraph - Get the graph specifying a parallel star forest

499:    Not Collective

501:    Input Arguments:
502: .  sf - star forest

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

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

513:    Level: intermediate

515: .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph()
516: @*/
517: PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
518: {

522:   if (nroots) *nroots = sf->nroots;
523:   if (nleaves) *nleaves = sf->nleaves;
524:   if (ilocal) *ilocal = sf->mine;
525:   if (iremote) *iremote = sf->remote;
526:   return(0);
527: }

529: /*@
530:    PetscSFGetLeafRange - Get the active leaf ranges

532:    Not Collective

534:    Input Arguments:
535: .  sf - star forest

537:    Output Arguments:
538: +  minleaf - minimum active leaf on this process
539: -  maxleaf - maximum active leaf on this process

541:    Level: developer

543: .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
544: @*/
545: PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
546: {

550:   PetscSFCheckGraphSet(sf,1);
551:   if (minleaf) *minleaf = sf->minleaf;
552:   if (maxleaf) *maxleaf = sf->maxleaf;
553:   return(0);
554: }

556: /*@C
557:    PetscSFView - view a star forest

559:    Collective

561:    Input Arguments:
562: +  sf - star forest
563: -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD

565:    Level: beginner

567: .seealso: PetscSFCreate(), PetscSFSetGraph()
568: @*/
569: PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
570: {
571:   PetscErrorCode    ierr;
572:   PetscBool         iascii;
573:   PetscViewerFormat format;

577:   if (!viewer) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);}
580:   if (sf->graphset) {PetscSFSetUp(sf);}
581:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
582:   if (iascii) {
583:     PetscMPIInt rank;
584:     PetscInt    ii,i,j;

586:     PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);
587:     PetscViewerASCIIPushTab(viewer);
588:     if (sf->ops->View) {(*sf->ops->View)(sf,viewer);}
589:     if (!sf->graphset) {
590:       PetscViewerASCIIPrintf(viewer,"PetscSFSetGraph() has not been called yet\n");
591:       PetscViewerASCIIPopTab(viewer);
592:       return(0);
593:     }
594:     MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);
595:     PetscViewerASCIIPushSynchronized(viewer);
596:     PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);
597:     for (i=0; i<sf->nleaves; i++) {
598:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- (%D,%D)\n",rank,sf->mine ? sf->mine[i] : i,sf->remote[i].rank,sf->remote[i].index);
599:     }
600:     PetscViewerFlush(viewer);
601:     PetscViewerGetFormat(viewer,&format);
602:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
603:       PetscMPIInt *tmpranks,*perm;
604:       PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm);
605:       PetscMemcpy(tmpranks,sf->ranks,sf->nranks*sizeof(tmpranks[0]));
606:       for (i=0; i<sf->nranks; i++) perm[i] = i;
607:       PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm);
608:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);
609:       for (ii=0; ii<sf->nranks; ii++) {
610:         i = perm[ii];
611:         PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);
612:         for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
613:           PetscViewerASCIISynchronizedPrintf(viewer,"[%d]    %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);
614:         }
615:       }
616:       PetscFree2(tmpranks,perm);
617:     }
618:     PetscViewerFlush(viewer);
619:     PetscViewerASCIIPopSynchronized(viewer);
620:     PetscViewerASCIIPopTab(viewer);
621:   }
622:   return(0);
623: }

625: /*@C
626:    PetscSFGetRanks - Get ranks and number of vertices referenced by leaves on this process

628:    Not Collective

630:    Input Arguments:
631: .  sf - star forest

633:    Output Arguments:
634: +  nranks - number of ranks referenced by local part
635: .  ranks - array of ranks
636: .  roffset - offset in rmine/rremote for each rank (length nranks+1)
637: .  rmine - concatenated array holding local indices referencing each remote rank
638: -  rremote - concatenated array holding remote indices referenced for each remote rank

640:    Level: developer

642: .seealso: PetscSFSetGraph()
643: @*/
644: PetscErrorCode PetscSFGetRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
645: {

649:   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
650:   if (nranks)  *nranks  = sf->nranks;
651:   if (ranks)   *ranks   = sf->ranks;
652:   if (roffset) *roffset = sf->roffset;
653:   if (rmine)   *rmine   = sf->rmine;
654:   if (rremote) *rremote = sf->rremote;
655:   return(0);
656: }

658: /*@C
659:    PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process

661:    Not Collective

663:    Input Arguments:
664: .  sf - star forest

666:    Output Arguments:
667: +  niranks - number of leaf ranks referencing roots on this process
668: .  iranks - array of ranks
669: .  ioffset - offset in irootloc for each rank (length niranks+1)
670: -  irootloc - concatenated array holding local indices of roots referenced by each leaf rank

672:    Level: developer

674: .seealso: PetscSFGetRanks()
675: @*/
676: PetscErrorCode PetscSFGetLeafRanks(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
677: {

682:   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
683:   if (sf->ops->GetLeafRanks) {
684:     (sf->ops->GetLeafRanks)(sf,niranks,iranks,ioffset,irootloc);
685:   } else {
686:     PetscSFType type;
687:     PetscSFGetType(sf,&type);
688:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type);
689:   }
690:   return(0);
691: }

693: static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) {
694:   PetscInt i;
695:   for (i=0; i<n; i++) {
696:     if (needle == list[i]) return PETSC_TRUE;
697:   }
698:   return PETSC_FALSE;
699: }

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

704:    Collective

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

710:    Level: developer

712: .seealso: PetscSFGetRanks()
713: @*/
714: PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup)
715: {
716:   PetscErrorCode     ierr;
717:   PetscTable         table;
718:   PetscTablePosition pos;
719:   PetscMPIInt        size,groupsize,*groupranks;
720:   PetscInt           i,*rcount,*ranks;

724:   PetscSFCheckGraphSet(sf,1);
725:   MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);
726:   PetscTableCreate(10,size,&table);
727:   for (i=0; i<sf->nleaves; i++) {
728:     /* Log 1-based rank */
729:     PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);
730:   }
731:   PetscTableGetCount(table,&sf->nranks);
732:   PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);
733:   PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);
734:   PetscTableGetHeadPosition(table,&pos);
735:   for (i=0; i<sf->nranks; i++) {
736:     PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);
737:     ranks[i]--;             /* Convert back to 0-based */
738:   }
739:   PetscTableDestroy(&table);

741:   /* We expect that dgroup is reliably "small" while nranks could be large */
742:   {
743:     MPI_Group group = MPI_GROUP_NULL;
744:     PetscMPIInt *dgroupranks;
745:     MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);
746:     MPI_Group_size(dgroup,&groupsize);
747:     PetscMalloc1(groupsize,&dgroupranks);
748:     PetscMalloc1(groupsize,&groupranks);
749:     for (i=0; i<groupsize; i++) dgroupranks[i] = i;
750:     if (groupsize) {MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);}
751:     MPI_Group_free(&group);
752:     PetscFree(dgroupranks);
753:   }

755:   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
756:   for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i; ) {
757:     for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */
758:       if (InList(ranks[i],groupsize,groupranks)) break;
759:     }
760:     for ( ; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
761:       if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break;
762:     }
763:     if (sf->ndranks < i) {                         /* Swap ranks[sf->ndranks] with ranks[i] */
764:       PetscInt    tmprank,tmpcount;
765:       tmprank = ranks[i];
766:       tmpcount = rcount[i];
767:       ranks[i] = ranks[sf->ndranks];
768:       rcount[i] = rcount[sf->ndranks];
769:       ranks[sf->ndranks] = tmprank;
770:       rcount[sf->ndranks] = tmpcount;
771:       sf->ndranks++;
772:     }
773:   }
774:   PetscFree(groupranks);
775:   PetscSortIntWithArray(sf->ndranks,ranks,rcount);
776:   PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);
777:   sf->roffset[0] = 0;
778:   for (i=0; i<sf->nranks; i++) {
779:     PetscMPIIntCast(ranks[i],sf->ranks+i);
780:     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
781:     rcount[i]        = 0;
782:   }
783:   for (i=0; i<sf->nleaves; i++) {
784:     PetscInt irank;
785:     /* Search for index of iremote[i].rank in sf->ranks */
786:     PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);
787:     if (irank < 0) {
788:       PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);
789:       if (irank >= 0) irank += sf->ndranks;
790:     }
791:     if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank);
792:     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
793:     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
794:     rcount[irank]++;
795:   }
796:   PetscFree2(rcount,ranks);
797:   return(0);
798: }

800: /*@C
801:    PetscSFGetGroups - gets incoming and outgoing process groups

803:    Collective

805:    Input Argument:
806: .  sf - star forest

808:    Output Arguments:
809: +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
810: -  outgoing - group of destination processes for outgoing edges (roots that I reference)

812:    Level: developer

814: .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
815: @*/
816: PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
817: {
819:   MPI_Group      group = MPI_GROUP_NULL;

822:   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining groups");
823:   if (sf->ingroup == MPI_GROUP_NULL) {
824:     PetscInt       i;
825:     const PetscInt *indegree;
826:     PetscMPIInt    rank,*outranks,*inranks;
827:     PetscSFNode    *remote;
828:     PetscSF        bgcount;

830:     /* Compute the number of incoming ranks */
831:     PetscMalloc1(sf->nranks,&remote);
832:     for (i=0; i<sf->nranks; i++) {
833:       remote[i].rank  = sf->ranks[i];
834:       remote[i].index = 0;
835:     }
836:     PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);
837:     PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);
838:     PetscSFComputeDegreeBegin(bgcount,&indegree);
839:     PetscSFComputeDegreeEnd(bgcount,&indegree);

841:     /* Enumerate the incoming ranks */
842:     PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);
843:     MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);
844:     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
845:     PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);
846:     PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);
847:     MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);
848:     MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);
849:     MPI_Group_free(&group);
850:     PetscFree2(inranks,outranks);
851:     PetscSFDestroy(&bgcount);
852:   }
853:   *incoming = sf->ingroup;

855:   if (sf->outgroup == MPI_GROUP_NULL) {
856:     MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);
857:     MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);
858:     MPI_Group_free(&group);
859:   }
860:   *outgoing = sf->outgroup;
861:   return(0);
862: }

864: /*@
865:    PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters

867:    Collective

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

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

875:    Level: developer

877:    Notes:

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

883: .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin(), PetscSFComputeMultiRootOriginalNumbering()
884: @*/
885: PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
886: {

892:   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
893:     PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);
894:     *multi = sf->multi;
895:     return(0);
896:   }
897:   if (!sf->multi) {
898:     const PetscInt *indegree;
899:     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
900:     PetscSFNode    *remote;
901:     maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
902:     PetscSFComputeDegreeBegin(sf,&indegree);
903:     PetscSFComputeDegreeEnd(sf,&indegree);
904:     PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);
905:     inoffset[0] = 0;
906:     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
907:     for (i=0; i<maxlocal; i++) outones[i] = 1;
908:     PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);
909:     PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);
910:     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
911: #if 0
912: #if defined(PETSC_USE_DEBUG)                                 /* Check that the expected number of increments occurred */
913:     for (i=0; i<sf->nroots; i++) {
914:       if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
915:     }
916: #endif
917: #endif
918:     PetscMalloc1(sf->nleaves,&remote);
919:     for (i=0; i<sf->nleaves; i++) {
920:       remote[i].rank  = sf->remote[i].rank;
921:       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
922:     }
923:     PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);
924:     PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);
925:     if (sf->rankorder) {        /* Sort the ranks */
926:       PetscMPIInt rank;
927:       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
928:       PetscSFNode *newremote;
929:       MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);
930:       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
931:       PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);
932:       for (i=0; i<maxlocal; i++) outranks[i] = rank;
933:       PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);
934:       PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);
935:       /* Sort the incoming ranks at each vertex, build the inverse map */
936:       for (i=0; i<sf->nroots; i++) {
937:         PetscInt j;
938:         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
939:         PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);
940:         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
941:       }
942:       PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);
943:       PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);
944:       PetscMalloc1(sf->nleaves,&newremote);
945:       for (i=0; i<sf->nleaves; i++) {
946:         newremote[i].rank  = sf->remote[i].rank;
947:         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
948:       }
949:       PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);
950:       PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);
951:     }
952:     PetscFree3(inoffset,outones,outoffset);
953:   }
954:   *multi = sf->multi;
955:   return(0);
956: }

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

961:    Collective

963:    Input Arguments:
964: +  sf - original star forest
965: .  nselected  - number of selected roots on this process
966: -  selected   - indices of the selected roots on this process

968:    Output Arguments:
969: .  newsf - new star forest

971:    Level: advanced

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

977: .seealso: PetscSFSetGraph(), PetscSFGetGraph()
978: @*/
979: PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
980: {
981:   PetscInt          *rootdata,*leafdata,*new_ilocal;
982:   PetscSFNode       *new_iremote;
983:   const PetscInt    *ilocal;
984:   const PetscSFNode *iremote;
985:   PetscInt          nleaves,nroots,n,i,new_nleaves = 0;
986:   PetscSF           tmpsf;
987:   PetscErrorCode    ierr;

991:   PetscSFCheckGraphSet(sf,1);

995:   /* Find out which leaves (not leaf data items) are still connected to roots in the embedded sf */
996:   PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);
997:   PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&tmpsf);
998:   PetscSFSetGraph(tmpsf,nroots,nleaves,NULL/*contiguous*/,PETSC_USE_POINTER,iremote,PETSC_USE_POINTER);
999:   PetscCalloc2(nroots,&rootdata,nleaves,&leafdata);
1000:   for (i=0; i<nselected; ++i) {
1001:     if (selected[i] < 0 || selected[i] >= nroots) SETERRQ2(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Root index %D is not in [0,%D)",selected[i],nroots);
1002:     rootdata[selected[i]] = 1;
1003:   }

1005:   PetscSFBcastBegin(tmpsf,MPIU_INT,rootdata,leafdata);
1006:   PetscSFBcastEnd(tmpsf,MPIU_INT,rootdata,leafdata);
1007:   PetscSFDestroy(&tmpsf);

1009:   /* Build newsf with leaves that are still connected */
1010:   for (i = 0; i < nleaves; ++i) new_nleaves += leafdata[i];
1011:   PetscMalloc1(new_nleaves,&new_ilocal);
1012:   PetscMalloc1(new_nleaves,&new_iremote);
1013:   for (i = 0, n = 0; i < nleaves; ++i) {
1014:     if (leafdata[i]) {
1015:       new_ilocal[n]        = sf->mine ? sf->mine[i] : i;
1016:       new_iremote[n].rank  = sf->remote[i].rank;
1017:       new_iremote[n].index = sf->remote[i].index;
1018:       ++n;
1019:     }
1020:   }
1021:   if (n != new_nleaves) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"There is a size mismatch in the SF embedding, %d != %d",n,new_nleaves);
1022:   PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,newsf);
1023:   PetscSFSetGraph(*newsf,nroots,new_nleaves,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);
1024:   PetscFree2(rootdata,leafdata);
1025:   return(0);
1026: }

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

1031:   Collective

1033:   Input Arguments:
1034: + sf - original star forest
1035: . nleaves - number of leaves to select on this process
1036: - selected - selected leaves on this process

1038:   Output Arguments:
1039: .  newsf - new star forest

1041:   Level: advanced

1043: .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph()
1044: @*/
1045: PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nleaves, const PetscInt *selected, PetscSF *newsf)
1046: {
1047:   PetscSFNode   *iremote;
1048:   PetscInt      *ilocal;
1049:   PetscInt       i;

1054:   PetscSFCheckGraphSet(sf, 1);
1057:   PetscMalloc1(nleaves, &ilocal);
1058:   PetscMalloc1(nleaves, &iremote);
1059:   for (i = 0; i < nleaves; ++i) {
1060:     const PetscInt l = selected[i];

1062:     ilocal[i]        = sf->mine ? sf->mine[l] : l;
1063:     iremote[i].rank  = sf->remote[l].rank;
1064:     iremote[i].index = sf->remote[l].index;
1065:   }
1066:   PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, newsf);
1067:   PetscSFSetGraph(*newsf, sf->nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);
1068:   return(0);
1069: }

1071: /*@C
1072:    PetscSFBcastBegin - begin pointwise broadcast to be concluded with call to PetscSFBcastEnd()

1074:    Collective on PetscSF

1076:    Input Arguments:
1077: +  sf - star forest on which to communicate
1078: .  unit - data type associated with each node
1079: -  rootdata - buffer to broadcast

1081:    Output Arguments:
1082: .  leafdata - buffer to update with values from each leaf's respective root

1084:    Level: intermediate

1086: .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFView(), PetscSFBcastEnd(), PetscSFReduceBegin()
1087: @*/
1088: PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
1089: {

1094:   PetscSFSetUp(sf);
1095:   PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);
1096:   (*sf->ops->BcastBegin)(sf,unit,rootdata,leafdata);
1097:   PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);
1098:   return(0);
1099: }

1101: /*@C
1102:    PetscSFBcastEnd - end a broadcast operation started with PetscSFBcastBegin()

1104:    Collective

1106:    Input Arguments:
1107: +  sf - star forest
1108: .  unit - data type
1109: -  rootdata - buffer to broadcast

1111:    Output Arguments:
1112: .  leafdata - buffer to update with values from each leaf's respective root

1114:    Level: intermediate

1116: .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
1117: @*/
1118: PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
1119: {

1124:   PetscSFSetUp(sf);
1125:   PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0);
1126:   (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata);
1127:   PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0);
1128:   return(0);
1129: }

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

1134:    Collective on PetscSF

1136:    Input Arguments:
1137: +  sf - star forest on which to communicate
1138: .  unit - data type associated with each node
1139: .  rootdata - buffer to broadcast
1140: -  op - operation to use for reduction

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

1145:    Level: intermediate

1147: .seealso: PetscSFBcastAndOpEnd(), PetscSFBcastBegin(), PetscSFBcastEnd()
1148: @*/
1149: PetscErrorCode PetscSFBcastAndOpBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
1150: {

1155:   PetscSFSetUp(sf);
1156:   PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);
1157:   (*sf->ops->BcastAndOpBegin)(sf,unit,rootdata,leafdata,op);
1158:   PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);
1159:   return(0);
1160: }

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

1165:    Collective

1167:    Input Arguments:
1168: +  sf - star forest
1169: .  unit - data type
1170: .  rootdata - buffer to broadcast
1171: -  op - operation to use for reduction

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

1176:    Level: intermediate

1178: .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
1179: @*/
1180: PetscErrorCode PetscSFBcastAndOpEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
1181: {

1186:   PetscSFSetUp(sf);
1187:   PetscLogEventBegin(PETSCSF_BcastAndOpEnd,sf,0,0,0);
1188:   (*sf->ops->BcastAndOpEnd)(sf,unit,rootdata,leafdata,op);
1189:   PetscLogEventEnd(PETSCSF_BcastAndOpEnd,sf,0,0,0);
1190:   return(0);
1191: }

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

1196:    Collective

1198:    Input Arguments:
1199: +  sf - star forest
1200: .  unit - data type
1201: .  leafdata - values to reduce
1202: -  op - reduction operation

1204:    Output Arguments:
1205: .  rootdata - result of reduction of values from all leaves of each root

1207:    Level: intermediate

1209: .seealso: PetscSFBcastBegin()
1210: @*/
1211: PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
1212: {

1217:   PetscSFSetUp(sf);
1218:   PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);
1219:   (sf->ops->ReduceBegin)(sf,unit,leafdata,rootdata,op);
1220:   PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);
1221:   return(0);
1222: }

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

1227:    Collective

1229:    Input Arguments:
1230: +  sf - star forest
1231: .  unit - data type
1232: .  leafdata - values to reduce
1233: -  op - reduction operation

1235:    Output Arguments:
1236: .  rootdata - result of reduction of values from all leaves of each root

1238:    Level: intermediate

1240: .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
1241: @*/
1242: PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
1243: {

1248:   PetscSFSetUp(sf);
1249:   PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);
1250:   (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);
1251:   PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);
1252:   return(0);
1253: }

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

1258:    Collective

1260:    Input Arguments:
1261: .  sf - star forest

1263:    Output Arguments:
1264: .  degree - degree of each root vertex

1266:    Level: advanced

1268: .seealso: PetscSFGatherBegin()
1269: @*/
1270: PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
1271: {

1276:   PetscSFCheckGraphSet(sf,1);
1278:   if (!sf->degreeknown) {
1279:     PetscInt i, nroots = sf->nroots, maxlocal = sf->maxleaf+1;  /* TODO: We should use PetscSFGetLeafRange() */
1280:     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
1281:     PetscMalloc1(nroots,&sf->degree);
1282:     PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
1283:     for (i=0; i<nroots; i++) sf->degree[i] = 0;
1284:     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
1285:     PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);
1286:   }
1287:   *degree = NULL;
1288:   return(0);
1289: }

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

1294:    Collective

1296:    Input Arguments:
1297: .  sf - star forest

1299:    Output Arguments:
1300: .  degree - degree of each root vertex

1302:    Level: developer

1304: .seealso:
1305: @*/
1306: PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
1307: {

1312:   PetscSFCheckGraphSet(sf,1);
1314:   if (!sf->degreeknown) {
1315:     if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
1316:     PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp,sf->degree,MPI_SUM);
1317:     PetscFree(sf->degreetmp);
1318:     sf->degreeknown = PETSC_TRUE;
1319:   }
1320:   *degree = sf->degree;
1321:   return(0);
1322: }


1325: /*@C
1326:    PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()). Each multi-root is assigned index of its original root.

1328:    Collective

1330:    Input Arguments:
1331: +  sf - star forest
1332: -  degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd()

1334:    Output Arguments:
1335: .  mRootsOrigNumbering - original indices of multi-roots; length of the array is equal to the number of multi-roots (roots of multi-SF)

1337:    Level: developer

1339: .seealso: PetscSFComputeDegreeBegin(), PetscSFComputeDegreeEnd(), PetscSFGetMultiSF()
1340: @*/
1341: PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *mRootsOrigNumbering[])
1342: {
1343:   PetscSF             msf;
1344:   PetscInt            i, j, k, nroots, nmroots;
1345:   PetscErrorCode      ierr;

1349:   PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);
1352:   PetscSFGetMultiSF(sf,&msf);
1353:   PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL);
1354:   PetscMalloc1(nmroots, mRootsOrigNumbering);
1355:   for (i=0,j=0,k=0; i<nroots; i++) {
1356:     if (!degree[i]) continue;
1357:     for (j=0; j<degree[i]; j++,k++) {
1358:       (*mRootsOrigNumbering)[k] = i;
1359:     }
1360:   }
1361: #if defined(PETSC_USE_DEBUG)
1362:   if (PetscUnlikely(k != nmroots)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail");
1363: #endif
1364:   return(0);
1365: }

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

1370:    Collective

1372:    Input Arguments:
1373: +  sf - star forest
1374: .  unit - data type
1375: .  leafdata - leaf values to use in reduction
1376: -  op - operation to use for reduction

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

1382:    Level: advanced

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

1390: .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
1391: @*/
1392: PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1393: {

1398:   PetscSFSetUp(sf);
1399:   PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);
1400:   (*sf->ops->FetchAndOpBegin)(sf,unit,rootdata,leafdata,leafupdate,op);
1401:   PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);
1402:   return(0);
1403: }

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

1408:    Collective

1410:    Input Arguments:
1411: +  sf - star forest
1412: .  unit - data type
1413: .  leafdata - leaf values to use in reduction
1414: -  op - operation to use for reduction

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

1420:    Level: advanced

1422: .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
1423: @*/
1424: PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1425: {

1430:   PetscSFSetUp(sf);
1431:   PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);
1432:   (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);
1433:   PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);
1434:   return(0);
1435: }

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

1440:    Collective

1442:    Input Arguments:
1443: +  sf - star forest
1444: .  unit - data type
1445: -  leafdata - leaf data to gather to roots

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

1450:    Level: intermediate

1452: .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
1453: @*/
1454: PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
1455: {
1457:   PetscSF        multi;

1461:   PetscSFSetUp(sf);
1462:   PetscSFGetMultiSF(sf,&multi);
1463:   PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);
1464:   return(0);
1465: }

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

1470:    Collective

1472:    Input Arguments:
1473: +  sf - star forest
1474: .  unit - data type
1475: -  leafdata - leaf data to gather to roots

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

1480:    Level: intermediate

1482: .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
1483: @*/
1484: PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
1485: {
1487:   PetscSF        multi;

1491:   PetscSFSetUp(sf);
1492:   PetscSFGetMultiSF(sf,&multi);
1493:   PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);
1494:   return(0);
1495: }

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

1500:    Collective

1502:    Input Arguments:
1503: +  sf - star forest
1504: .  unit - data type
1505: -  multirootdata - root buffer to send to each leaf, one unit of data per leaf

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

1510:    Level: intermediate

1512: .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
1513: @*/
1514: PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
1515: {
1517:   PetscSF        multi;

1521:   PetscSFSetUp(sf);
1522:   PetscSFGetMultiSF(sf,&multi);
1523:   PetscSFBcastBegin(multi,unit,multirootdata,leafdata);
1524:   return(0);
1525: }

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

1530:    Collective

1532:    Input Arguments:
1533: +  sf - star forest
1534: .  unit - data type
1535: -  multirootdata - root buffer to send to each leaf, one unit of data per leaf

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

1540:    Level: intermediate

1542: .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
1543: @*/
1544: PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
1545: {
1547:   PetscSF        multi;

1551:   PetscSFSetUp(sf);
1552:   PetscSFGetMultiSF(sf,&multi);
1553:   PetscSFBcastEnd(multi,unit,multirootdata,leafdata);
1554:   return(0);
1555: }

1557: /*@
1558:   PetscSFCompose - Compose a new PetscSF equivalent to action to PetscSFs

1560:   Input Parameters:
1561: + sfA - The first PetscSF
1562: - sfB - The second PetscSF

1564:   Output Parameters:
1565: . sfBA - equvalent PetscSF for applying A then B

1567:   Level: developer

1569: .seealso: PetscSF, PetscSFGetGraph(), PetscSFSetGraph()
1570: @*/
1571: PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA)
1572: {
1573:   MPI_Comm           comm;
1574:   const PetscSFNode *remotePointsA, *remotePointsB;
1575:   PetscSFNode       *remotePointsBA;
1576:   const PetscInt    *localPointsA, *localPointsB;
1577:   PetscInt           numRootsA, numLeavesA, numRootsB, numLeavesB;
1578:   PetscErrorCode     ierr;

1582:   PetscSFCheckGraphSet(sfA, 1);
1584:   PetscSFCheckGraphSet(sfB, 2);
1586:   PetscObjectGetComm((PetscObject) sfA, &comm);
1587:   PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);
1588:   PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);
1589:   PetscMalloc1(numLeavesB, &remotePointsBA);
1590:   PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, remotePointsBA);
1591:   PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, remotePointsBA);
1592:   PetscSFCreate(comm, sfBA);
1593:   PetscSFSetGraph(*sfBA, numRootsA, numLeavesB, localPointsB, PETSC_COPY_VALUES, remotePointsBA, PETSC_OWN_POINTER);
1594:   return(0);
1595: }