/*@C PetscSFGetGroups - gets incoming and outgoing process groups Collective Input Argument: . sf - star forest Output Arguments: + incoming - group of origin processes for incoming edges (leaves that reference my roots) - outgoing - group of destination processes for outgoing edges (roots that I reference) Level: developer .seealso: PetscSFGetWindow(), PetscSFRestoreWindow() @*/ PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing) { PetscErrorCode ierr; MPI_Group group; PetscFunctionBegin; if (sf->ingroup == MPI_GROUP_NULL) { PetscInt i; const PetscInt *indegree; PetscMPIInt rank,*outranks,*inranks; PetscSFNode *remote; PetscSF bgcount; /* Compute the number of incoming ranks */ ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr); for (i=0; i<sf->nranks; i++) { remote[i].rank = sf->ranks[i]; remote[i].index = 0; } ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr); ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr); ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr); /* Enumerate the incoming ranks */ ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr); ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); for (i=0; i<sf->nranks; i++) outranks[i] = rank; ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr); ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr); ierr = MPI_Group_free(&group);CHKERRQ(ierr); ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr); ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr); } *incoming = sf->ingroup; if (sf->outgroup == MPI_GROUP_NULL) { ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr); ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr); ierr = MPI_Group_free(&group);CHKERRQ(ierr); } *outgoing = sf->outgroup; PetscFunctionReturn(0); }
/*@C PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map Collective Input Arguments: . sf - star forest to invert Output Arguments: . isf - inverse of sf Level: advanced Notes: All roots must have degree 1. The local space may be a permutation, but cannot be sparse. .seealso: PetscSFSetGraph() @*/ PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf) { PetscErrorCode ierr; PetscMPIInt rank; PetscInt i,nroots,nleaves,maxlocal,count,*newilocal; const PetscInt *ilocal; PetscSFNode *roots,*leaves; PetscFunctionBegin; ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr); for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,(ilocal ? ilocal[i] : i)+1); ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr); for (i=0; i<maxlocal; i++) { leaves[i].rank = rank; leaves[i].index = i; } for (i=0; i <nroots; i++) { roots[i].rank = -1; roots[i].index = -1; } ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr); ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr); /* Check whether our leaves are sparse */ for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++; if (count == nroots) newilocal = NULL; else { /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */ ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr); for (i=0,count=0; i<nroots; i++) { if (roots[i].rank >= 0) { newilocal[count] = i; roots[count].rank = roots[i].rank; roots[count].index = roots[i].index; count++; } } } ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr); ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr); ierr = PetscFree2(roots,leaves);CHKERRQ(ierr); PetscFunctionReturn(0); }
/*@C PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices Collective Input Arguments: + sf - original star forest . nroots - number of roots to select on this process - selected - selected roots on this process Output Arguments: . newsf - new star forest Level: advanced Note: To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can be done by calling PetscSFGetGraph(). .seealso: PetscSFSetGraph(), PetscSFGetGraph() @*/ PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nroots,const PetscInt *selected,PetscSF *newsf) { PetscInt *rootdata, *leafdata, *ilocal; PetscSFNode *iremote; PetscInt leafsize = 0, nleaves = 0, n, i; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); if (nroots) PetscValidPointer(selected,3); PetscValidPointer(newsf,4); if (sf->mine) for (i = 0; i < sf->nleaves; ++i) {leafsize = PetscMax(leafsize, sf->mine[i]+1);} else leafsize = sf->nleaves; ierr = PetscCalloc2(sf->nroots,&rootdata,leafsize,&leafdata);CHKERRQ(ierr); for (i=0; i<nroots; ++i) rootdata[selected[i]] = 1; ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); for (i = 0; i < leafsize; ++i) nleaves += leafdata[i]; ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr); ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr); for (i = 0, n = 0; i < sf->nleaves; ++i) { const PetscInt lidx = sf->mine ? sf->mine[i] : i; if (leafdata[lidx]) { ilocal[n] = lidx; iremote[n].rank = sf->remote[i].rank; iremote[n].index = sf->remote[i].index; ++n; } } if (n != nleaves) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "There is a size mismatch in the SF embedding, %d != %d", n, nleaves); ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,newsf);CHKERRQ(ierr); ierr = PetscSFSetGraph(*newsf,sf->nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); PetscFunctionReturn(0); }
/*@C PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices Collective Input Arguments: + sf - original star forest . nleaves - number of leaves to select on this process - selected - selected leaves on this process Output Arguments: . newsf - new star forest Level: advanced .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph() @*/ PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nleaves, const PetscInt *selected, PetscSF *newsf) { PetscSFNode *iremote; PetscInt *ilocal; PetscInt i; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); if (nleaves) PetscValidPointer(selected, 3); PetscValidPointer(newsf, 4); ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr); ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr); for (i = 0; i < nleaves; ++i) { const PetscInt l = selected[i]; ilocal[i] = sf->mine ? sf->mine[l] : l; iremote[i].rank = sf->remote[l].rank; iremote[i].index = sf->remote[l].index; } ierr = PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, newsf);CHKERRQ(ierr); ierr = PetscSFSetGraph(*newsf, sf->nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); PetscFunctionReturn(0); }
PETSC_EXTERN void PETSC_STDCALL petscsfduplicate_(PetscSF sf,PetscSFDuplicateOption *opt,PetscSF *newsf, int *__ierr ){ *__ierr = PetscSFDuplicate( (PetscSF)PetscToPointer((sf) ),*opt,newsf); }
/*@C PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters Collective Input Argument: . sf - star forest that may contain roots with 0 or with more than 1 vertex Output Arguments: . multi - star forest with split roots, such that each root has degree exactly 1 Level: developer Notes: In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming edge, it is a candidate for future optimization that might involve its removal. .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin() @*/ PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); PetscValidPointer(multi,2); if (sf->nroots < 0) { /* Graph has not been set yet; why do we need this? */ ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); *multi = sf->multi; PetscFunctionReturn(0); } if (!sf->multi) { const PetscInt *indegree; PetscInt i,*inoffset,*outones,*outoffset,maxlocal; PetscSFNode *remote; ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr); ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr); for (i=0,maxlocal=0; i<sf->nleaves; i++) maxlocal = PetscMax(maxlocal,(sf->mine ? sf->mine[i] : i)+1); ierr = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr); inoffset[0] = 0; for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i]; for (i=0; i<maxlocal; i++) outones[i] = 1; ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr); ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr); for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */ #if 0 #if defined(PETSC_USE_DEBUG) /* Check that the expected number of increments occurred */ for (i=0; i<sf->nroots; i++) { if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp"); } #endif #endif ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr); for (i=0; i<sf->nleaves; i++) { remote[i].rank = sf->remote[i].rank; remote[i].index = outoffset[sf->mine ? sf->mine[i] : i]; } ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr); ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); if (sf->rankorder) { /* Sort the ranks */ PetscMPIInt rank; PetscInt *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree; PetscSFNode *newremote; ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr); for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]); ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr); for (i=0; i<maxlocal; i++) outranks[i] = rank; ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr); ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr); /* Sort the incoming ranks at each vertex, build the inverse map */ for (i=0; i<sf->nroots; i++) { PetscInt j; for (j=0; j<indegree[i]; j++) tmpoffset[j] = j; ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr); for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j; } ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr); ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr); for (i=0; i<sf->nleaves; i++) { newremote[i].rank = sf->remote[i].rank; newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i]; } ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr); ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr); } ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr); } *multi = sf->multi; 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; }