/*@ PetscSFDestroy - destroy star forest Collective Input Arguments: . sf - address of star forest Level: intermediate .seealso: PetscSFCreate(), PetscSFReset() @*/ PetscErrorCode PetscSFDestroy(PetscSF *sf) { PetscErrorCode ierr; PetscFunctionBegin; if (!*sf) PetscFunctionReturn(0); PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1); if (--((PetscObject)(*sf))->refct > 0) {*sf = 0; PetscFunctionReturn(0);} ierr = PetscSFReset(*sf);CHKERRQ(ierr); if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);} ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr); PetscFunctionReturn(0); }
/*@C PetscSFSetGraph - Set a parallel star forest Collective Input Arguments: + sf - star forest . nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) . nleaves - number of leaf vertices on the current process, each of these references a root on any process . ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage . localmode - copy mode for ilocal . iremote - remote locations of root vertices for each leaf on the current process - remotemode - copy mode for iremote Level: intermediate .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph() @*/ PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode) { PetscErrorCode ierr; PetscTable table; PetscTablePosition pos; PetscMPIInt size; PetscInt i,*rcount,*ranks; PetscFunctionBegin; PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr); if (nleaves && ilocal) PetscValidIntPointer(ilocal,4); if (nleaves) PetscValidPointer(iremote,6); if (nroots < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"roots %D, cannot be negative",nroots); if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves); ierr = PetscSFReset(sf);CHKERRQ(ierr); sf->nroots = nroots; sf->nleaves = nleaves; if (ilocal) { switch (localmode) { case PETSC_COPY_VALUES: ierr = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr); sf->mine = sf->mine_alloc; ierr = PetscMemcpy(sf->mine,ilocal,nleaves*sizeof(*sf->mine));CHKERRQ(ierr); sf->minleaf = PETSC_MAX_INT; sf->maxleaf = PETSC_MIN_INT; for (i=0; i<nleaves; i++) { sf->minleaf = PetscMin(sf->minleaf,ilocal[i]); sf->maxleaf = PetscMax(sf->maxleaf,ilocal[i]); } break; case PETSC_OWN_POINTER: sf->mine_alloc = (PetscInt*)ilocal; sf->mine = sf->mine_alloc; break; case PETSC_USE_POINTER: sf->mine = (PetscInt*)ilocal; break; default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode"); } } if (!ilocal || nleaves > 0) { sf->minleaf = 0; sf->maxleaf = nleaves - 1; } switch (remotemode) { case PETSC_COPY_VALUES: ierr = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr); sf->remote = sf->remote_alloc; ierr = PetscMemcpy(sf->remote,iremote,nleaves*sizeof(*sf->remote));CHKERRQ(ierr); break; case PETSC_OWN_POINTER: sf->remote_alloc = (PetscSFNode*)iremote; sf->remote = sf->remote_alloc; break; case PETSC_USE_POINTER: sf->remote = (PetscSFNode*)iremote; break; default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode"); } ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr); ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr); for (i=0; i<nleaves; i++) { /* Log 1-based rank */ ierr = PetscTableAdd(table,iremote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr); } ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr); ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,nleaves,&sf->rmine,nleaves,&sf->rremote);CHKERRQ(ierr); ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr); ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr); for (i=0; i<sf->nranks; i++) { ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr); ranks[i]--; /* Convert back to 0-based */ } ierr = PetscTableDestroy(&table);CHKERRQ(ierr); ierr = PetscSortIntWithArray(sf->nranks,ranks,rcount);CHKERRQ(ierr); sf->roffset[0] = 0; for (i=0; i<sf->nranks; i++) { ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr); sf->roffset[i+1] = sf->roffset[i] + rcount[i]; rcount[i] = 0; } for (i=0; i<nleaves; i++) { PetscInt lo,hi,irank; /* Search for index of iremote[i].rank in sf->ranks */ lo = 0; hi = sf->nranks; while (hi - lo > 1) { PetscInt mid = lo + (hi - lo)/2; if (iremote[i].rank < sf->ranks[mid]) hi = mid; else lo = mid; } if (hi - lo == 1 && iremote[i].rank == sf->ranks[lo]) irank = lo; else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",iremote[i].rank); sf->rmine[sf->roffset[irank] + rcount[irank]] = ilocal ? ilocal[i] : i; sf->rremote[sf->roffset[irank] + rcount[irank]] = iremote[i].index; rcount[irank]++; } ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr); #if !defined(PETSC_USE_64BIT_INDICES) if (nroots == PETSC_DETERMINE) { /* Jed, if you have a better way to do this, put it in */ PetscInt *numRankLeaves, *leafOff, *leafIndices, *numRankRoots, *rootOff, *rootIndices, maxRoots = 0; /* All to all to determine number of leaf indices from each (you can do this using Scan and asynch messages) */ ierr = PetscMalloc4(size,&numRankLeaves,size+1,&leafOff,size,&numRankRoots,size+1,&rootOff);CHKERRQ(ierr); ierr = PetscMemzero(numRankLeaves, size * sizeof(PetscInt));CHKERRQ(ierr); for (i = 0; i < nleaves; ++i) ++numRankLeaves[iremote[i].rank]; ierr = MPI_Alltoall(numRankLeaves, 1, MPIU_INT, numRankRoots, 1, MPIU_INT, PetscObjectComm((PetscObject)sf));CHKERRQ(ierr); /* Could set nroots to this maximum */ for (i = 0; i < size; ++i) maxRoots += numRankRoots[i]; /* Gather all indices */ ierr = PetscMalloc2(nleaves,&leafIndices,maxRoots,&rootIndices);CHKERRQ(ierr); leafOff[0] = 0; for (i = 0; i < size; ++i) leafOff[i+1] = leafOff[i] + numRankLeaves[i]; for (i = 0; i < nleaves; ++i) leafIndices[leafOff[iremote[i].rank]++] = iremote[i].index; leafOff[0] = 0; for (i = 0; i < size; ++i) leafOff[i+1] = leafOff[i] + numRankLeaves[i]; rootOff[0] = 0; for (i = 0; i < size; ++i) rootOff[i+1] = rootOff[i] + numRankRoots[i]; ierr = MPI_Alltoallv(leafIndices, numRankLeaves, leafOff, MPIU_INT, rootIndices, numRankRoots, rootOff, MPIU_INT, PetscObjectComm((PetscObject)sf));CHKERRQ(ierr); /* Sort and reduce */ ierr = PetscSortRemoveDupsInt(&maxRoots, rootIndices);CHKERRQ(ierr); ierr = PetscFree2(leafIndices,rootIndices);CHKERRQ(ierr); ierr = PetscFree4(numRankLeaves,leafOff,numRankRoots,rootOff);CHKERRQ(ierr); sf->nroots = maxRoots; } #endif sf->graphset = PETSC_TRUE; ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); }
int main(int argc,char **argv) { PetscSF sf,sfDup,sfInv,sfEmbed,sfA,sfB,sfBA; const PetscInt *degree; PetscErrorCode ierr; ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = CheckGraphNotSet(sf);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = CheckGraphNotSet(sf);CHKERRQ(ierr); ierr = PetscSFReset(sf);CHKERRQ(ierr); ierr = CheckGraphNotSet(sf);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = CheckGraphNotSet(sf);CHKERRQ(ierr); ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); ierr = CheckGraphNotSet(sf);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = CheckGraphNotSet(sf);CHKERRQ(ierr); ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); ierr = CheckGraphNotSet(sf);CHKERRQ(ierr); ierr = PetscSFReset(sf);CHKERRQ(ierr); ierr = CheckGraphNotSet(sf);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_COPY_VALUES,NULL,PETSC_COPY_VALUES);CHKERRQ(ierr); ierr = CheckGraphEmpty(sf);CHKERRQ(ierr); ierr = PetscSFReset(sf);CHKERRQ(ierr); ierr = CheckGraphNotSet(sf);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_COPY_VALUES,NULL,PETSC_COPY_VALUES);CHKERRQ(ierr); ierr = CheckGraphEmpty(sf);CHKERRQ(ierr); ierr = PetscSFReset(sf);CHKERRQ(ierr); ierr = CheckGraphNotSet(sf);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test setup */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = CheckRanksNotSet(sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_COPY_VALUES,NULL,PETSC_COPY_VALUES);CHKERRQ(ierr); ierr = CheckRanksNotSet(sf);CHKERRQ(ierr); ierr = PetscSFSetUp(sf);CHKERRQ(ierr); ierr = CheckRanksEmpty(sf);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test setup then reset */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_COPY_VALUES,NULL,PETSC_COPY_VALUES);CHKERRQ(ierr); ierr = PetscSFSetUp(sf);CHKERRQ(ierr); ierr = PetscSFReset(sf);CHKERRQ(ierr); ierr = CheckRanksNotSet(sf);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test view (no graph set, no type set) */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFView(sf,NULL);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test set graph then view (no type set) */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_COPY_VALUES,NULL,PETSC_COPY_VALUES);CHKERRQ(ierr); ierr = PetscSFView(sf,NULL);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test set type then view (no graph set) */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); ierr = PetscSFView(sf,NULL);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test set type then graph then view */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_COPY_VALUES,NULL,PETSC_COPY_VALUES);CHKERRQ(ierr); ierr = PetscSFView(sf,NULL);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test set graph then type */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_COPY_VALUES,NULL,PETSC_COPY_VALUES);CHKERRQ(ierr); ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); ierr = CheckGraphEmpty(sf);CHKERRQ(ierr); ierr = PetscSFReset(sf);CHKERRQ(ierr); ierr = CheckGraphNotSet(sf);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test Bcast */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_COPY_VALUES,NULL,PETSC_COPY_VALUES);CHKERRQ(ierr); ierr = PetscSFBcastBegin(sf,MPI_INT,NULL,NULL);CHKERRQ(ierr); ierr = PetscSFBcastEnd (sf,MPI_INT,NULL,NULL);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test Reduce */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_COPY_VALUES,NULL,PETSC_COPY_VALUES);CHKERRQ(ierr); ierr = PetscSFReduceBegin(sf,MPI_INT,NULL,NULL,MPIU_REPLACE);CHKERRQ(ierr); ierr = PetscSFReduceEnd (sf,MPI_INT,NULL,NULL,MPIU_REPLACE);CHKERRQ(ierr); ierr = PetscSFReduceBegin(sf,MPI_INT,NULL,NULL,MPI_SUM);CHKERRQ(ierr); ierr = PetscSFReduceEnd (sf,MPI_INT,NULL,NULL,MPI_SUM);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test FetchAndOp */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_COPY_VALUES,NULL,PETSC_COPY_VALUES);CHKERRQ(ierr); ierr = PetscSFFetchAndOpBegin(sf,MPI_INT,NULL,NULL,NULL,MPI_SUM);CHKERRQ(ierr); ierr = PetscSFFetchAndOpEnd (sf,MPI_INT,NULL,NULL,NULL,MPI_SUM);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test ComputeDegree */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_COPY_VALUES,NULL,PETSC_COPY_VALUES);CHKERRQ(ierr); ierr = PetscSFComputeDegreeBegin(sf,°ree);CHKERRQ(ierr); ierr = PetscSFComputeDegreeEnd(sf,°ree);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test PetscSFDuplicate() */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_USE_POINTER,NULL,PETSC_USE_POINTER);CHKERRQ(ierr); ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_GRAPH,&sfDup);CHKERRQ(ierr); ierr = CheckGraphEmpty(sfDup);CHKERRQ(ierr); ierr = PetscSFDestroy(&sfDup);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test PetscSFCreateInverseSF() */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_USE_POINTER,NULL,PETSC_USE_POINTER);CHKERRQ(ierr); ierr = PetscSFCreateInverseSF(sf,&sfInv);CHKERRQ(ierr); ierr = CheckGraphEmpty(sfInv);CHKERRQ(ierr); ierr = PetscSFDestroy(&sfInv);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test PetscSFCreateEmbeddedSF() */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_USE_POINTER,NULL,PETSC_USE_POINTER);CHKERRQ(ierr); ierr = PetscSFCreateEmbeddedSF(sf,0,NULL,&sfEmbed);CHKERRQ(ierr); ierr = CheckGraphEmpty(sfEmbed);CHKERRQ(ierr); ierr = PetscSFDestroy(&sfEmbed);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test PetscSFCreateEmbeddedLeafSF() */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf,0,0,NULL,PETSC_USE_POINTER,NULL,PETSC_USE_POINTER);CHKERRQ(ierr); ierr = PetscSFCreateEmbeddedLeafSF(sf,0,NULL,&sfEmbed);CHKERRQ(ierr); ierr = CheckGraphEmpty(sfEmbed);CHKERRQ(ierr); ierr = PetscSFDestroy(&sfEmbed);CHKERRQ(ierr); ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); /* Test PetscSFCompose() */ ierr = PetscSFCreate(PETSC_COMM_WORLD,&sfA);CHKERRQ(ierr); ierr = PetscSFSetGraph(sfA,0,0,NULL,PETSC_USE_POINTER,NULL,PETSC_USE_POINTER);CHKERRQ(ierr); ierr = PetscSFCreate(PETSC_COMM_WORLD,&sfB);CHKERRQ(ierr); ierr = PetscSFSetGraph(sfB,0,0,NULL,PETSC_USE_POINTER,NULL,PETSC_USE_POINTER);CHKERRQ(ierr); ierr = PetscSFCompose(sfA,sfB,&sfBA);CHKERRQ(ierr); ierr = CheckGraphEmpty(sfBA);CHKERRQ(ierr); ierr = PetscSFDestroy(&sfBA);CHKERRQ(ierr); ierr = PetscSFDestroy(&sfA);CHKERRQ(ierr); ierr = PetscSFDestroy(&sfB);CHKERRQ(ierr); ierr = PetscFinalize(); return ierr; }