xref: /petsc/src/vec/is/sf/interface/sf.c (revision 59af0bd3658d6c64d35e37f76ad6a8a026fa611f)
1af0996ceSBarry Smith #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/
2c4e6a40aSLawrence Mitchell #include <petsc/private/hashseti.h>
395fce210SBarry Smith #include <petscctable.h>
495fce210SBarry Smith 
57fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
67fd2d3dbSJunchao Zhang   #include <cuda_runtime.h>
77fd2d3dbSJunchao Zhang #endif
87fd2d3dbSJunchao Zhang 
97fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_HIP)
107fd2d3dbSJunchao Zhang   #include <hip/hip_runtime.h>
117fd2d3dbSJunchao Zhang #endif
127fd2d3dbSJunchao Zhang 
1395fce210SBarry Smith #if defined(PETSC_USE_DEBUG)
1495fce210SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) do {                          \
1595fce210SBarry Smith     if (PetscUnlikely(!(sf)->graphset))                              \
16dd5b3ca6SJunchao Zhang       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() or PetscSFSetGraphWithPattern() on argument %D \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME); \
1795fce210SBarry Smith   } while (0)
1895fce210SBarry Smith #else
1995fce210SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) do {} while (0)
2095fce210SBarry Smith #endif
2195fce210SBarry Smith 
224c8fdceaSLisandro Dalcin const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",NULL};
2395fce210SBarry Smith 
247fd2d3dbSJunchao Zhang PETSC_STATIC_INLINE PetscErrorCode PetscGetMemType(const void *data,PetscMemType *type)
257fd2d3dbSJunchao Zhang {
267fd2d3dbSJunchao Zhang   PetscFunctionBegin;
277fd2d3dbSJunchao Zhang   PetscValidPointer(type,2);
287fd2d3dbSJunchao Zhang   *type = PETSC_MEMTYPE_HOST;
297fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
307fd2d3dbSJunchao Zhang   if (PetscCUDAInitialized && data) {
317fd2d3dbSJunchao Zhang     cudaError_t                  cerr;
327fd2d3dbSJunchao Zhang     struct cudaPointerAttributes attr;
337fd2d3dbSJunchao Zhang     enum cudaMemoryType          mtype;
347fd2d3dbSJunchao Zhang     cerr = cudaPointerGetAttributes(&attr,data); /* Do not check error since before CUDA 11.0, passing a host pointer returns cudaErrorInvalidValue */
357fd2d3dbSJunchao Zhang     cudaGetLastError(); /* Reset the last error */
367fd2d3dbSJunchao Zhang     #if (CUDART_VERSION < 10000)
377fd2d3dbSJunchao Zhang       mtype = attr.memoryType;
387fd2d3dbSJunchao Zhang     #else
397fd2d3dbSJunchao Zhang       mtype = attr.type;
407fd2d3dbSJunchao Zhang     #endif
417fd2d3dbSJunchao Zhang     if (cerr == cudaSuccess && mtype == cudaMemoryTypeDevice) *type = PETSC_MEMTYPE_DEVICE;
427fd2d3dbSJunchao Zhang   }
437fd2d3dbSJunchao Zhang #endif
447fd2d3dbSJunchao Zhang 
457fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_HIP)
467fd2d3dbSJunchao Zhang   if (PetscHIPInitialized && data) {
477fd2d3dbSJunchao Zhang     hipError_t                   cerr;
487fd2d3dbSJunchao Zhang     struct hipPointerAttribute_t attr;
497fd2d3dbSJunchao Zhang     enum hipMemoryType           mtype;
50*59af0bd3SScott Kruger     cerr = hipPointerGetAttributes(&attr,data);
517fd2d3dbSJunchao Zhang     hipGetLastError(); /* Reset the last error */
527fd2d3dbSJunchao Zhang     mtype = attr.memoryType;
537fd2d3dbSJunchao Zhang     if (cerr == hipSuccess && mtype == hipMemoryTypeDevice) *type = PETSC_MEMTYPE_DEVICE;
547fd2d3dbSJunchao Zhang   }
557fd2d3dbSJunchao Zhang #endif
567fd2d3dbSJunchao Zhang   PetscFunctionReturn(0);
577fd2d3dbSJunchao Zhang }
587fd2d3dbSJunchao Zhang 
598af6ec1cSBarry Smith /*@
6095fce210SBarry Smith    PetscSFCreate - create a star forest communication context
6195fce210SBarry Smith 
62d083f849SBarry Smith    Collective
6395fce210SBarry Smith 
6495fce210SBarry Smith    Input Arguments:
6595fce210SBarry Smith .  comm - communicator on which the star forest will operate
6695fce210SBarry Smith 
6795fce210SBarry Smith    Output Arguments:
6895fce210SBarry Smith .  sf - new star forest context
6995fce210SBarry Smith 
70dd5b3ca6SJunchao Zhang    Options Database Keys:
71dd5b3ca6SJunchao Zhang +  -sf_type basic     -Use MPI persistent Isend/Irecv for communication (Default)
72dd5b3ca6SJunchao Zhang .  -sf_type window    -Use MPI-3 one-sided window for communication
73dd5b3ca6SJunchao Zhang -  -sf_type neighbor  -Use MPI-3 neighborhood collectives for communication
74dd5b3ca6SJunchao Zhang 
7595fce210SBarry Smith    Level: intermediate
7695fce210SBarry Smith 
77dd5b3ca6SJunchao Zhang    Notes:
78dd5b3ca6SJunchao Zhang    When one knows the communication graph is one of the predefined graph, such as MPI_Alltoall, MPI_Allgatherv,
79dd5b3ca6SJunchao Zhang    MPI_Gatherv, one can create a PetscSF and then set its graph with PetscSFSetGraphWithPattern(). These special
80dd5b3ca6SJunchao Zhang    SFs are optimized and they have better performance than general SFs.
81dd5b3ca6SJunchao Zhang 
82dd5b3ca6SJunchao Zhang .seealso: PetscSFSetGraph(), PetscSFSetGraphWithPattern(), PetscSFDestroy()
8395fce210SBarry Smith @*/
8495fce210SBarry Smith PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf)
8595fce210SBarry Smith {
8695fce210SBarry Smith   PetscErrorCode ierr;
8795fce210SBarry Smith   PetscSF        b;
8895fce210SBarry Smith 
8995fce210SBarry Smith   PetscFunctionBegin;
9095fce210SBarry Smith   PetscValidPointer(sf,2);
91607a6623SBarry Smith   ierr = PetscSFInitializePackage();CHKERRQ(ierr);
9295fce210SBarry Smith 
9373107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(b,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);CHKERRQ(ierr);
9495fce210SBarry Smith 
9595fce210SBarry Smith   b->nroots    = -1;
9695fce210SBarry Smith   b->nleaves   = -1;
9729046d53SLisandro Dalcin   b->minleaf   = PETSC_MAX_INT;
9829046d53SLisandro Dalcin   b->maxleaf   = PETSC_MIN_INT;
9995fce210SBarry Smith   b->nranks    = -1;
10095fce210SBarry Smith   b->rankorder = PETSC_TRUE;
10195fce210SBarry Smith   b->ingroup   = MPI_GROUP_NULL;
10295fce210SBarry Smith   b->outgroup  = MPI_GROUP_NULL;
10395fce210SBarry Smith   b->graphset  = PETSC_FALSE;
10420c24465SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
10520c24465SJunchao Zhang   b->use_gpu_aware_mpi    = use_gpu_aware_mpi;
10620c24465SJunchao Zhang   b->use_stream_aware_mpi = PETSC_FALSE;
10720c24465SJunchao Zhang   b->use_default_stream   = PETSC_TRUE; /* The assumption is true for PETSc internal use of SF */
108*59af0bd3SScott Kruger   /* Set the default */
10927f636e8SJunchao Zhang   #if defined(PETSC_HAVE_KOKKOS) /* Prefer kokkos over cuda*/
11020c24465SJunchao Zhang     b->backend = PETSCSF_BACKEND_KOKKOS;
11127f636e8SJunchao Zhang   #elif defined(PETSC_HAVE_CUDA)
11227f636e8SJunchao Zhang     b->backend = PETSCSF_BACKEND_CUDA;
113*59af0bd3SScott Kruger   #elif defined(PETSC_HAVE_HIP)
114*59af0bd3SScott Kruger     b->backend = PETSCSF_BACKEND_HIP;
11520c24465SJunchao Zhang   #endif
11620c24465SJunchao Zhang #endif
11795fce210SBarry Smith   *sf = b;
11895fce210SBarry Smith   PetscFunctionReturn(0);
11995fce210SBarry Smith }
12095fce210SBarry Smith 
12129046d53SLisandro Dalcin /*@
12295fce210SBarry Smith    PetscSFReset - Reset a star forest so that different sizes or neighbors can be used
12395fce210SBarry Smith 
12495fce210SBarry Smith    Collective
12595fce210SBarry Smith 
12695fce210SBarry Smith    Input Arguments:
12795fce210SBarry Smith .  sf - star forest
12895fce210SBarry Smith 
12995fce210SBarry Smith    Level: advanced
13095fce210SBarry Smith 
13195fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy()
13295fce210SBarry Smith @*/
13395fce210SBarry Smith PetscErrorCode PetscSFReset(PetscSF sf)
13495fce210SBarry Smith {
13595fce210SBarry Smith   PetscErrorCode ierr;
13695fce210SBarry Smith 
13795fce210SBarry Smith   PetscFunctionBegin;
13895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
13979715d56SJed Brown   if (sf->ops->Reset) {ierr = (*sf->ops->Reset)(sf);CHKERRQ(ierr);}
14029046d53SLisandro Dalcin   sf->nroots   = -1;
14129046d53SLisandro Dalcin   sf->nleaves  = -1;
14229046d53SLisandro Dalcin   sf->minleaf  = PETSC_MAX_INT;
14329046d53SLisandro Dalcin   sf->maxleaf  = PETSC_MIN_INT;
14495fce210SBarry Smith   sf->mine     = NULL;
14595fce210SBarry Smith   sf->remote   = NULL;
14629046d53SLisandro Dalcin   sf->graphset = PETSC_FALSE;
14729046d53SLisandro Dalcin   ierr = PetscFree(sf->mine_alloc);CHKERRQ(ierr);
14895fce210SBarry Smith   ierr = PetscFree(sf->remote_alloc);CHKERRQ(ierr);
14921c688dcSJed Brown   sf->nranks = -1;
15029046d53SLisandro Dalcin   ierr = PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);CHKERRQ(ierr);
1517fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
15220c24465SJunchao Zhang   for (PetscInt i=0; i<2; i++) {ierr = PetscSFFree(sf,PETSC_MEMTYPE_DEVICE,sf->rmine_d[i]);CHKERRQ(ierr);}
153eb02082bSJunchao Zhang #endif
15429046d53SLisandro Dalcin   sf->degreeknown = PETSC_FALSE;
15595fce210SBarry Smith   ierr = PetscFree(sf->degree);CHKERRQ(ierr);
15695fce210SBarry Smith   if (sf->ingroup  != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRQ(ierr);}
15795fce210SBarry Smith   if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRQ(ierr);}
15895fce210SBarry Smith   ierr = PetscSFDestroy(&sf->multi);CHKERRQ(ierr);
159dd5b3ca6SJunchao Zhang   ierr = PetscLayoutDestroy(&sf->map);CHKERRQ(ierr);
16095fce210SBarry Smith   sf->setupcalled = PETSC_FALSE;
16195fce210SBarry Smith   PetscFunctionReturn(0);
16295fce210SBarry Smith }
16395fce210SBarry Smith 
16495fce210SBarry Smith /*@C
16529046d53SLisandro Dalcin    PetscSFSetType - Set the PetscSF communication implementation
16695fce210SBarry Smith 
16795fce210SBarry Smith    Collective on PetscSF
16895fce210SBarry Smith 
16995fce210SBarry Smith    Input Parameters:
17095fce210SBarry Smith +  sf - the PetscSF context
17195fce210SBarry Smith -  type - a known method
17295fce210SBarry Smith 
17395fce210SBarry Smith    Options Database Key:
17495fce210SBarry Smith .  -sf_type <type> - Sets the method; use -help for a list
17570616304SStefano Zampini    of available methods (for instance, window, basic, neighbor)
17695fce210SBarry Smith 
17795fce210SBarry Smith    Notes:
17895fce210SBarry Smith    See "include/petscsf.h" for available methods (for instance)
17995fce210SBarry Smith +    PETSCSFWINDOW - MPI-2/3 one-sided
18095fce210SBarry Smith -    PETSCSFBASIC - basic implementation using MPI-1 two-sided
18195fce210SBarry Smith 
18295fce210SBarry Smith   Level: intermediate
18395fce210SBarry Smith 
18495fce210SBarry Smith .seealso: PetscSFType, PetscSFCreate()
18595fce210SBarry Smith @*/
18695fce210SBarry Smith PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type)
18795fce210SBarry Smith {
18895fce210SBarry Smith   PetscErrorCode ierr,(*r)(PetscSF);
18995fce210SBarry Smith   PetscBool      match;
19095fce210SBarry Smith 
19195fce210SBarry Smith   PetscFunctionBegin;
19295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
19395fce210SBarry Smith   PetscValidCharPointer(type,2);
19495fce210SBarry Smith 
19595fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)sf,type,&match);CHKERRQ(ierr);
19695fce210SBarry Smith   if (match) PetscFunctionReturn(0);
19795fce210SBarry Smith 
198adc40e5bSBarry Smith   ierr = PetscFunctionListFind(PetscSFList,type,&r);CHKERRQ(ierr);
19995fce210SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type);
20029046d53SLisandro Dalcin   /* Destroy the previous PetscSF implementation context */
20129046d53SLisandro Dalcin   if (sf->ops->Destroy) {ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr);}
20295fce210SBarry Smith   ierr = PetscMemzero(sf->ops,sizeof(*sf->ops));CHKERRQ(ierr);
20395fce210SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)sf,type);CHKERRQ(ierr);
20495fce210SBarry Smith   ierr = (*r)(sf);CHKERRQ(ierr);
20595fce210SBarry Smith   PetscFunctionReturn(0);
20695fce210SBarry Smith }
20795fce210SBarry Smith 
20829046d53SLisandro Dalcin /*@C
20929046d53SLisandro Dalcin   PetscSFGetType - Get the PetscSF communication implementation
21029046d53SLisandro Dalcin 
21129046d53SLisandro Dalcin   Not Collective
21229046d53SLisandro Dalcin 
21329046d53SLisandro Dalcin   Input Parameter:
21429046d53SLisandro Dalcin . sf  - the PetscSF context
21529046d53SLisandro Dalcin 
21629046d53SLisandro Dalcin   Output Parameter:
21729046d53SLisandro Dalcin . type - the PetscSF type name
21829046d53SLisandro Dalcin 
21929046d53SLisandro Dalcin   Level: intermediate
22029046d53SLisandro Dalcin 
22129046d53SLisandro Dalcin .seealso: PetscSFSetType(), PetscSFCreate()
22229046d53SLisandro Dalcin @*/
22329046d53SLisandro Dalcin PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type)
22429046d53SLisandro Dalcin {
22529046d53SLisandro Dalcin   PetscFunctionBegin;
22629046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID,1);
22729046d53SLisandro Dalcin   PetscValidPointer(type,2);
22829046d53SLisandro Dalcin   *type = ((PetscObject)sf)->type_name;
22929046d53SLisandro Dalcin   PetscFunctionReturn(0);
23029046d53SLisandro Dalcin }
23129046d53SLisandro Dalcin 
232d36d33e4SMatthew G. Knepley /*@
23395fce210SBarry Smith    PetscSFDestroy - destroy star forest
23495fce210SBarry Smith 
23595fce210SBarry Smith    Collective
23695fce210SBarry Smith 
23795fce210SBarry Smith    Input Arguments:
23895fce210SBarry Smith .  sf - address of star forest
23995fce210SBarry Smith 
24095fce210SBarry Smith    Level: intermediate
24195fce210SBarry Smith 
24295fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFReset()
24395fce210SBarry Smith @*/
24495fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf)
24595fce210SBarry Smith {
24695fce210SBarry Smith   PetscErrorCode ierr;
24795fce210SBarry Smith 
24895fce210SBarry Smith   PetscFunctionBegin;
24995fce210SBarry Smith   if (!*sf) PetscFunctionReturn(0);
25095fce210SBarry Smith   PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1);
25129046d53SLisandro Dalcin   if (--((PetscObject)(*sf))->refct > 0) {*sf = NULL; PetscFunctionReturn(0);}
25295fce210SBarry Smith   ierr = PetscSFReset(*sf);CHKERRQ(ierr);
25395fce210SBarry Smith   if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);}
25497929ea7SJunchao Zhang   ierr = PetscSFDestroy(&(*sf)->vscat.lsf);CHKERRQ(ierr);
25597929ea7SJunchao Zhang   if ((*sf)->vscat.bs > 1) {ierr = MPI_Type_free(&(*sf)->vscat.unit);CHKERRQ(ierr);}
25695fce210SBarry Smith   ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr);
25795fce210SBarry Smith   PetscFunctionReturn(0);
25895fce210SBarry Smith }
25995fce210SBarry Smith 
260c4e6a40aSLawrence Mitchell static PetscErrorCode PetscSFCheckGraphValid_Private(PetscSF sf)
261c4e6a40aSLawrence Mitchell {
262c4e6a40aSLawrence Mitchell   PetscInt           i, nleaves;
263c4e6a40aSLawrence Mitchell   PetscMPIInt        size;
264c4e6a40aSLawrence Mitchell   const PetscInt    *ilocal;
265c4e6a40aSLawrence Mitchell   const PetscSFNode *iremote;
266c4e6a40aSLawrence Mitchell   PetscErrorCode     ierr;
267c4e6a40aSLawrence Mitchell 
268c4e6a40aSLawrence Mitchell   PetscFunctionBegin;
26976bd3646SJed Brown   if (!sf->graphset || !PetscDefined(USE_DEBUG)) PetscFunctionReturn(0);
270c4e6a40aSLawrence Mitchell   ierr = PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
271c4e6a40aSLawrence Mitchell   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
272c4e6a40aSLawrence Mitchell   for (i = 0; i < nleaves; i++) {
273c4e6a40aSLawrence Mitchell     const PetscInt rank = iremote[i].rank;
274c4e6a40aSLawrence Mitchell     const PetscInt remote = iremote[i].index;
275c4e6a40aSLawrence Mitchell     const PetscInt leaf = ilocal ? ilocal[i] : i;
276c4e6a40aSLawrence Mitchell     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);
277c4e6a40aSLawrence Mitchell     if (remote < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided index (%D) for remote %D is invalid, should be >= 0",remote,i);
278c4e6a40aSLawrence Mitchell     if (leaf < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided location (%D) for leaf %D is invalid, should be >= 0",leaf,i);
279c4e6a40aSLawrence Mitchell   }
280c4e6a40aSLawrence Mitchell   PetscFunctionReturn(0);
281c4e6a40aSLawrence Mitchell }
282c4e6a40aSLawrence Mitchell 
28320c24465SJunchao Zhang 
28420c24465SJunchao Zhang 
28520c24465SJunchao Zhang 
28695fce210SBarry Smith /*@
28795fce210SBarry Smith    PetscSFSetUp - set up communication structures
28895fce210SBarry Smith 
28995fce210SBarry Smith    Collective
29095fce210SBarry Smith 
29195fce210SBarry Smith    Input Arguments:
29295fce210SBarry Smith .  sf - star forest communication object
29395fce210SBarry Smith 
29495fce210SBarry Smith    Level: beginner
29595fce210SBarry Smith 
29695fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFSetType()
29795fce210SBarry Smith @*/
29895fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf)
29995fce210SBarry Smith {
30095fce210SBarry Smith   PetscErrorCode ierr;
30195fce210SBarry Smith 
30295fce210SBarry Smith   PetscFunctionBegin;
30329046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
30429046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
30595fce210SBarry Smith   if (sf->setupcalled) PetscFunctionReturn(0);
30629046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr);
30720c24465SJunchao Zhang   ierr = PetscSFCheckGraphValid_Private(sf);CHKERRQ(ierr);
30820c24465SJunchao Zhang   if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);} /* Zero all sf->ops */
30995fce210SBarry Smith   if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);}
31020c24465SJunchao Zhang #if defined(PETSC_HAVE_CUDA)
31120c24465SJunchao Zhang   if (sf->backend == PETSCSF_BACKEND_CUDA) {
31220c24465SJunchao Zhang     sf->ops->Malloc = PetscSFMalloc_Cuda;
31320c24465SJunchao Zhang     sf->ops->Free   = PetscSFFree_Cuda;
31420c24465SJunchao Zhang   }
31520c24465SJunchao Zhang #endif
316*59af0bd3SScott Kruger #if defined(PETSC_HAVE_HIP)
317*59af0bd3SScott Kruger   if (sf->backend == PETSCSF_BACKEND_HIP) {
318*59af0bd3SScott Kruger     sf->ops->Malloc = PetscSFMalloc_HIP;
319*59af0bd3SScott Kruger     sf->ops->Free   = PetscSFFree_HIP;
320*59af0bd3SScott Kruger   }
321*59af0bd3SScott Kruger #endif
32220c24465SJunchao Zhang 
323*59af0bd3SScott Kruger #
32420c24465SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
32520c24465SJunchao Zhang   if (sf->backend == PETSCSF_BACKEND_KOKKOS) {
32620c24465SJunchao Zhang     sf->ops->Malloc = PetscSFMalloc_Kokkos;
32720c24465SJunchao Zhang     sf->ops->Free   = PetscSFFree_Kokkos;
32820c24465SJunchao Zhang   }
32920c24465SJunchao Zhang #endif
33029046d53SLisandro Dalcin   ierr = PetscLogEventEnd(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr);
33195fce210SBarry Smith   sf->setupcalled = PETSC_TRUE;
33295fce210SBarry Smith   PetscFunctionReturn(0);
33395fce210SBarry Smith }
33495fce210SBarry Smith 
3358af6ec1cSBarry Smith /*@
33695fce210SBarry Smith    PetscSFSetFromOptions - set PetscSF options using the options database
33795fce210SBarry Smith 
33895fce210SBarry Smith    Logically Collective
33995fce210SBarry Smith 
34095fce210SBarry Smith    Input Arguments:
34195fce210SBarry Smith .  sf - star forest
34295fce210SBarry Smith 
34395fce210SBarry Smith    Options Database Keys:
34460263706SJed Brown +  -sf_type               - implementation type, see PetscSFSetType()
34551ccb202SJunchao Zhang .  -sf_rank_order         - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
346b85e67b7SJunchao Zhang .  -sf_use_default_stream - Assume callers of SF computed the input root/leafdata with the default cuda stream. SF will also
347c2a741eeSJunchao Zhang                             use the default stream to process data. Therefore, no stream synchronization is needed between SF and its caller (default: true).
348c06a8e02SRichard Tran Mills                             If true, this option only works with -use_gpu_aware_mpi 1.
34920c24465SJunchao Zhang .  -sf_use_stream_aware_mpi  - Assume the underlying MPI is cuda-stream aware and SF won't sync streams for send/recv buffers passed to MPI (default: false).
350c06a8e02SRichard Tran Mills                                If true, this option only works with -use_gpu_aware_mpi 1.
35195fce210SBarry Smith 
352*59af0bd3SScott Kruger -  -sf_backend cuda | hip | kokkos -Select the device backend SF uses. Currently SF has these backends: cuda, hip and Kokkos.
353*59af0bd3SScott Kruger                               On CUDA (HIP) devices, one can choose cuda (hip) or kokkos with the default being kokkos. On other devices,
35420c24465SJunchao Zhang                               the only available is kokkos.
35520c24465SJunchao Zhang 
35695fce210SBarry Smith    Level: intermediate
35795fce210SBarry Smith @*/
35895fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
35995fce210SBarry Smith {
36095fce210SBarry Smith   PetscSFType    deft;
36195fce210SBarry Smith   char           type[256];
36295fce210SBarry Smith   PetscErrorCode ierr;
36395fce210SBarry Smith   PetscBool      flg;
36495fce210SBarry Smith 
36595fce210SBarry Smith   PetscFunctionBegin;
36695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
36795fce210SBarry Smith   ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr);
36895fce210SBarry Smith   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
36929046d53SLisandro Dalcin   ierr = PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,sizeof(type),&flg);CHKERRQ(ierr);
37095fce210SBarry Smith   ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr);
37195fce210SBarry Smith   ierr = 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);CHKERRQ(ierr);
3727fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
37320c24465SJunchao Zhang   {
37420c24465SJunchao Zhang     char        backendstr[32] = {0};
375*59af0bd3SScott Kruger     PetscBool   isCuda = PETSC_FALSE,isHip = PETSC_FALSE,isKokkos = PETSC_FALSE,set;
37620c24465SJunchao Zhang     /* Change the defaults set in PetscSFCreate() with command line options */
377cd620004SJunchao Zhang     ierr = PetscOptionsBool("-sf_use_default_stream","Assume SF's input and output root/leafdata is computed on the default stream","PetscSFSetFromOptions",sf->use_default_stream,&sf->use_default_stream,NULL);CHKERRQ(ierr);
378b85e67b7SJunchao Zhang     ierr = PetscOptionsBool("-sf_use_stream_aware_mpi","Assume the underlying MPI is cuda-stream aware","PetscSFSetFromOptions",sf->use_stream_aware_mpi,&sf->use_stream_aware_mpi,NULL);CHKERRQ(ierr);
37920c24465SJunchao Zhang     ierr = PetscOptionsString("-sf_backend","Select the device backend SF uses","PetscSFSetFromOptions",NULL,backendstr,sizeof(backendstr),&set);CHKERRQ(ierr);
38020c24465SJunchao Zhang     ierr = PetscStrcasecmp("cuda",backendstr,&isCuda);CHKERRQ(ierr);
38120c24465SJunchao Zhang     ierr = PetscStrcasecmp("kokkos",backendstr,&isKokkos);CHKERRQ(ierr);
382*59af0bd3SScott Kruger     ierr = PetscStrcasecmp("hip",backendstr,&isHip);CHKERRQ(ierr);
383*59af0bd3SScott Kruger   #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
38420c24465SJunchao Zhang     if (isCuda) sf->backend = PETSCSF_BACKEND_CUDA;
38520c24465SJunchao Zhang     else if (isKokkos) sf->backend = PETSCSF_BACKEND_KOKKOS;
386*59af0bd3SScott Kruger     else if (isHip) sf->backend = PETSCSF_BACKEND_HIP;
387*59af0bd3SScott Kruger     else if (set) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"-sf_backend %s is not supported. You may choose cuda, hip or kokkos (if installed)", backendstr);
38820c24465SJunchao Zhang   #elif defined(PETSC_HAVE_KOKKOS)
38920c24465SJunchao Zhang     if (set && !isKokkos) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"-sf_backend %s is not supported. You can only choose kokkos", backendstr);
39020c24465SJunchao Zhang   #endif
39120c24465SJunchao Zhang   }
392c2a741eeSJunchao Zhang #endif
393e55864a3SBarry Smith   if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(PetscOptionsObject,sf);CHKERRQ(ierr);}
39495fce210SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
39595fce210SBarry Smith   PetscFunctionReturn(0);
39695fce210SBarry Smith }
39795fce210SBarry Smith 
39829046d53SLisandro Dalcin /*@
39995fce210SBarry Smith    PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order
40095fce210SBarry Smith 
40195fce210SBarry Smith    Logically Collective
40295fce210SBarry Smith 
40395fce210SBarry Smith    Input Arguments:
40495fce210SBarry Smith +  sf - star forest
40595fce210SBarry Smith -  flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic)
40695fce210SBarry Smith 
40795fce210SBarry Smith    Level: advanced
40895fce210SBarry Smith 
40995fce210SBarry Smith .seealso: PetscSFGatherBegin(), PetscSFScatterBegin()
41095fce210SBarry Smith @*/
41195fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg)
41295fce210SBarry Smith {
41395fce210SBarry Smith   PetscFunctionBegin;
41495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
41595fce210SBarry Smith   PetscValidLogicalCollectiveBool(sf,flg,2);
41695fce210SBarry Smith   if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
41795fce210SBarry Smith   sf->rankorder = flg;
41895fce210SBarry Smith   PetscFunctionReturn(0);
41995fce210SBarry Smith }
42095fce210SBarry Smith 
4218af6ec1cSBarry Smith /*@
42295fce210SBarry Smith    PetscSFSetGraph - Set a parallel star forest
42395fce210SBarry Smith 
42495fce210SBarry Smith    Collective
42595fce210SBarry Smith 
42695fce210SBarry Smith    Input Arguments:
42795fce210SBarry Smith +  sf - star forest
42895fce210SBarry Smith .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
42995fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
430c4e6a40aSLawrence Mitchell .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage (locations must be >= 0, enforced
431c4e6a40aSLawrence Mitchell during setup in debug mode)
43295fce210SBarry Smith .  localmode - copy mode for ilocal
433c4e6a40aSLawrence Mitchell .  iremote - remote locations of root vertices for each leaf on the current process (locations must be >= 0, enforced
434c4e6a40aSLawrence Mitchell during setup in debug mode)
43595fce210SBarry Smith -  remotemode - copy mode for iremote
43695fce210SBarry Smith 
43795fce210SBarry Smith    Level: intermediate
43895fce210SBarry Smith 
43995452b02SPatrick Sanan    Notes:
44095452b02SPatrick Sanan     In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode
44138ab3f8aSBarry Smith 
4422ad1e87fSLisandro Dalcin    Developers Note: Local indices which are the identity permutation in the range [0,nleaves) are discarded as they
4432ad1e87fSLisandro Dalcin    encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not
4442ad1e87fSLisandro Dalcin    needed
4452ad1e87fSLisandro Dalcin 
446c4e6a40aSLawrence Mitchell    Developers Note: This object does not necessarily encode a true star forest in the graph theoretic sense, since leaf
447c4e6a40aSLawrence Mitchell    indices are not required to be unique. Some functions, however, rely on unique leaf indices (checked in debug mode).
448c4e6a40aSLawrence Mitchell 
44995fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph()
45095fce210SBarry Smith @*/
45195fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode)
45295fce210SBarry Smith {
45395fce210SBarry Smith   PetscErrorCode ierr;
45495fce210SBarry Smith 
45595fce210SBarry Smith   PetscFunctionBegin;
45695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
45729046d53SLisandro Dalcin   if (nleaves > 0 && ilocal) PetscValidIntPointer(ilocal,4);
45829046d53SLisandro Dalcin   if (nleaves > 0) PetscValidPointer(iremote,6);
45929046d53SLisandro Dalcin   if (nroots  < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nroots %D, cannot be negative",nroots);
46095fce210SBarry Smith   if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %D, cannot be negative",nleaves);
46129046d53SLisandro Dalcin 
4622a67d2daSStefano Zampini   if (sf->nroots >= 0) { /* Reset only if graph already set */
46395fce210SBarry Smith     ierr = PetscSFReset(sf);CHKERRQ(ierr);
4642a67d2daSStefano Zampini   }
4652a67d2daSStefano Zampini 
46629046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
46729046d53SLisandro Dalcin 
46895fce210SBarry Smith   sf->nroots  = nroots;
46995fce210SBarry Smith   sf->nleaves = nleaves;
47029046d53SLisandro Dalcin 
47129046d53SLisandro Dalcin   if (nleaves && ilocal) {
47221c688dcSJed Brown     PetscInt i;
47329046d53SLisandro Dalcin     PetscInt minleaf = PETSC_MAX_INT;
47429046d53SLisandro Dalcin     PetscInt maxleaf = PETSC_MIN_INT;
4752ad1e87fSLisandro Dalcin     int      contiguous = 1;
47629046d53SLisandro Dalcin     for (i=0; i<nleaves; i++) {
47729046d53SLisandro Dalcin       minleaf = PetscMin(minleaf,ilocal[i]);
47829046d53SLisandro Dalcin       maxleaf = PetscMax(maxleaf,ilocal[i]);
4792ad1e87fSLisandro Dalcin       contiguous &= (ilocal[i] == i);
48029046d53SLisandro Dalcin     }
48129046d53SLisandro Dalcin     sf->minleaf = minleaf;
48229046d53SLisandro Dalcin     sf->maxleaf = maxleaf;
4832ad1e87fSLisandro Dalcin     if (contiguous) {
4842ad1e87fSLisandro Dalcin       if (localmode == PETSC_OWN_POINTER) {
4852ad1e87fSLisandro Dalcin         ierr = PetscFree(ilocal);CHKERRQ(ierr);
4862ad1e87fSLisandro Dalcin       }
4872ad1e87fSLisandro Dalcin       ilocal = NULL;
4882ad1e87fSLisandro Dalcin     }
48929046d53SLisandro Dalcin   } else {
49029046d53SLisandro Dalcin     sf->minleaf = 0;
49129046d53SLisandro Dalcin     sf->maxleaf = nleaves - 1;
49229046d53SLisandro Dalcin   }
49329046d53SLisandro Dalcin 
49429046d53SLisandro Dalcin   if (ilocal) {
49595fce210SBarry Smith     switch (localmode) {
49695fce210SBarry Smith     case PETSC_COPY_VALUES:
497785e854fSJed Brown       ierr = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr);
498580bdb30SBarry Smith       ierr = PetscArraycpy(sf->mine_alloc,ilocal,nleaves);CHKERRQ(ierr);
49995fce210SBarry Smith       sf->mine = sf->mine_alloc;
50095fce210SBarry Smith       break;
50195fce210SBarry Smith     case PETSC_OWN_POINTER:
50295fce210SBarry Smith       sf->mine_alloc = (PetscInt*)ilocal;
50395fce210SBarry Smith       sf->mine       = sf->mine_alloc;
50495fce210SBarry Smith       break;
50595fce210SBarry Smith     case PETSC_USE_POINTER:
50629046d53SLisandro Dalcin       sf->mine_alloc = NULL;
50795fce210SBarry Smith       sf->mine       = (PetscInt*)ilocal;
50895fce210SBarry Smith       break;
50995fce210SBarry Smith     default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode");
51095fce210SBarry Smith     }
51195fce210SBarry Smith   }
51229046d53SLisandro Dalcin 
51395fce210SBarry Smith   switch (remotemode) {
51495fce210SBarry Smith   case PETSC_COPY_VALUES:
515785e854fSJed Brown     ierr = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr);
516580bdb30SBarry Smith     ierr = PetscArraycpy(sf->remote_alloc,iremote,nleaves);CHKERRQ(ierr);
51795fce210SBarry Smith     sf->remote = sf->remote_alloc;
51895fce210SBarry Smith     break;
51995fce210SBarry Smith   case PETSC_OWN_POINTER:
52095fce210SBarry Smith     sf->remote_alloc = (PetscSFNode*)iremote;
52195fce210SBarry Smith     sf->remote       = sf->remote_alloc;
52295fce210SBarry Smith     break;
52395fce210SBarry Smith   case PETSC_USE_POINTER:
52429046d53SLisandro Dalcin     sf->remote_alloc = NULL;
52595fce210SBarry Smith     sf->remote       = (PetscSFNode*)iremote;
52695fce210SBarry Smith     break;
52795fce210SBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode");
52895fce210SBarry Smith   }
52995fce210SBarry Smith 
530563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
53129046d53SLisandro Dalcin   sf->graphset = PETSC_TRUE;
53295fce210SBarry Smith   PetscFunctionReturn(0);
53395fce210SBarry Smith }
53495fce210SBarry Smith 
53529046d53SLisandro Dalcin /*@
536dd5b3ca6SJunchao Zhang   PetscSFSetGraphWithPattern - Sets the graph of an SF with a specific pattern
537dd5b3ca6SJunchao Zhang 
538dd5b3ca6SJunchao Zhang   Collective
539dd5b3ca6SJunchao Zhang 
540dd5b3ca6SJunchao Zhang   Input Parameters:
541dd5b3ca6SJunchao Zhang + sf      - The PetscSF
542dd5b3ca6SJunchao Zhang . map     - Layout of roots over all processes (insignificant when pattern is PETSCSF_PATTERN_ALLTOALL)
543dd5b3ca6SJunchao Zhang - pattern - One of PETSCSF_PATTERN_ALLGATHER, PETSCSF_PATTERN_GATHER, PETSCSF_PATTERN_ALLTOALL
544dd5b3ca6SJunchao Zhang 
545dd5b3ca6SJunchao Zhang   Notes:
546dd5b3ca6SJunchao Zhang   It is easier to explain PetscSFPattern using vectors. Suppose we have an MPI vector x and its layout is map.
547dd5b3ca6SJunchao Zhang   n and N are local and global sizes of x respectively.
548dd5b3ca6SJunchao Zhang 
549dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLGATHER, the routine creates a graph that if one does Bcast on it, it will copy x to
550dd5b3ca6SJunchao Zhang   sequential vectors y on all ranks.
551dd5b3ca6SJunchao Zhang 
552dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_GATHER, the routine creates a graph that if one does Bcast on it, it will copy x to a
553dd5b3ca6SJunchao Zhang   sequential vector y on rank 0.
554dd5b3ca6SJunchao Zhang 
555dd5b3ca6SJunchao Zhang   In above cases, entries of x are roots and entries of y are leaves.
556dd5b3ca6SJunchao Zhang 
557dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLTOALL, map is insignificant. Suppose NP is size of sf's communicator. The routine
558dd5b3ca6SJunchao Zhang   creates a graph that every rank has NP leaves and NP roots. On rank i, its leaf j is connected to root i
559dd5b3ca6SJunchao Zhang   of rank j. Here 0 <=i,j<NP. It is a kind of MPI_Alltoall with sendcount/recvcount being 1. Note that it does
560dd5b3ca6SJunchao Zhang   not mean one can not send multiple items. One just needs to create a new MPI datatype for the mulptiple data
561dd5b3ca6SJunchao Zhang   items with MPI_Type_contiguous() and use that as the <unit> argument in SF routines.
562dd5b3ca6SJunchao Zhang 
563dd5b3ca6SJunchao Zhang   In this case, roots and leaves are symmetric.
564dd5b3ca6SJunchao Zhang 
565dd5b3ca6SJunchao Zhang   Level: intermediate
566dd5b3ca6SJunchao Zhang  @*/
567dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFSetGraphWithPattern(PetscSF sf,PetscLayout map,PetscSFPattern pattern)
568dd5b3ca6SJunchao Zhang {
569dd5b3ca6SJunchao Zhang   MPI_Comm       comm;
570dd5b3ca6SJunchao Zhang   PetscInt       n,N,res[2];
571dd5b3ca6SJunchao Zhang   PetscMPIInt    rank,size;
572dd5b3ca6SJunchao Zhang   PetscSFType    type;
573dd5b3ca6SJunchao Zhang   PetscErrorCode ierr;
574dd5b3ca6SJunchao Zhang 
575dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
576dd5b3ca6SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
577dd5b3ca6SJunchao Zhang   if (pattern < PETSCSF_PATTERN_ALLGATHER || pattern > PETSCSF_PATTERN_ALLTOALL) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Unsupported PetscSFPattern %D\n",pattern);
578dd5b3ca6SJunchao Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
579dd5b3ca6SJunchao Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
580dd5b3ca6SJunchao Zhang 
581dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLTOALL) {
582dd5b3ca6SJunchao Zhang     type = PETSCSFALLTOALL;
583dd5b3ca6SJunchao Zhang     ierr = PetscLayoutCreate(comm,&sf->map);CHKERRQ(ierr);
584dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetLocalSize(sf->map,size);CHKERRQ(ierr);
585dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetSize(sf->map,((PetscInt)size)*size);CHKERRQ(ierr);
586dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetUp(sf->map);CHKERRQ(ierr);
587dd5b3ca6SJunchao Zhang   } else {
588dd5b3ca6SJunchao Zhang     ierr   = PetscLayoutGetLocalSize(map,&n);CHKERRQ(ierr);
589dd5b3ca6SJunchao Zhang     ierr   = PetscLayoutGetSize(map,&N);CHKERRQ(ierr);
590dd5b3ca6SJunchao Zhang     res[0] = n;
591dd5b3ca6SJunchao Zhang     res[1] = -n;
592dd5b3ca6SJunchao Zhang     /* Check if n are same over all ranks so that we can optimize it */
593dd5b3ca6SJunchao Zhang     ierr   = MPIU_Allreduce(MPI_IN_PLACE,res,2,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
594dd5b3ca6SJunchao Zhang     if (res[0] == -res[1]) { /* same n */
595dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHER  : PETSCSFGATHER;
596dd5b3ca6SJunchao Zhang     } else {
597dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHERV : PETSCSFGATHERV;
598dd5b3ca6SJunchao Zhang     }
599dd5b3ca6SJunchao Zhang     ierr = PetscLayoutReference(map,&sf->map);CHKERRQ(ierr);
600dd5b3ca6SJunchao Zhang   }
601dd5b3ca6SJunchao Zhang   ierr = PetscSFSetType(sf,type);CHKERRQ(ierr);
602dd5b3ca6SJunchao Zhang 
603dd5b3ca6SJunchao Zhang   sf->pattern = pattern;
604dd5b3ca6SJunchao Zhang   sf->mine    = NULL; /* Contiguous */
605dd5b3ca6SJunchao Zhang 
606dd5b3ca6SJunchao Zhang   /* Set nleaves, nroots here in case user calls PetscSFGetGraph, which is legal to call even before PetscSFSetUp is called.
607dd5b3ca6SJunchao Zhang      Also set other easy stuff.
608dd5b3ca6SJunchao Zhang    */
609dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLGATHER) {
610dd5b3ca6SJunchao Zhang     sf->nleaves      = N;
611dd5b3ca6SJunchao Zhang     sf->nroots       = n;
612dd5b3ca6SJunchao Zhang     sf->nranks       = size;
613dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
614dd5b3ca6SJunchao Zhang     sf->maxleaf      = N - 1;
615dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_GATHER) {
616dd5b3ca6SJunchao Zhang     sf->nleaves      = rank ? 0 : N;
617dd5b3ca6SJunchao Zhang     sf->nroots       = n;
618dd5b3ca6SJunchao Zhang     sf->nranks       = rank ? 0 : size;
619dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
620dd5b3ca6SJunchao Zhang     sf->maxleaf      = rank ? -1 : N - 1;
621dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_ALLTOALL) {
622dd5b3ca6SJunchao Zhang     sf->nleaves      = size;
623dd5b3ca6SJunchao Zhang     sf->nroots       = size;
624dd5b3ca6SJunchao Zhang     sf->nranks       = size;
625dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
626dd5b3ca6SJunchao Zhang     sf->maxleaf      = size - 1;
627dd5b3ca6SJunchao Zhang   }
628dd5b3ca6SJunchao Zhang   sf->ndranks  = 0; /* We do not need to separate out distinguished ranks for patterned graphs to improve communication performance */
629dd5b3ca6SJunchao Zhang   sf->graphset = PETSC_TRUE;
630dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
631dd5b3ca6SJunchao Zhang }
632dd5b3ca6SJunchao Zhang 
633dd5b3ca6SJunchao Zhang /*@
63495fce210SBarry Smith    PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map
63595fce210SBarry Smith 
63695fce210SBarry Smith    Collective
63795fce210SBarry Smith 
63895fce210SBarry Smith    Input Arguments:
639dd5b3ca6SJunchao Zhang 
64095fce210SBarry Smith .  sf - star forest to invert
64195fce210SBarry Smith 
64295fce210SBarry Smith    Output Arguments:
64395fce210SBarry Smith .  isf - inverse of sf
64495fce210SBarry Smith    Level: advanced
64595fce210SBarry Smith 
64695fce210SBarry Smith    Notes:
64795fce210SBarry Smith    All roots must have degree 1.
64895fce210SBarry Smith 
64995fce210SBarry Smith    The local space may be a permutation, but cannot be sparse.
65095fce210SBarry Smith 
65195fce210SBarry Smith .seealso: PetscSFSetGraph()
65295fce210SBarry Smith @*/
65395fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
65495fce210SBarry Smith {
65595fce210SBarry Smith   PetscErrorCode ierr;
65695fce210SBarry Smith   PetscMPIInt    rank;
65795fce210SBarry Smith   PetscInt       i,nroots,nleaves,maxlocal,count,*newilocal;
65895fce210SBarry Smith   const PetscInt *ilocal;
65995fce210SBarry Smith   PetscSFNode    *roots,*leaves;
66095fce210SBarry Smith 
66195fce210SBarry Smith   PetscFunctionBegin;
66229046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
66329046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
66429046d53SLisandro Dalcin   PetscValidPointer(isf,2);
66529046d53SLisandro Dalcin 
66695fce210SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
66729046d53SLisandro Dalcin   maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
66829046d53SLisandro Dalcin 
66929046d53SLisandro Dalcin   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
670ae9aee6dSMatthew G. Knepley   ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr);
671ae9aee6dSMatthew G. Knepley   for (i=0; i<maxlocal; i++) {
67295fce210SBarry Smith     leaves[i].rank  = rank;
67395fce210SBarry Smith     leaves[i].index = i;
67495fce210SBarry Smith   }
67595fce210SBarry Smith   for (i=0; i <nroots; i++) {
67695fce210SBarry Smith     roots[i].rank  = -1;
67795fce210SBarry Smith     roots[i].index = -1;
67895fce210SBarry Smith   }
6798bfbc91cSJed Brown   ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
6808bfbc91cSJed Brown   ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPIU_REPLACE);CHKERRQ(ierr);
68195fce210SBarry Smith 
68295fce210SBarry Smith   /* Check whether our leaves are sparse */
68395fce210SBarry Smith   for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
68495fce210SBarry Smith   if (count == nroots) newilocal = NULL;
68595fce210SBarry Smith   else {                        /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
686785e854fSJed Brown     ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr);
68795fce210SBarry Smith     for (i=0,count=0; i<nroots; i++) {
68895fce210SBarry Smith       if (roots[i].rank >= 0) {
68995fce210SBarry Smith         newilocal[count]   = i;
69095fce210SBarry Smith         roots[count].rank  = roots[i].rank;
69195fce210SBarry Smith         roots[count].index = roots[i].index;
69295fce210SBarry Smith         count++;
69395fce210SBarry Smith       }
69495fce210SBarry Smith     }
69595fce210SBarry Smith   }
69695fce210SBarry Smith 
69795fce210SBarry Smith   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr);
69895fce210SBarry Smith   ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr);
69995fce210SBarry Smith   ierr = PetscFree2(roots,leaves);CHKERRQ(ierr);
70095fce210SBarry Smith   PetscFunctionReturn(0);
70195fce210SBarry Smith }
70295fce210SBarry Smith 
70395fce210SBarry Smith /*@
70495fce210SBarry Smith    PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph
70595fce210SBarry Smith 
70695fce210SBarry Smith    Collective
70795fce210SBarry Smith 
70895fce210SBarry Smith    Input Arguments:
70995fce210SBarry Smith +  sf - communication object to duplicate
71095fce210SBarry Smith -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)
71195fce210SBarry Smith 
71295fce210SBarry Smith    Output Arguments:
71395fce210SBarry Smith .  newsf - new communication object
71495fce210SBarry Smith 
71595fce210SBarry Smith    Level: beginner
71695fce210SBarry Smith 
71795fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph()
71895fce210SBarry Smith @*/
71995fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
72095fce210SBarry Smith {
72129046d53SLisandro Dalcin   PetscSFType    type;
72295fce210SBarry Smith   PetscErrorCode ierr;
72397929ea7SJunchao Zhang   MPI_Datatype   dtype=MPIU_SCALAR;
72495fce210SBarry Smith 
72595fce210SBarry Smith   PetscFunctionBegin;
72629046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
72729046d53SLisandro Dalcin   PetscValidLogicalCollectiveEnum(sf,opt,2);
72829046d53SLisandro Dalcin   PetscValidPointer(newsf,3);
72995fce210SBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr);
73029046d53SLisandro Dalcin   ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr);
73129046d53SLisandro Dalcin   if (type) {ierr = PetscSFSetType(*newsf,type);CHKERRQ(ierr);}
73295fce210SBarry Smith   if (opt == PETSCSF_DUPLICATE_GRAPH) {
733dd5b3ca6SJunchao Zhang     PetscSFCheckGraphSet(sf,1);
734dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
73595fce210SBarry Smith       PetscInt          nroots,nleaves;
73695fce210SBarry Smith       const PetscInt    *ilocal;
73795fce210SBarry Smith       const PetscSFNode *iremote;
73895fce210SBarry Smith       ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
73995fce210SBarry Smith       ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr);
740dd5b3ca6SJunchao Zhang     } else {
741dd5b3ca6SJunchao Zhang       ierr = PetscSFSetGraphWithPattern(*newsf,sf->map,sf->pattern);CHKERRQ(ierr);
742dd5b3ca6SJunchao Zhang     }
74395fce210SBarry Smith   }
74497929ea7SJunchao Zhang   /* Since oldtype is committed, so is newtype, according to MPI */
74597929ea7SJunchao Zhang   if (sf->vscat.bs > 1) {ierr = MPI_Type_dup(sf->vscat.unit,&dtype);CHKERRQ(ierr);}
74697929ea7SJunchao Zhang   (*newsf)->vscat.bs     = sf->vscat.bs;
74797929ea7SJunchao Zhang   (*newsf)->vscat.unit   = dtype;
74897929ea7SJunchao Zhang   (*newsf)->vscat.to_n   = sf->vscat.to_n;
74997929ea7SJunchao Zhang   (*newsf)->vscat.from_n = sf->vscat.from_n;
75097929ea7SJunchao Zhang   /* Do not copy lsf. Build it on demand since it is rarely used */
75197929ea7SJunchao Zhang 
75220c24465SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
75320c24465SJunchao Zhang   (*newsf)->backend              = sf->backend;
75420c24465SJunchao Zhang   (*newsf)->use_default_stream   = sf->use_default_stream;
75520c24465SJunchao Zhang   (*newsf)->use_gpu_aware_mpi    = sf->use_gpu_aware_mpi;
75620c24465SJunchao Zhang   (*newsf)->use_stream_aware_mpi = sf->use_stream_aware_mpi;
75720c24465SJunchao Zhang #endif
75829046d53SLisandro Dalcin   if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);}
75920c24465SJunchao Zhang   /* Don't do PetscSFSetUp() since the new sf's graph might have not been set. */
76095fce210SBarry Smith   PetscFunctionReturn(0);
76195fce210SBarry Smith }
76295fce210SBarry Smith 
76395fce210SBarry Smith /*@C
76495fce210SBarry Smith    PetscSFGetGraph - Get the graph specifying a parallel star forest
76595fce210SBarry Smith 
76695fce210SBarry Smith    Not Collective
76795fce210SBarry Smith 
76895fce210SBarry Smith    Input Arguments:
76995fce210SBarry Smith .  sf - star forest
77095fce210SBarry Smith 
77195fce210SBarry Smith    Output Arguments:
77295fce210SBarry Smith +  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
77395fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
77495fce210SBarry Smith .  ilocal - locations of leaves in leafdata buffers
77595fce210SBarry Smith -  iremote - remote locations of root vertices for each leaf on the current process
77695fce210SBarry Smith 
777373e0d91SLisandro Dalcin    Notes:
778373e0d91SLisandro Dalcin    We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set yet
779373e0d91SLisandro Dalcin 
780245d9833Sprj-    When called from Fortran, the returned iremote array is a copy and must be deallocated after use. Consequently, if you
781ca797d7aSLawrence Mitchell    want to update the graph, you must call PetscSFSetGraph after modifying the iremote array.
782ca797d7aSLawrence Mitchell 
78395fce210SBarry Smith    Level: intermediate
78495fce210SBarry Smith 
78595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph()
78695fce210SBarry Smith @*/
78795fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
78895fce210SBarry Smith {
789b8dee149SJunchao Zhang   PetscErrorCode ierr;
79095fce210SBarry Smith 
79195fce210SBarry Smith   PetscFunctionBegin;
79295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
793b8dee149SJunchao Zhang   if (sf->ops->GetGraph) {
794b8dee149SJunchao Zhang     ierr = (sf->ops->GetGraph)(sf,nroots,nleaves,ilocal,iremote);CHKERRQ(ierr);
795b8dee149SJunchao Zhang   } else {
79695fce210SBarry Smith     if (nroots) *nroots = sf->nroots;
79795fce210SBarry Smith     if (nleaves) *nleaves = sf->nleaves;
79895fce210SBarry Smith     if (ilocal) *ilocal = sf->mine;
79995fce210SBarry Smith     if (iremote) *iremote = sf->remote;
800b8dee149SJunchao Zhang   }
80195fce210SBarry Smith   PetscFunctionReturn(0);
80295fce210SBarry Smith }
80395fce210SBarry Smith 
80429046d53SLisandro Dalcin /*@
80595fce210SBarry Smith    PetscSFGetLeafRange - Get the active leaf ranges
80695fce210SBarry Smith 
80795fce210SBarry Smith    Not Collective
80895fce210SBarry Smith 
80995fce210SBarry Smith    Input Arguments:
81095fce210SBarry Smith .  sf - star forest
81195fce210SBarry Smith 
81295fce210SBarry Smith    Output Arguments:
813dd5b3ca6SJunchao Zhang +  minleaf - minimum active leaf on this process. Return 0 if there are no leaves.
814dd5b3ca6SJunchao Zhang -  maxleaf - maximum active leaf on this process. Return -1 if there are no leaves.
81595fce210SBarry Smith 
81695fce210SBarry Smith    Level: developer
81795fce210SBarry Smith 
81895fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
81995fce210SBarry Smith @*/
82095fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
82195fce210SBarry Smith {
82295fce210SBarry Smith   PetscFunctionBegin;
82395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
82429046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
82595fce210SBarry Smith   if (minleaf) *minleaf = sf->minleaf;
82695fce210SBarry Smith   if (maxleaf) *maxleaf = sf->maxleaf;
82795fce210SBarry Smith   PetscFunctionReturn(0);
82895fce210SBarry Smith }
82995fce210SBarry Smith 
83095fce210SBarry Smith /*@C
831fe2efc57SMark    PetscSFViewFromOptions - View from Options
832fe2efc57SMark 
833fe2efc57SMark    Collective on PetscSF
834fe2efc57SMark 
835fe2efc57SMark    Input Parameters:
836fe2efc57SMark +  A - the star forest
837736c3998SJose E. Roman .  obj - Optional object
838736c3998SJose E. Roman -  name - command line option
839fe2efc57SMark 
840fe2efc57SMark    Level: intermediate
841fe2efc57SMark .seealso:  PetscSF, PetscSFView, PetscObjectViewFromOptions(), PetscSFCreate()
842fe2efc57SMark @*/
843fe2efc57SMark PetscErrorCode  PetscSFViewFromOptions(PetscSF A,PetscObject obj,const char name[])
844fe2efc57SMark {
845fe2efc57SMark   PetscErrorCode ierr;
846fe2efc57SMark 
847fe2efc57SMark   PetscFunctionBegin;
848fe2efc57SMark   PetscValidHeaderSpecific(A,PETSCSF_CLASSID,1);
849fe2efc57SMark   ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr);
850fe2efc57SMark   PetscFunctionReturn(0);
851fe2efc57SMark }
852fe2efc57SMark 
853fe2efc57SMark /*@C
85495fce210SBarry Smith    PetscSFView - view a star forest
85595fce210SBarry Smith 
85695fce210SBarry Smith    Collective
85795fce210SBarry Smith 
85895fce210SBarry Smith    Input Arguments:
85995fce210SBarry Smith +  sf - star forest
86095fce210SBarry Smith -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD
86195fce210SBarry Smith 
86295fce210SBarry Smith    Level: beginner
86395fce210SBarry Smith 
86495fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph()
86595fce210SBarry Smith @*/
86695fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
86795fce210SBarry Smith {
86895fce210SBarry Smith   PetscErrorCode    ierr;
86995fce210SBarry Smith   PetscBool         iascii;
87095fce210SBarry Smith   PetscViewerFormat format;
87195fce210SBarry Smith 
87295fce210SBarry Smith   PetscFunctionBegin;
87395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
87495fce210SBarry Smith   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);}
87595fce210SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
87695fce210SBarry Smith   PetscCheckSameComm(sf,1,viewer,2);
87780153354SVaclav Hapla   if (sf->graphset) {ierr = PetscSFSetUp(sf);CHKERRQ(ierr);}
87895fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
87995fce210SBarry Smith   if (iascii) {
88095fce210SBarry Smith     PetscMPIInt rank;
88181bfa7aaSJed Brown     PetscInt    ii,i,j;
88295fce210SBarry Smith 
883dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr);
88495fce210SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
88595fce210SBarry Smith     if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);}
886dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
88780153354SVaclav Hapla       if (!sf->graphset) {
88880153354SVaclav Hapla         ierr = PetscViewerASCIIPrintf(viewer,"PetscSFSetGraph() has not been called yet\n");CHKERRQ(ierr);
88980153354SVaclav Hapla         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
89080153354SVaclav Hapla         PetscFunctionReturn(0);
89180153354SVaclav Hapla       }
89295fce210SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
8931575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
89495fce210SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr);
89595fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
89695fce210SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- (%D,%D)\n",rank,sf->mine ? sf->mine[i] : i,sf->remote[i].rank,sf->remote[i].index);CHKERRQ(ierr);
89795fce210SBarry Smith       }
89895fce210SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
89995fce210SBarry Smith       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
90095fce210SBarry Smith       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
90181bfa7aaSJed Brown         PetscMPIInt *tmpranks,*perm;
90281bfa7aaSJed Brown         ierr = PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm);CHKERRQ(ierr);
903580bdb30SBarry Smith         ierr = PetscArraycpy(tmpranks,sf->ranks,sf->nranks);CHKERRQ(ierr);
90481bfa7aaSJed Brown         for (i=0; i<sf->nranks; i++) perm[i] = i;
90581bfa7aaSJed Brown         ierr = PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm);CHKERRQ(ierr);
90695fce210SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr);
90781bfa7aaSJed Brown         for (ii=0; ii<sf->nranks; ii++) {
90881bfa7aaSJed Brown           i = perm[ii];
9097904a332SBarry Smith           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %D edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr);
91095fce210SBarry Smith           for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
91195fce210SBarry Smith             ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d]    %D <- %D\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr);
91295fce210SBarry Smith           }
91395fce210SBarry Smith         }
91481bfa7aaSJed Brown         ierr = PetscFree2(tmpranks,perm);CHKERRQ(ierr);
91595fce210SBarry Smith       }
91695fce210SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
9171575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
918dd5b3ca6SJunchao Zhang     }
91995fce210SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
92095fce210SBarry Smith   }
92195fce210SBarry Smith   PetscFunctionReturn(0);
92295fce210SBarry Smith }
92395fce210SBarry Smith 
92495fce210SBarry Smith /*@C
925dec1416fSJunchao Zhang    PetscSFGetRootRanks - Get root ranks and number of vertices referenced by leaves on this process
92695fce210SBarry Smith 
92795fce210SBarry Smith    Not Collective
92895fce210SBarry Smith 
92995fce210SBarry Smith    Input Arguments:
93095fce210SBarry Smith .  sf - star forest
93195fce210SBarry Smith 
93295fce210SBarry Smith    Output Arguments:
93395fce210SBarry Smith +  nranks - number of ranks referenced by local part
93495fce210SBarry Smith .  ranks - array of ranks
93595fce210SBarry Smith .  roffset - offset in rmine/rremote for each rank (length nranks+1)
93695fce210SBarry Smith .  rmine - concatenated array holding local indices referencing each remote rank
93795fce210SBarry Smith -  rremote - concatenated array holding remote indices referenced for each remote rank
93895fce210SBarry Smith 
93995fce210SBarry Smith    Level: developer
94095fce210SBarry Smith 
941dec1416fSJunchao Zhang .seealso: PetscSFGetLeafRanks()
94295fce210SBarry Smith @*/
943dec1416fSJunchao Zhang PetscErrorCode PetscSFGetRootRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
94495fce210SBarry Smith {
945dec1416fSJunchao Zhang   PetscErrorCode ierr;
94695fce210SBarry Smith 
94795fce210SBarry Smith   PetscFunctionBegin;
94895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
94929046d53SLisandro Dalcin   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
950dec1416fSJunchao Zhang   if (sf->ops->GetRootRanks) {
951dec1416fSJunchao Zhang     ierr = (sf->ops->GetRootRanks)(sf,nranks,ranks,roffset,rmine,rremote);CHKERRQ(ierr);
952dec1416fSJunchao Zhang   } else {
953dec1416fSJunchao Zhang     /* The generic implementation */
95495fce210SBarry Smith     if (nranks)  *nranks  = sf->nranks;
95595fce210SBarry Smith     if (ranks)   *ranks   = sf->ranks;
95695fce210SBarry Smith     if (roffset) *roffset = sf->roffset;
95795fce210SBarry Smith     if (rmine)   *rmine   = sf->rmine;
95895fce210SBarry Smith     if (rremote) *rremote = sf->rremote;
959dec1416fSJunchao Zhang   }
96095fce210SBarry Smith   PetscFunctionReturn(0);
96195fce210SBarry Smith }
96295fce210SBarry Smith 
9638750ddebSJunchao Zhang /*@C
9648750ddebSJunchao Zhang    PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process
9658750ddebSJunchao Zhang 
9668750ddebSJunchao Zhang    Not Collective
9678750ddebSJunchao Zhang 
9688750ddebSJunchao Zhang    Input Arguments:
9698750ddebSJunchao Zhang .  sf - star forest
9708750ddebSJunchao Zhang 
9718750ddebSJunchao Zhang    Output Arguments:
9728750ddebSJunchao Zhang +  niranks - number of leaf ranks referencing roots on this process
9738750ddebSJunchao Zhang .  iranks - array of ranks
9748750ddebSJunchao Zhang .  ioffset - offset in irootloc for each rank (length niranks+1)
9758750ddebSJunchao Zhang -  irootloc - concatenated array holding local indices of roots referenced by each leaf rank
9768750ddebSJunchao Zhang 
9778750ddebSJunchao Zhang    Level: developer
9788750ddebSJunchao Zhang 
979dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks()
9808750ddebSJunchao Zhang @*/
9818750ddebSJunchao Zhang PetscErrorCode PetscSFGetLeafRanks(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
9828750ddebSJunchao Zhang {
9838750ddebSJunchao Zhang   PetscErrorCode ierr;
9848750ddebSJunchao Zhang 
9858750ddebSJunchao Zhang   PetscFunctionBegin;
9868750ddebSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
9878750ddebSJunchao Zhang   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
9888750ddebSJunchao Zhang   if (sf->ops->GetLeafRanks) {
9898750ddebSJunchao Zhang     ierr = (sf->ops->GetLeafRanks)(sf,niranks,iranks,ioffset,irootloc);CHKERRQ(ierr);
9908750ddebSJunchao Zhang   } else {
9918750ddebSJunchao Zhang     PetscSFType type;
9928750ddebSJunchao Zhang     ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr);
9938750ddebSJunchao Zhang     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type);
9948750ddebSJunchao Zhang   }
9958750ddebSJunchao Zhang   PetscFunctionReturn(0);
9968750ddebSJunchao Zhang }
9978750ddebSJunchao Zhang 
998b5a8e515SJed Brown static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) {
999b5a8e515SJed Brown   PetscInt i;
1000b5a8e515SJed Brown   for (i=0; i<n; i++) {
1001b5a8e515SJed Brown     if (needle == list[i]) return PETSC_TRUE;
1002b5a8e515SJed Brown   }
1003b5a8e515SJed Brown   return PETSC_FALSE;
1004b5a8e515SJed Brown }
1005b5a8e515SJed Brown 
100695fce210SBarry Smith /*@C
100721c688dcSJed Brown    PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations.
100821c688dcSJed Brown 
100921c688dcSJed Brown    Collective
101021c688dcSJed Brown 
101121c688dcSJed Brown    Input Arguments:
1012b5a8e515SJed Brown +  sf - PetscSF to set up; PetscSFSetGraph() must have been called
1013b5a8e515SJed Brown -  dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange)
101421c688dcSJed Brown 
101521c688dcSJed Brown    Level: developer
101621c688dcSJed Brown 
1017dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks()
101821c688dcSJed Brown @*/
1019b5a8e515SJed Brown PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup)
102021c688dcSJed Brown {
102121c688dcSJed Brown   PetscErrorCode     ierr;
102221c688dcSJed Brown   PetscTable         table;
102321c688dcSJed Brown   PetscTablePosition pos;
1024b5a8e515SJed Brown   PetscMPIInt        size,groupsize,*groupranks;
1025247e8311SStefano Zampini   PetscInt           *rcount,*ranks;
1026247e8311SStefano Zampini   PetscInt           i, irank = -1,orank = -1;
102721c688dcSJed Brown 
102821c688dcSJed Brown   PetscFunctionBegin;
102921c688dcSJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
103029046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
103121c688dcSJed Brown   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRQ(ierr);
103221c688dcSJed Brown   ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr);
103321c688dcSJed Brown   for (i=0; i<sf->nleaves; i++) {
103421c688dcSJed Brown     /* Log 1-based rank */
103521c688dcSJed Brown     ierr = PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr);
103621c688dcSJed Brown   }
103721c688dcSJed Brown   ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr);
103821c688dcSJed Brown   ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr);
103921c688dcSJed Brown   ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr);
104021c688dcSJed Brown   ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr);
104121c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
104221c688dcSJed Brown     ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr);
104321c688dcSJed Brown     ranks[i]--;             /* Convert back to 0-based */
104421c688dcSJed Brown   }
104521c688dcSJed Brown   ierr = PetscTableDestroy(&table);CHKERRQ(ierr);
1046b5a8e515SJed Brown 
1047b5a8e515SJed Brown   /* We expect that dgroup is reliably "small" while nranks could be large */
1048b5a8e515SJed Brown   {
10497fb8a5e4SKarl Rupp     MPI_Group group = MPI_GROUP_NULL;
1050b5a8e515SJed Brown     PetscMPIInt *dgroupranks;
1051b5a8e515SJed Brown     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
1052b5a8e515SJed Brown     ierr = MPI_Group_size(dgroup,&groupsize);CHKERRQ(ierr);
1053b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&dgroupranks);CHKERRQ(ierr);
1054b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&groupranks);CHKERRQ(ierr);
1055b5a8e515SJed Brown     for (i=0; i<groupsize; i++) dgroupranks[i] = i;
1056f3fc9a17SJed Brown     if (groupsize) {ierr = MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);CHKERRQ(ierr);}
1057b5a8e515SJed Brown     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
1058b5a8e515SJed Brown     ierr = PetscFree(dgroupranks);CHKERRQ(ierr);
1059b5a8e515SJed Brown   }
1060b5a8e515SJed Brown 
1061b5a8e515SJed Brown   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
1062b5a8e515SJed Brown   for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i;) {
1063b5a8e515SJed Brown     for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */
1064b5a8e515SJed Brown       if (InList(ranks[i],groupsize,groupranks)) break;
1065b5a8e515SJed Brown     }
1066b5a8e515SJed Brown     for (; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
1067b5a8e515SJed Brown       if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break;
1068b5a8e515SJed Brown     }
1069b5a8e515SJed Brown     if (sf->ndranks < i) {                         /* Swap ranks[sf->ndranks] with ranks[i] */
1070b5a8e515SJed Brown       PetscInt tmprank,tmpcount;
1071247e8311SStefano Zampini 
1072b5a8e515SJed Brown       tmprank             = ranks[i];
1073b5a8e515SJed Brown       tmpcount            = rcount[i];
1074b5a8e515SJed Brown       ranks[i]            = ranks[sf->ndranks];
1075b5a8e515SJed Brown       rcount[i]           = rcount[sf->ndranks];
1076b5a8e515SJed Brown       ranks[sf->ndranks]  = tmprank;
1077b5a8e515SJed Brown       rcount[sf->ndranks] = tmpcount;
1078b5a8e515SJed Brown       sf->ndranks++;
1079b5a8e515SJed Brown     }
1080b5a8e515SJed Brown   }
1081b5a8e515SJed Brown   ierr = PetscFree(groupranks);CHKERRQ(ierr);
1082b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->ndranks,ranks,rcount);CHKERRQ(ierr);
1083b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);CHKERRQ(ierr);
108421c688dcSJed Brown   sf->roffset[0] = 0;
108521c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
108621c688dcSJed Brown     ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr);
108721c688dcSJed Brown     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
108821c688dcSJed Brown     rcount[i]        = 0;
108921c688dcSJed Brown   }
1090247e8311SStefano Zampini   for (i=0, irank = -1, orank = -1; i<sf->nleaves; i++) {
1091247e8311SStefano Zampini     /* short circuit */
1092247e8311SStefano Zampini     if (orank != sf->remote[i].rank) {
109321c688dcSJed Brown       /* Search for index of iremote[i].rank in sf->ranks */
1094b5a8e515SJed Brown       ierr = PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);CHKERRQ(ierr);
1095b5a8e515SJed Brown       if (irank < 0) {
1096b5a8e515SJed Brown         ierr = PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);CHKERRQ(ierr);
1097b5a8e515SJed Brown         if (irank >= 0) irank += sf->ndranks;
109821c688dcSJed Brown       }
1099247e8311SStefano Zampini       orank = sf->remote[i].rank;
1100247e8311SStefano Zampini     }
1101b5a8e515SJed Brown     if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %D in array",sf->remote[i].rank);
110221c688dcSJed Brown     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
110321c688dcSJed Brown     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
110421c688dcSJed Brown     rcount[irank]++;
110521c688dcSJed Brown   }
110621c688dcSJed Brown   ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr);
110721c688dcSJed Brown   PetscFunctionReturn(0);
110821c688dcSJed Brown }
110921c688dcSJed Brown 
111021c688dcSJed Brown /*@C
111195fce210SBarry Smith    PetscSFGetGroups - gets incoming and outgoing process groups
111295fce210SBarry Smith 
111395fce210SBarry Smith    Collective
111495fce210SBarry Smith 
111595fce210SBarry Smith    Input Argument:
111695fce210SBarry Smith .  sf - star forest
111795fce210SBarry Smith 
111895fce210SBarry Smith    Output Arguments:
111995fce210SBarry Smith +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
112095fce210SBarry Smith -  outgoing - group of destination processes for outgoing edges (roots that I reference)
112195fce210SBarry Smith 
112295fce210SBarry Smith    Level: developer
112395fce210SBarry Smith 
112495fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
112595fce210SBarry Smith @*/
112695fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
112795fce210SBarry Smith {
112895fce210SBarry Smith   PetscErrorCode ierr;
11297fb8a5e4SKarl Rupp   MPI_Group      group = MPI_GROUP_NULL;
113095fce210SBarry Smith 
113195fce210SBarry Smith   PetscFunctionBegin;
113244ee17edSStefano Zampini   if (sf->nranks < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUpRanks() before obtaining groups");
113395fce210SBarry Smith   if (sf->ingroup == MPI_GROUP_NULL) {
113495fce210SBarry Smith     PetscInt       i;
113595fce210SBarry Smith     const PetscInt *indegree;
113695fce210SBarry Smith     PetscMPIInt    rank,*outranks,*inranks;
113795fce210SBarry Smith     PetscSFNode    *remote;
113895fce210SBarry Smith     PetscSF        bgcount;
113995fce210SBarry Smith 
114095fce210SBarry Smith     /* Compute the number of incoming ranks */
1141785e854fSJed Brown     ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr);
114295fce210SBarry Smith     for (i=0; i<sf->nranks; i++) {
114395fce210SBarry Smith       remote[i].rank  = sf->ranks[i];
114495fce210SBarry Smith       remote[i].index = 0;
114595fce210SBarry Smith     }
114695fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr);
114795fce210SBarry Smith     ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
114895fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr);
114995fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr);
115095fce210SBarry Smith     /* Enumerate the incoming ranks */
1151dcca6d9dSJed Brown     ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr);
115295fce210SBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
115395fce210SBarry Smith     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
115495fce210SBarry Smith     ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
115595fce210SBarry Smith     ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
115695fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
115795fce210SBarry Smith     ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRQ(ierr);
115895fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
115995fce210SBarry Smith     ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr);
116095fce210SBarry Smith     ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr);
116195fce210SBarry Smith   }
116295fce210SBarry Smith   *incoming = sf->ingroup;
116395fce210SBarry Smith 
116495fce210SBarry Smith   if (sf->outgroup == MPI_GROUP_NULL) {
116595fce210SBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRQ(ierr);
116695fce210SBarry Smith     ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRQ(ierr);
116795fce210SBarry Smith     ierr = MPI_Group_free(&group);CHKERRQ(ierr);
116895fce210SBarry Smith   }
116995fce210SBarry Smith   *outgoing = sf->outgroup;
117095fce210SBarry Smith   PetscFunctionReturn(0);
117195fce210SBarry Smith }
117295fce210SBarry Smith 
117329046d53SLisandro Dalcin /*@
117495fce210SBarry Smith    PetscSFGetMultiSF - gets the inner SF implemeting gathers and scatters
117595fce210SBarry Smith 
117695fce210SBarry Smith    Collective
117795fce210SBarry Smith 
117895fce210SBarry Smith    Input Argument:
117995fce210SBarry Smith .  sf - star forest that may contain roots with 0 or with more than 1 vertex
118095fce210SBarry Smith 
118195fce210SBarry Smith    Output Arguments:
118295fce210SBarry Smith .  multi - star forest with split roots, such that each root has degree exactly 1
118395fce210SBarry Smith 
118495fce210SBarry Smith    Level: developer
118595fce210SBarry Smith 
118695fce210SBarry Smith    Notes:
118795fce210SBarry Smith 
118895fce210SBarry Smith    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
118995fce210SBarry Smith    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
119095fce210SBarry Smith    edge, it is a candidate for future optimization that might involve its removal.
119195fce210SBarry Smith 
1192673100f5SVaclav Hapla .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin(), PetscSFComputeMultiRootOriginalNumbering()
119395fce210SBarry Smith @*/
119495fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
119595fce210SBarry Smith {
119695fce210SBarry Smith   PetscErrorCode ierr;
119795fce210SBarry Smith 
119895fce210SBarry Smith   PetscFunctionBegin;
119995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
120095fce210SBarry Smith   PetscValidPointer(multi,2);
120195fce210SBarry Smith   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
120295fce210SBarry Smith     ierr   = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
120395fce210SBarry Smith     *multi = sf->multi;
120495fce210SBarry Smith     PetscFunctionReturn(0);
120595fce210SBarry Smith   }
120695fce210SBarry Smith   if (!sf->multi) {
120795fce210SBarry Smith     const PetscInt *indegree;
12089837ea96SMatthew G. Knepley     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
120995fce210SBarry Smith     PetscSFNode    *remote;
121029046d53SLisandro Dalcin     maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
121195fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr);
121295fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr);
12139837ea96SMatthew G. Knepley     ierr = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr);
121495fce210SBarry Smith     inoffset[0] = 0;
121595fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
12169837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) outones[i] = 1;
1217dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
1218dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
121995fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
122076bd3646SJed Brown     if (PetscDefined(USE_DEBUG)) { /* Check that the expected number of increments occurred */
122195fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
122295fce210SBarry Smith         if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
122395fce210SBarry Smith       }
122476bd3646SJed Brown     }
1225785e854fSJed Brown     ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr);
122695fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
122795fce210SBarry Smith       remote[i].rank  = sf->remote[i].rank;
122838e7336fSToby Isaac       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
122995fce210SBarry Smith     }
123095fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
123101365b40SToby Isaac     ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
123295fce210SBarry Smith     if (sf->rankorder) {        /* Sort the ranks */
123395fce210SBarry Smith       PetscMPIInt rank;
123495fce210SBarry Smith       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
123595fce210SBarry Smith       PetscSFNode *newremote;
123695fce210SBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRQ(ierr);
123795fce210SBarry Smith       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
12389837ea96SMatthew G. Knepley       ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr);
12399837ea96SMatthew G. Knepley       for (i=0; i<maxlocal; i++) outranks[i] = rank;
12408bfbc91cSJed Brown       ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
12418bfbc91cSJed Brown       ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPIU_REPLACE);CHKERRQ(ierr);
124295fce210SBarry Smith       /* Sort the incoming ranks at each vertex, build the inverse map */
124395fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
124495fce210SBarry Smith         PetscInt j;
124595fce210SBarry Smith         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
124695fce210SBarry Smith         ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr);
124795fce210SBarry Smith         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
124895fce210SBarry Smith       }
124995fce210SBarry Smith       ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
125095fce210SBarry Smith       ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset);CHKERRQ(ierr);
1251785e854fSJed Brown       ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr);
125295fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
125395fce210SBarry Smith         newremote[i].rank  = sf->remote[i].rank;
125401365b40SToby Isaac         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
125595fce210SBarry Smith       }
125601365b40SToby Isaac       ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
125795fce210SBarry Smith       ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr);
125895fce210SBarry Smith     }
125995fce210SBarry Smith     ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr);
126095fce210SBarry Smith   }
126195fce210SBarry Smith   *multi = sf->multi;
126295fce210SBarry Smith   PetscFunctionReturn(0);
126395fce210SBarry Smith }
126495fce210SBarry Smith 
126595fce210SBarry Smith /*@C
126695fce210SBarry Smith    PetscSFCreateEmbeddedSF - removes edges from all but the selected roots, does not remap indices
126795fce210SBarry Smith 
126895fce210SBarry Smith    Collective
126995fce210SBarry Smith 
127095fce210SBarry Smith    Input Arguments:
127195fce210SBarry Smith +  sf - original star forest
1272ba2a7774SJunchao Zhang .  nselected  - number of selected roots on this process
1273ba2a7774SJunchao Zhang -  selected   - indices of the selected roots on this process
127495fce210SBarry Smith 
127595fce210SBarry Smith    Output Arguments:
1276cd620004SJunchao Zhang .  esf - new star forest
127795fce210SBarry Smith 
127895fce210SBarry Smith    Level: advanced
127995fce210SBarry Smith 
128095fce210SBarry Smith    Note:
128195fce210SBarry Smith    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
128295fce210SBarry Smith    be done by calling PetscSFGetGraph().
128395fce210SBarry Smith 
128495fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph()
128595fce210SBarry Smith @*/
1286cd620004SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *esf)
128795fce210SBarry Smith {
1288cd620004SJunchao Zhang   PetscInt          i,j,n,nroots,nleaves,esf_nleaves,*new_ilocal,minleaf,maxleaf,maxlocal;
1289cd620004SJunchao Zhang   const PetscInt    *ilocal;
1290cd620004SJunchao Zhang   signed char       *rootdata,*leafdata,*leafmem;
1291ba2a7774SJunchao Zhang   const PetscSFNode *iremote;
1292f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1293f659e5c7SJunchao Zhang   MPI_Comm          comm;
12940511a646SMatthew G. Knepley   PetscErrorCode    ierr;
129595fce210SBarry Smith 
129695fce210SBarry Smith   PetscFunctionBegin;
129795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
129829046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1299ba2a7774SJunchao Zhang   if (nselected) PetscValidPointer(selected,3);
1300cd620004SJunchao Zhang   PetscValidPointer(esf,4);
13010511a646SMatthew G. Knepley 
1302f659e5c7SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1303140a1472SStefano Zampini   ierr = PetscLogEventBegin(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr);
1304f659e5c7SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1305cd620004SJunchao Zhang   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
1306cd620004SJunchao Zhang 
130776bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {  /* Error out if selected[] has dups or  out of range indices */
1308cd620004SJunchao Zhang     PetscBool dups;
1309cd620004SJunchao Zhang     ierr = PetscCheckDupsInt(nselected,selected,&dups);CHKERRQ(ierr);
1310cd620004SJunchao Zhang     if (dups) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"selected[] has dups");
1311cd620004SJunchao Zhang     for (i=0; i<nselected; i++)
1312cd620004SJunchao Zhang       if (selected[i] < 0 || selected[i] >= nroots) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"selected root indice %D is out of [0,%D)",selected[i],nroots);
1313cd620004SJunchao Zhang   }
1314f659e5c7SJunchao Zhang 
1315f659e5c7SJunchao Zhang   if (sf->ops->CreateEmbeddedSF) {
1316cd620004SJunchao Zhang     ierr = (*sf->ops->CreateEmbeddedSF)(sf,nselected,selected,esf);CHKERRQ(ierr);
1317f659e5c7SJunchao Zhang   } else {
1318cd620004SJunchao Zhang     /* A generic version of creating embedded sf */
1319cd620004SJunchao Zhang     ierr = PetscSFGetLeafRange(sf,&minleaf,&maxleaf);CHKERRQ(ierr);
1320cd620004SJunchao Zhang     maxlocal = maxleaf - minleaf + 1;
1321cd620004SJunchao Zhang     ierr = PetscCalloc2(nroots,&rootdata,maxlocal,&leafmem);CHKERRQ(ierr);
1322cd620004SJunchao Zhang     leafdata = leafmem - minleaf;
1323cd620004SJunchao Zhang     /* Tag selected roots and bcast to leaves */
1324cd620004SJunchao Zhang     for (i=0; i<nselected; i++) rootdata[selected[i]] = 1;
1325cd620004SJunchao Zhang     ierr = PetscSFBcastBegin(sf,MPI_SIGNED_CHAR,rootdata,leafdata);CHKERRQ(ierr);
1326cd620004SJunchao Zhang     ierr = PetscSFBcastEnd(sf,MPI_SIGNED_CHAR,rootdata,leafdata);CHKERRQ(ierr);
1327ba2a7774SJunchao Zhang 
1328cd620004SJunchao Zhang     /* Build esf with leaves that are still connected */
1329cd620004SJunchao Zhang     esf_nleaves = 0;
1330cd620004SJunchao Zhang     for (i=0; i<nleaves; i++) {
1331cd620004SJunchao Zhang       j = ilocal ? ilocal[i] : i;
1332cd620004SJunchao Zhang       /* esf_nleaves += leafdata[j] should work in theory, but failed with SFWindow bugs
1333cd620004SJunchao Zhang          with PetscSFBcast. See https://gitlab.com/petsc/petsc/issues/555
1334cd620004SJunchao Zhang       */
1335cd620004SJunchao Zhang       esf_nleaves += (leafdata[j] ? 1 : 0);
1336cd620004SJunchao Zhang     }
1337cd620004SJunchao Zhang     ierr = PetscMalloc1(esf_nleaves,&new_ilocal);CHKERRQ(ierr);
1338cd620004SJunchao Zhang     ierr = PetscMalloc1(esf_nleaves,&new_iremote);CHKERRQ(ierr);
1339cd620004SJunchao Zhang     for (i=n=0; i<nleaves; i++) {
1340cd620004SJunchao Zhang       j = ilocal ? ilocal[i] : i;
1341cd620004SJunchao Zhang       if (leafdata[j]) {
1342cd620004SJunchao Zhang         new_ilocal[n]        = j;
1343cd620004SJunchao Zhang         new_iremote[n].rank  = iremote[i].rank;
1344cd620004SJunchao Zhang         new_iremote[n].index = iremote[i].index;
1345fc1ede2bSMatthew G. Knepley         ++n;
134695fce210SBarry Smith       }
134795fce210SBarry Smith     }
1348cd620004SJunchao Zhang     ierr = PetscSFCreate(comm,esf);CHKERRQ(ierr);
1349cd620004SJunchao Zhang     ierr = PetscSFSetFromOptions(*esf);CHKERRQ(ierr);
1350cd620004SJunchao Zhang     ierr = PetscSFSetGraph(*esf,nroots,esf_nleaves,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1351cd620004SJunchao Zhang     ierr = PetscFree2(rootdata,leafmem);CHKERRQ(ierr);
1352f659e5c7SJunchao Zhang   }
1353140a1472SStefano Zampini   ierr = PetscLogEventEnd(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr);
135495fce210SBarry Smith   PetscFunctionReturn(0);
135595fce210SBarry Smith }
135695fce210SBarry Smith 
13572f5fb4c2SMatthew G. Knepley /*@C
13582f5fb4c2SMatthew G. Knepley   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
13592f5fb4c2SMatthew G. Knepley 
13602f5fb4c2SMatthew G. Knepley   Collective
13612f5fb4c2SMatthew G. Knepley 
13622f5fb4c2SMatthew G. Knepley   Input Arguments:
13632f5fb4c2SMatthew G. Knepley + sf - original star forest
1364f659e5c7SJunchao Zhang . nselected  - number of selected leaves on this process
1365f659e5c7SJunchao Zhang - selected   - indices of the selected leaves on this process
13662f5fb4c2SMatthew G. Knepley 
13672f5fb4c2SMatthew G. Knepley   Output Arguments:
13682f5fb4c2SMatthew G. Knepley .  newsf - new star forest
13692f5fb4c2SMatthew G. Knepley 
13702f5fb4c2SMatthew G. Knepley   Level: advanced
13712f5fb4c2SMatthew G. Knepley 
13722f5fb4c2SMatthew G. Knepley .seealso: PetscSFCreateEmbeddedSF(), PetscSFSetGraph(), PetscSFGetGraph()
13732f5fb4c2SMatthew G. Knepley @*/
1374f659e5c7SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
13752f5fb4c2SMatthew G. Knepley {
1376f659e5c7SJunchao Zhang   const PetscSFNode *iremote;
1377f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1378f659e5c7SJunchao Zhang   const PetscInt    *ilocal;
1379f659e5c7SJunchao Zhang   PetscInt          i,nroots,*leaves,*new_ilocal;
1380f659e5c7SJunchao Zhang   MPI_Comm          comm;
13812f5fb4c2SMatthew G. Knepley   PetscErrorCode    ierr;
13822f5fb4c2SMatthew G. Knepley 
13832f5fb4c2SMatthew G. Knepley   PetscFunctionBegin;
13842f5fb4c2SMatthew G. Knepley   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
138529046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1386f659e5c7SJunchao Zhang   if (nselected) PetscValidPointer(selected,3);
13872f5fb4c2SMatthew G. Knepley   PetscValidPointer(newsf,4);
13882f5fb4c2SMatthew G. Knepley 
1389f659e5c7SJunchao Zhang   /* Uniq selected[] and put results in leaves[] */
1390f659e5c7SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1391f659e5c7SJunchao Zhang   ierr = PetscMalloc1(nselected,&leaves);CHKERRQ(ierr);
1392dd5b3ca6SJunchao Zhang   ierr = PetscArraycpy(leaves,selected,nselected);CHKERRQ(ierr);
1393f659e5c7SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(&nselected,leaves);CHKERRQ(ierr);
1394f659e5c7SJunchao Zhang   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);
1395f659e5c7SJunchao Zhang 
1396f659e5c7SJunchao Zhang   /* Optimize the routine only when sf is setup and hence we can reuse sf's communication pattern */
1397f659e5c7SJunchao Zhang   if (sf->setupcalled && sf->ops->CreateEmbeddedLeafSF) {
1398f659e5c7SJunchao Zhang     ierr = (*sf->ops->CreateEmbeddedLeafSF)(sf,nselected,leaves,newsf);CHKERRQ(ierr);
1399f659e5c7SJunchao Zhang   } else {
1400f659e5c7SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,NULL,&ilocal,&iremote);CHKERRQ(ierr);
1401f659e5c7SJunchao Zhang     ierr = PetscMalloc1(nselected,&new_ilocal);CHKERRQ(ierr);
1402f659e5c7SJunchao Zhang     ierr = PetscMalloc1(nselected,&new_iremote);CHKERRQ(ierr);
1403f659e5c7SJunchao Zhang     for (i=0; i<nselected; ++i) {
1404f659e5c7SJunchao Zhang       const PetscInt l     = leaves[i];
1405f659e5c7SJunchao Zhang       new_ilocal[i]        = ilocal ? ilocal[l] : l;
1406f659e5c7SJunchao Zhang       new_iremote[i].rank  = iremote[l].rank;
1407f659e5c7SJunchao Zhang       new_iremote[i].index = iremote[l].index;
14082f5fb4c2SMatthew G. Knepley     }
14094cc19a0cSStefano Zampini     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,newsf);CHKERRQ(ierr);
1410f659e5c7SJunchao Zhang     ierr = PetscSFSetGraph(*newsf,nroots,nselected,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1411f659e5c7SJunchao Zhang   }
1412f659e5c7SJunchao Zhang   ierr = PetscFree(leaves);CHKERRQ(ierr);
14132f5fb4c2SMatthew G. Knepley   PetscFunctionReturn(0);
14142f5fb4c2SMatthew G. Knepley }
14152f5fb4c2SMatthew G. Knepley 
141695fce210SBarry Smith /*@C
14173482bfa8SJunchao Zhang    PetscSFBcastAndOpBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastAndOpEnd()
14183482bfa8SJunchao Zhang 
14193482bfa8SJunchao Zhang    Collective on PetscSF
14203482bfa8SJunchao Zhang 
14213482bfa8SJunchao Zhang    Input Arguments:
14223482bfa8SJunchao Zhang +  sf - star forest on which to communicate
14233482bfa8SJunchao Zhang .  unit - data type associated with each node
14243482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
14253482bfa8SJunchao Zhang -  op - operation to use for reduction
14263482bfa8SJunchao Zhang 
14273482bfa8SJunchao Zhang    Output Arguments:
14283482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
14293482bfa8SJunchao Zhang 
14303482bfa8SJunchao Zhang    Level: intermediate
14313482bfa8SJunchao Zhang 
1432d0295fc0SJunchao Zhang    Notes:
1433d0295fc0SJunchao Zhang     When petsc is configured with device support, it will use its own mechanism to figure out whether the given data pointers
1434d0295fc0SJunchao Zhang     are host pointers or device pointers, which may incur a noticable cost. If you already knew the info, you should
1435d0295fc0SJunchao Zhang     use PetscSFBcastAndOpWithMemTypeBegin() instead.
1436d0295fc0SJunchao Zhang .seealso: PetscSFBcastAndOpEnd(), PetscSFBcastBegin(), PetscSFBcastEnd(), PetscSFBcastAndOpWithMemTypeBegin()
14373482bfa8SJunchao Zhang @*/
14383482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
14393482bfa8SJunchao Zhang {
14403482bfa8SJunchao Zhang   PetscErrorCode ierr;
1441eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
14423482bfa8SJunchao Zhang 
14433482bfa8SJunchao Zhang   PetscFunctionBegin;
14443482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
14453482bfa8SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
144697929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);}
1447eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1448eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1449eb02082bSJunchao Zhang   ierr = (*sf->ops->BcastAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op);CHKERRQ(ierr);
145097929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);}
14513482bfa8SJunchao Zhang   PetscFunctionReturn(0);
14523482bfa8SJunchao Zhang }
14533482bfa8SJunchao Zhang 
14543482bfa8SJunchao Zhang /*@C
1455d0295fc0SJunchao Zhang    PetscSFBcastAndOpWithMemTypeBegin - begin pointwise broadcast with root value being reduced to leaf value with explicit memory types, to be concluded with call to PetscSFBcastAndOpEnd()
1456d0295fc0SJunchao Zhang 
1457d0295fc0SJunchao Zhang    Collective on PetscSF
1458d0295fc0SJunchao Zhang 
1459d0295fc0SJunchao Zhang    Input Arguments:
1460d0295fc0SJunchao Zhang +  sf - star forest on which to communicate
1461d0295fc0SJunchao Zhang .  unit - data type associated with each node
1462d0295fc0SJunchao Zhang .  rootmtype - memory type of rootdata
1463d0295fc0SJunchao Zhang .  rootdata - buffer to broadcast
1464d0295fc0SJunchao Zhang .  leafmtype - memory type of leafdata
1465d0295fc0SJunchao Zhang -  op - operation to use for reduction
1466d0295fc0SJunchao Zhang 
1467d0295fc0SJunchao Zhang    Output Arguments:
1468d0295fc0SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
1469d0295fc0SJunchao Zhang 
1470d0295fc0SJunchao Zhang    Level: intermediate
1471d0295fc0SJunchao Zhang 
1472d0295fc0SJunchao Zhang .seealso: PetscSFBcastAndOpEnd(), PetscSFBcastBegin(), PetscSFBcastEnd(),PetscSFBcastAndOpBegin()
1473d0295fc0SJunchao Zhang @*/
1474d0295fc0SJunchao Zhang PetscErrorCode PetscSFBcastAndOpWithMemTypeBegin(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,const void *rootdata,PetscMemType leafmtype,void *leafdata,MPI_Op op)
1475d0295fc0SJunchao Zhang {
1476d0295fc0SJunchao Zhang   PetscErrorCode ierr;
1477d0295fc0SJunchao Zhang 
1478d0295fc0SJunchao Zhang   PetscFunctionBegin;
1479d0295fc0SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1480d0295fc0SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
148197929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);}
1482d0295fc0SJunchao Zhang   ierr = (*sf->ops->BcastAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op);CHKERRQ(ierr);
148397929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);}
1484d0295fc0SJunchao Zhang   PetscFunctionReturn(0);
1485d0295fc0SJunchao Zhang }
1486d0295fc0SJunchao Zhang 
1487d0295fc0SJunchao Zhang /*@C
14883482bfa8SJunchao Zhang    PetscSFBcastAndOpEnd - end a broadcast & reduce operation started with PetscSFBcastAndOpBegin()
14893482bfa8SJunchao Zhang 
14903482bfa8SJunchao Zhang    Collective
14913482bfa8SJunchao Zhang 
14923482bfa8SJunchao Zhang    Input Arguments:
14933482bfa8SJunchao Zhang +  sf - star forest
14943482bfa8SJunchao Zhang .  unit - data type
14953482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
14963482bfa8SJunchao Zhang -  op - operation to use for reduction
14973482bfa8SJunchao Zhang 
14983482bfa8SJunchao Zhang    Output Arguments:
14993482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
15003482bfa8SJunchao Zhang 
15013482bfa8SJunchao Zhang    Level: intermediate
15023482bfa8SJunchao Zhang 
15033482bfa8SJunchao Zhang .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
15043482bfa8SJunchao Zhang @*/
15053482bfa8SJunchao Zhang PetscErrorCode PetscSFBcastAndOpEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
15063482bfa8SJunchao Zhang {
15073482bfa8SJunchao Zhang   PetscErrorCode ierr;
15083482bfa8SJunchao Zhang 
15093482bfa8SJunchao Zhang   PetscFunctionBegin;
15103482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
151197929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventBegin(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);}
151200816365SJunchao Zhang   ierr = (*sf->ops->BcastAndOpEnd)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr);
151397929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventEnd(PETSCSF_BcastAndOpEnd,sf,0,0,0);CHKERRQ(ierr);}
15143482bfa8SJunchao Zhang   PetscFunctionReturn(0);
15153482bfa8SJunchao Zhang }
15163482bfa8SJunchao Zhang 
15173482bfa8SJunchao Zhang /*@C
151895fce210SBarry Smith    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
151995fce210SBarry Smith 
152095fce210SBarry Smith    Collective
152195fce210SBarry Smith 
152295fce210SBarry Smith    Input Arguments:
152395fce210SBarry Smith +  sf - star forest
152495fce210SBarry Smith .  unit - data type
152595fce210SBarry Smith .  leafdata - values to reduce
152695fce210SBarry Smith -  op - reduction operation
152795fce210SBarry Smith 
152895fce210SBarry Smith    Output Arguments:
152995fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
153095fce210SBarry Smith 
153195fce210SBarry Smith    Level: intermediate
153295fce210SBarry Smith 
1533d0295fc0SJunchao Zhang    Notes:
1534d0295fc0SJunchao Zhang     When petsc is configured with device support, it will use its own mechanism to figure out whether the given data pointers
1535d0295fc0SJunchao Zhang     are host pointers or device pointers, which may incur a noticable cost. If you already knew the info, you should
1536d0295fc0SJunchao Zhang     use PetscSFReduceWithMemTypeBegin() instead.
1537d0295fc0SJunchao Zhang 
1538d0295fc0SJunchao Zhang .seealso: PetscSFBcastBegin(), PetscSFReduceWithMemTypeBegin()
153995fce210SBarry Smith @*/
154095fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
154195fce210SBarry Smith {
154295fce210SBarry Smith   PetscErrorCode ierr;
1543eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
154495fce210SBarry Smith 
154595fce210SBarry Smith   PetscFunctionBegin;
154695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
154795fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
154897929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);}
1549eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1550eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1551eb02082bSJunchao Zhang   ierr = (sf->ops->ReduceBegin)(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op);CHKERRQ(ierr);
155297929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);}
155395fce210SBarry Smith   PetscFunctionReturn(0);
155495fce210SBarry Smith }
155595fce210SBarry Smith 
155695fce210SBarry Smith /*@C
1557d0295fc0SJunchao Zhang    PetscSFReduceWithMemTypeBegin - begin reduction of leafdata into rootdata with explicit memory types, to be completed with call to PetscSFReduceEnd()
1558d0295fc0SJunchao Zhang 
1559d0295fc0SJunchao Zhang    Collective
1560d0295fc0SJunchao Zhang 
1561d0295fc0SJunchao Zhang    Input Arguments:
1562d0295fc0SJunchao Zhang +  sf - star forest
1563d0295fc0SJunchao Zhang .  unit - data type
1564d0295fc0SJunchao Zhang .  leafmtype - memory type of leafdata
1565d0295fc0SJunchao Zhang .  leafdata - values to reduce
1566d0295fc0SJunchao Zhang .  rootmtype - memory type of rootdata
1567d0295fc0SJunchao Zhang -  op - reduction operation
1568d0295fc0SJunchao Zhang 
1569d0295fc0SJunchao Zhang    Output Arguments:
1570d0295fc0SJunchao Zhang .  rootdata - result of reduction of values from all leaves of each root
1571d0295fc0SJunchao Zhang 
1572d0295fc0SJunchao Zhang    Level: intermediate
1573d0295fc0SJunchao Zhang 
1574d0295fc0SJunchao Zhang .seealso: PetscSFBcastBegin(), PetscSFReduceBegin()
1575d0295fc0SJunchao Zhang @*/
1576d0295fc0SJunchao Zhang PetscErrorCode PetscSFReduceWithMemTypeBegin(PetscSF sf,MPI_Datatype unit,PetscMemType leafmtype,const void *leafdata,PetscMemType rootmtype,void *rootdata,MPI_Op op)
1577d0295fc0SJunchao Zhang {
1578d0295fc0SJunchao Zhang   PetscErrorCode ierr;
1579d0295fc0SJunchao Zhang 
1580d0295fc0SJunchao Zhang   PetscFunctionBegin;
1581d0295fc0SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1582d0295fc0SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
158397929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);}
1584d0295fc0SJunchao Zhang   ierr = (sf->ops->ReduceBegin)(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op);CHKERRQ(ierr);
158597929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);}
1586d0295fc0SJunchao Zhang   PetscFunctionReturn(0);
1587d0295fc0SJunchao Zhang }
1588d0295fc0SJunchao Zhang 
1589d0295fc0SJunchao Zhang /*@C
159095fce210SBarry Smith    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
159195fce210SBarry Smith 
159295fce210SBarry Smith    Collective
159395fce210SBarry Smith 
159495fce210SBarry Smith    Input Arguments:
159595fce210SBarry Smith +  sf - star forest
159695fce210SBarry Smith .  unit - data type
159795fce210SBarry Smith .  leafdata - values to reduce
159895fce210SBarry Smith -  op - reduction operation
159995fce210SBarry Smith 
160095fce210SBarry Smith    Output Arguments:
160195fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
160295fce210SBarry Smith 
160395fce210SBarry Smith    Level: intermediate
160495fce210SBarry Smith 
160595fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
160695fce210SBarry Smith @*/
160795fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
160895fce210SBarry Smith {
160995fce210SBarry Smith   PetscErrorCode ierr;
161095fce210SBarry Smith 
161195fce210SBarry Smith   PetscFunctionBegin;
161295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
161397929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);}
161400816365SJunchao Zhang   ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
161597929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);}
161695fce210SBarry Smith   PetscFunctionReturn(0);
161795fce210SBarry Smith }
161895fce210SBarry Smith 
161995fce210SBarry Smith /*@C
1620a1729e3fSJunchao Zhang    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
1621a1729e3fSJunchao Zhang 
1622a1729e3fSJunchao Zhang    Collective
1623a1729e3fSJunchao Zhang 
1624a1729e3fSJunchao Zhang    Input Arguments:
1625a1729e3fSJunchao Zhang +  sf - star forest
1626a1729e3fSJunchao Zhang .  unit - data type
1627a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1628a1729e3fSJunchao Zhang -  op - operation to use for reduction
1629a1729e3fSJunchao Zhang 
1630a1729e3fSJunchao Zhang    Output Arguments:
1631a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1632a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1633a1729e3fSJunchao Zhang 
1634a1729e3fSJunchao Zhang    Level: advanced
1635a1729e3fSJunchao Zhang 
1636a1729e3fSJunchao Zhang    Note:
1637a1729e3fSJunchao Zhang    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
1638a1729e3fSJunchao Zhang    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
1639a1729e3fSJunchao Zhang    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
1640a1729e3fSJunchao Zhang    integers.
1641a1729e3fSJunchao Zhang 
1642a1729e3fSJunchao Zhang .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
1643a1729e3fSJunchao Zhang @*/
1644a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1645a1729e3fSJunchao Zhang {
1646a1729e3fSJunchao Zhang   PetscErrorCode ierr;
1647eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype,leafupdatemtype;
1648a1729e3fSJunchao Zhang 
1649a1729e3fSJunchao Zhang   PetscFunctionBegin;
1650a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1651a1729e3fSJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1652a1729e3fSJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1653eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1654eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1655eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafupdate,&leafupdatemtype);CHKERRQ(ierr);
1656eb02082bSJunchao Zhang   if (leafmtype != leafupdatemtype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for leafdata and leafupdate in different memory types");
1657eb02082bSJunchao Zhang   ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,leafupdate,op);CHKERRQ(ierr);
1658a1729e3fSJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1659a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1660a1729e3fSJunchao Zhang }
1661a1729e3fSJunchao Zhang 
1662a1729e3fSJunchao Zhang /*@C
1663a1729e3fSJunchao Zhang    PetscSFFetchAndOpEnd - end operation started in matching call to PetscSFFetchAndOpBegin() to fetch values from roots and update atomically by applying operation using my leaf value
1664a1729e3fSJunchao Zhang 
1665a1729e3fSJunchao Zhang    Collective
1666a1729e3fSJunchao Zhang 
1667a1729e3fSJunchao Zhang    Input Arguments:
1668a1729e3fSJunchao Zhang +  sf - star forest
1669a1729e3fSJunchao Zhang .  unit - data type
1670a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1671a1729e3fSJunchao Zhang -  op - operation to use for reduction
1672a1729e3fSJunchao Zhang 
1673a1729e3fSJunchao Zhang    Output Arguments:
1674a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1675a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1676a1729e3fSJunchao Zhang 
1677a1729e3fSJunchao Zhang    Level: advanced
1678a1729e3fSJunchao Zhang 
1679a1729e3fSJunchao Zhang .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
1680a1729e3fSJunchao Zhang @*/
1681a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1682a1729e3fSJunchao Zhang {
1683a1729e3fSJunchao Zhang   PetscErrorCode ierr;
1684a1729e3fSJunchao Zhang 
1685a1729e3fSJunchao Zhang   PetscFunctionBegin;
1686a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1687a1729e3fSJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
168800816365SJunchao Zhang   ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1689a1729e3fSJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1690a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1691a1729e3fSJunchao Zhang }
1692a1729e3fSJunchao Zhang 
1693a1729e3fSJunchao Zhang /*@C
169495fce210SBarry Smith    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
169595fce210SBarry Smith 
169695fce210SBarry Smith    Collective
169795fce210SBarry Smith 
169895fce210SBarry Smith    Input Arguments:
169995fce210SBarry Smith .  sf - star forest
170095fce210SBarry Smith 
170195fce210SBarry Smith    Output Arguments:
170295fce210SBarry Smith .  degree - degree of each root vertex
170395fce210SBarry Smith 
170495fce210SBarry Smith    Level: advanced
170595fce210SBarry Smith 
1706ffe67aa5SVáclav Hapla    Notes:
1707ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1708ffe67aa5SVáclav Hapla 
170995fce210SBarry Smith .seealso: PetscSFGatherBegin()
171095fce210SBarry Smith @*/
171195fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
171295fce210SBarry Smith {
171395fce210SBarry Smith   PetscErrorCode ierr;
171495fce210SBarry Smith 
171595fce210SBarry Smith   PetscFunctionBegin;
171695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
171795fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
171895fce210SBarry Smith   PetscValidPointer(degree,2);
1719803bd9e8SMatthew G. Knepley   if (!sf->degreeknown) {
17205b0d146aSStefano Zampini     PetscInt i, nroots = sf->nroots, maxlocal;
1721803bd9e8SMatthew G. Knepley     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
17225b0d146aSStefano Zampini     maxlocal = sf->maxleaf-sf->minleaf+1;
172329046d53SLisandro Dalcin     ierr = PetscMalloc1(nroots,&sf->degree);CHKERRQ(ierr);
172429046d53SLisandro Dalcin     ierr = PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp);CHKERRQ(ierr); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
172529046d53SLisandro Dalcin     for (i=0; i<nroots; i++) sf->degree[i] = 0;
17269837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
17275b0d146aSStefano Zampini     ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp-sf->minleaf,sf->degree,MPI_SUM);CHKERRQ(ierr);
172895fce210SBarry Smith   }
172995fce210SBarry Smith   *degree = NULL;
173095fce210SBarry Smith   PetscFunctionReturn(0);
173195fce210SBarry Smith }
173295fce210SBarry Smith 
173395fce210SBarry Smith /*@C
173495fce210SBarry Smith    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
173595fce210SBarry Smith 
173695fce210SBarry Smith    Collective
173795fce210SBarry Smith 
173895fce210SBarry Smith    Input Arguments:
173995fce210SBarry Smith .  sf - star forest
174095fce210SBarry Smith 
174195fce210SBarry Smith    Output Arguments:
174295fce210SBarry Smith .  degree - degree of each root vertex
174395fce210SBarry Smith 
174495fce210SBarry Smith    Level: developer
174595fce210SBarry Smith 
1746ffe67aa5SVáclav Hapla    Notes:
1747ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1748ffe67aa5SVáclav Hapla 
174995fce210SBarry Smith .seealso:
175095fce210SBarry Smith @*/
175195fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
175295fce210SBarry Smith {
175395fce210SBarry Smith   PetscErrorCode ierr;
175495fce210SBarry Smith 
175595fce210SBarry Smith   PetscFunctionBegin;
175695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
175795fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
175829046d53SLisandro Dalcin   PetscValidPointer(degree,2);
175995fce210SBarry Smith   if (!sf->degreeknown) {
176029046d53SLisandro Dalcin     if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
17615b0d146aSStefano Zampini     ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp-sf->minleaf,sf->degree,MPI_SUM);CHKERRQ(ierr);
176295fce210SBarry Smith     ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr);
176395fce210SBarry Smith     sf->degreeknown = PETSC_TRUE;
176495fce210SBarry Smith   }
176595fce210SBarry Smith   *degree = sf->degree;
176695fce210SBarry Smith   PetscFunctionReturn(0);
176795fce210SBarry Smith }
176895fce210SBarry Smith 
1769673100f5SVaclav Hapla 
1770673100f5SVaclav Hapla /*@C
177166dfcd1aSVaclav Hapla    PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()).
177266dfcd1aSVaclav Hapla    Each multi-root is assigned index of the corresponding original root.
1773673100f5SVaclav Hapla 
1774673100f5SVaclav Hapla    Collective
1775673100f5SVaclav Hapla 
1776673100f5SVaclav Hapla    Input Arguments:
1777673100f5SVaclav Hapla +  sf - star forest
1778673100f5SVaclav Hapla -  degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd()
1779673100f5SVaclav Hapla 
1780673100f5SVaclav Hapla    Output Arguments:
178166dfcd1aSVaclav Hapla +  nMultiRoots - (optional) number of multi-roots (roots of multi-SF)
178266dfcd1aSVaclav Hapla -  multiRootsOrigNumbering - original indices of multi-roots; length of this array is nMultiRoots
1783673100f5SVaclav Hapla 
1784673100f5SVaclav Hapla    Level: developer
1785673100f5SVaclav Hapla 
1786ffe67aa5SVáclav Hapla    Notes:
1787ffe67aa5SVáclav Hapla    The returned array multiRootsOrigNumbering is newly allocated and should be destroyed with PetscFree() when no longer needed.
1788ffe67aa5SVáclav Hapla 
1789673100f5SVaclav Hapla .seealso: PetscSFComputeDegreeBegin(), PetscSFComputeDegreeEnd(), PetscSFGetMultiSF()
1790673100f5SVaclav Hapla @*/
179166dfcd1aSVaclav Hapla PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[])
1792673100f5SVaclav Hapla {
1793673100f5SVaclav Hapla   PetscSF             msf;
1794673100f5SVaclav Hapla   PetscInt            i, j, k, nroots, nmroots;
1795673100f5SVaclav Hapla   PetscErrorCode      ierr;
1796673100f5SVaclav Hapla 
1797673100f5SVaclav Hapla   PetscFunctionBegin;
1798673100f5SVaclav Hapla   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1799673100f5SVaclav Hapla   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
180029328920SVaclav Hapla   if (nroots) PetscValidIntPointer(degree,2);
180166dfcd1aSVaclav Hapla   if (nMultiRoots) PetscValidIntPointer(nMultiRoots,3);
180266dfcd1aSVaclav Hapla   PetscValidPointer(multiRootsOrigNumbering,4);
1803673100f5SVaclav Hapla   ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr);
1804673100f5SVaclav Hapla   ierr = PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL);CHKERRQ(ierr);
180566dfcd1aSVaclav Hapla   ierr = PetscMalloc1(nmroots, multiRootsOrigNumbering);CHKERRQ(ierr);
1806673100f5SVaclav Hapla   for (i=0,j=0,k=0; i<nroots; i++) {
1807673100f5SVaclav Hapla     if (!degree[i]) continue;
1808673100f5SVaclav Hapla     for (j=0; j<degree[i]; j++,k++) {
180966dfcd1aSVaclav Hapla       (*multiRootsOrigNumbering)[k] = i;
1810673100f5SVaclav Hapla     }
1811673100f5SVaclav Hapla   }
1812673100f5SVaclav Hapla   if (PetscUnlikely(k != nmroots)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail");
181366dfcd1aSVaclav Hapla   if (nMultiRoots) *nMultiRoots = nmroots;
1814673100f5SVaclav Hapla   PetscFunctionReturn(0);
1815673100f5SVaclav Hapla }
1816673100f5SVaclav Hapla 
181795fce210SBarry Smith /*@C
181895fce210SBarry Smith    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
181995fce210SBarry Smith 
182095fce210SBarry Smith    Collective
182195fce210SBarry Smith 
182295fce210SBarry Smith    Input Arguments:
182395fce210SBarry Smith +  sf - star forest
182495fce210SBarry Smith .  unit - data type
182595fce210SBarry Smith -  leafdata - leaf data to gather to roots
182695fce210SBarry Smith 
182795fce210SBarry Smith    Output Argument:
182895fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
182995fce210SBarry Smith 
183095fce210SBarry Smith    Level: intermediate
183195fce210SBarry Smith 
183295fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
183395fce210SBarry Smith @*/
183495fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
183595fce210SBarry Smith {
183695fce210SBarry Smith   PetscErrorCode ierr;
1837a5526d50SJunchao Zhang   PetscSF        multi = NULL;
183895fce210SBarry Smith 
183995fce210SBarry Smith   PetscFunctionBegin;
184095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
184129046d53SLisandro Dalcin   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
184295fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
18438bfbc91cSJed Brown   ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
184495fce210SBarry Smith   PetscFunctionReturn(0);
184595fce210SBarry Smith }
184695fce210SBarry Smith 
184795fce210SBarry Smith /*@C
184895fce210SBarry Smith    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
184995fce210SBarry Smith 
185095fce210SBarry Smith    Collective
185195fce210SBarry Smith 
185295fce210SBarry Smith    Input Arguments:
185395fce210SBarry Smith +  sf - star forest
185495fce210SBarry Smith .  unit - data type
185595fce210SBarry Smith -  leafdata - leaf data to gather to roots
185695fce210SBarry Smith 
185795fce210SBarry Smith    Output Argument:
185895fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
185995fce210SBarry Smith 
186095fce210SBarry Smith    Level: intermediate
186195fce210SBarry Smith 
186295fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
186395fce210SBarry Smith @*/
186495fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
186595fce210SBarry Smith {
186695fce210SBarry Smith   PetscErrorCode ierr;
1867a5526d50SJunchao Zhang   PetscSF        multi = NULL;
186895fce210SBarry Smith 
186995fce210SBarry Smith   PetscFunctionBegin;
187095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
187195fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
18728bfbc91cSJed Brown   ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPIU_REPLACE);CHKERRQ(ierr);
187395fce210SBarry Smith   PetscFunctionReturn(0);
187495fce210SBarry Smith }
187595fce210SBarry Smith 
187695fce210SBarry Smith /*@C
187795fce210SBarry Smith    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
187895fce210SBarry Smith 
187995fce210SBarry Smith    Collective
188095fce210SBarry Smith 
188195fce210SBarry Smith    Input Arguments:
188295fce210SBarry Smith +  sf - star forest
188395fce210SBarry Smith .  unit - data type
188495fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
188595fce210SBarry Smith 
188695fce210SBarry Smith    Output Argument:
188795fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
188895fce210SBarry Smith 
188995fce210SBarry Smith    Level: intermediate
189095fce210SBarry Smith 
189195fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
189295fce210SBarry Smith @*/
189395fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
189495fce210SBarry Smith {
189595fce210SBarry Smith   PetscErrorCode ierr;
1896a5526d50SJunchao Zhang   PetscSF        multi = NULL;
189795fce210SBarry Smith 
189895fce210SBarry Smith   PetscFunctionBegin;
189995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
190095fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
190195fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
190295fce210SBarry Smith   ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
190395fce210SBarry Smith   PetscFunctionReturn(0);
190495fce210SBarry Smith }
190595fce210SBarry Smith 
190695fce210SBarry Smith /*@C
190795fce210SBarry Smith    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
190895fce210SBarry Smith 
190995fce210SBarry Smith    Collective
191095fce210SBarry Smith 
191195fce210SBarry Smith    Input Arguments:
191295fce210SBarry Smith +  sf - star forest
191395fce210SBarry Smith .  unit - data type
191495fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
191595fce210SBarry Smith 
191695fce210SBarry Smith    Output Argument:
191795fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
191895fce210SBarry Smith 
191995fce210SBarry Smith    Level: intermediate
192095fce210SBarry Smith 
192195fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
192295fce210SBarry Smith @*/
192395fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
192495fce210SBarry Smith {
192595fce210SBarry Smith   PetscErrorCode ierr;
1926a5526d50SJunchao Zhang   PetscSF        multi = NULL;
192795fce210SBarry Smith 
192895fce210SBarry Smith   PetscFunctionBegin;
192995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
193095fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
193195fce210SBarry Smith   ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata);CHKERRQ(ierr);
193295fce210SBarry Smith   PetscFunctionReturn(0);
193395fce210SBarry Smith }
1934a7b3aa13SAta Mesgarnejad 
1935a072220fSLawrence Mitchell static PetscErrorCode PetscSFCheckLeavesUnique_Private(PetscSF sf)
1936a072220fSLawrence Mitchell {
1937a072220fSLawrence Mitchell   PetscInt        i, n, nleaves;
1938a072220fSLawrence Mitchell   const PetscInt *ilocal = NULL;
1939a072220fSLawrence Mitchell   PetscHSetI      seen;
1940a072220fSLawrence Mitchell   PetscErrorCode  ierr;
1941a072220fSLawrence Mitchell 
1942a072220fSLawrence Mitchell   PetscFunctionBegin;
1943b458e8f1SJose E. Roman   if (PetscDefined(USE_DEBUG)) {
1944a072220fSLawrence Mitchell     ierr = PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
1945a072220fSLawrence Mitchell     ierr = PetscHSetICreate(&seen);CHKERRQ(ierr);
1946a072220fSLawrence Mitchell     for (i = 0; i < nleaves; i++) {
1947a072220fSLawrence Mitchell       const PetscInt leaf = ilocal ? ilocal[i] : i;
1948a072220fSLawrence Mitchell       ierr = PetscHSetIAdd(seen,leaf);CHKERRQ(ierr);
1949a072220fSLawrence Mitchell     }
1950a072220fSLawrence Mitchell     ierr = PetscHSetIGetSize(seen,&n);CHKERRQ(ierr);
1951a072220fSLawrence Mitchell     if (n != nleaves) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided leaves have repeated values: all leaves must be unique");
1952a072220fSLawrence Mitchell     ierr = PetscHSetIDestroy(&seen);CHKERRQ(ierr);
1953b458e8f1SJose E. Roman   }
1954a072220fSLawrence Mitchell   PetscFunctionReturn(0);
1955a072220fSLawrence Mitchell }
195654729392SStefano Zampini 
1957a7b3aa13SAta Mesgarnejad /*@
195804c0ada0SJunchao Zhang   PetscSFCompose - Compose a new PetscSF by putting the second SF under the first one in a top (roots) down (leaves) view
1959a7b3aa13SAta Mesgarnejad 
1960a7b3aa13SAta Mesgarnejad   Input Parameters:
196154729392SStefano Zampini + sfA - The first PetscSF
196254729392SStefano Zampini - sfB - The second PetscSF
1963a7b3aa13SAta Mesgarnejad 
1964a7b3aa13SAta Mesgarnejad   Output Parameters:
196554729392SStefano Zampini . sfBA - The composite SF
1966a7b3aa13SAta Mesgarnejad 
1967a7b3aa13SAta Mesgarnejad   Level: developer
1968a7b3aa13SAta Mesgarnejad 
1969a072220fSLawrence Mitchell   Notes:
197054729392SStefano Zampini   Currently, the two SFs must be defined on congruent communicators and they must be true star
197154729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots.
197254729392SStefano Zampini 
197354729392SStefano Zampini   sfA's leaf space and sfB's root space might be partially overlapped. The composition builds
197454729392SStefano Zampini   a graph with sfA's roots and sfB's leaves only when there is a path between them. Unconnected
197554729392SStefano Zampini   nodes (roots or leaves) are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a
197654729392SStefano Zampini   Bcast on sfA, then a Bcast on sfB, on connected nodes.
1977a072220fSLawrence Mitchell 
197804c0ada0SJunchao Zhang .seealso: PetscSF, PetscSFComposeInverse(), PetscSFGetGraph(), PetscSFSetGraph()
1979a7b3aa13SAta Mesgarnejad @*/
1980a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
1981a7b3aa13SAta Mesgarnejad {
198204c0ada0SJunchao Zhang   PetscErrorCode    ierr;
1983a7b3aa13SAta Mesgarnejad   const PetscSFNode *remotePointsA,*remotePointsB;
1984d41018fbSJunchao Zhang   PetscSFNode       *remotePointsBA=NULL,*reorderedRemotePointsA = NULL,*leafdataB;
198554729392SStefano Zampini   const PetscInt    *localPointsA,*localPointsB;
198654729392SStefano Zampini   PetscInt          *localPointsBA;
198754729392SStefano Zampini   PetscInt          i,numRootsA,numLeavesA,numRootsB,numLeavesB,minleaf,maxleaf,numLeavesBA;
198854729392SStefano Zampini   PetscBool         denseB;
1989a7b3aa13SAta Mesgarnejad 
1990a7b3aa13SAta Mesgarnejad   PetscFunctionBegin;
1991a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1);
199229046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfA,1);
199329046d53SLisandro Dalcin   PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2);
199429046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfB,2);
199554729392SStefano Zampini   PetscCheckSameComm(sfA,1,sfB,2);
199629046d53SLisandro Dalcin   PetscValidPointer(sfBA,3);
199754729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfA);CHKERRQ(ierr);
199854729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfB);CHKERRQ(ierr);
199954729392SStefano Zampini 
2000a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfA,&numRootsA,&numLeavesA,&localPointsA,&remotePointsA);CHKERRQ(ierr);
2001a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfB,&numRootsB,&numLeavesB,&localPointsB,&remotePointsB);CHKERRQ(ierr);
2002d41018fbSJunchao Zhang   if (localPointsA) {
200354729392SStefano Zampini     ierr = PetscMalloc1(numRootsB,&reorderedRemotePointsA);CHKERRQ(ierr);
200454729392SStefano Zampini     for (i=0; i<numRootsB; i++) {
200554729392SStefano Zampini       reorderedRemotePointsA[i].rank = -1;
200654729392SStefano Zampini       reorderedRemotePointsA[i].index = -1;
200754729392SStefano Zampini     }
200854729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
200954729392SStefano Zampini       if (localPointsA[i] >= numRootsB) continue;
201054729392SStefano Zampini       reorderedRemotePointsA[localPointsA[i]] = remotePointsA[i];
201154729392SStefano Zampini     }
2012d41018fbSJunchao Zhang     remotePointsA = reorderedRemotePointsA;
2013d41018fbSJunchao Zhang   }
2014d41018fbSJunchao Zhang   ierr = PetscSFGetLeafRange(sfB,&minleaf,&maxleaf);CHKERRQ(ierr);
2015d41018fbSJunchao Zhang   ierr = PetscMalloc1(maxleaf-minleaf+1,&leafdataB);CHKERRQ(ierr);
2016d41018fbSJunchao Zhang   ierr = PetscSFBcastBegin(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf);CHKERRQ(ierr);
2017d41018fbSJunchao Zhang   ierr = PetscSFBcastEnd(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf);CHKERRQ(ierr);
2018d41018fbSJunchao Zhang   ierr = PetscFree(reorderedRemotePointsA);CHKERRQ(ierr);
2019d41018fbSJunchao Zhang 
202054729392SStefano Zampini   denseB = (PetscBool)!localPointsB;
202154729392SStefano Zampini   for (i=0,numLeavesBA=0; i<numLeavesB; i++) {
202254729392SStefano Zampini     if (leafdataB[localPointsB ? localPointsB[i]-minleaf : i].rank == -1) denseB = PETSC_FALSE;
202354729392SStefano Zampini     else numLeavesBA++;
202454729392SStefano Zampini   }
202554729392SStefano Zampini   if (denseB) {
2026d41018fbSJunchao Zhang     localPointsBA  = NULL;
2027d41018fbSJunchao Zhang     remotePointsBA = leafdataB;
2028d41018fbSJunchao Zhang   } else {
202954729392SStefano Zampini     ierr = PetscMalloc1(numLeavesBA,&localPointsBA);CHKERRQ(ierr);
203054729392SStefano Zampini     ierr = PetscMalloc1(numLeavesBA,&remotePointsBA);CHKERRQ(ierr);
203154729392SStefano Zampini     for (i=0,numLeavesBA=0; i<numLeavesB; i++) {
203254729392SStefano Zampini       const PetscInt l = localPointsB ? localPointsB[i] : i;
203354729392SStefano Zampini 
203454729392SStefano Zampini       if (leafdataB[l-minleaf].rank == -1) continue;
203554729392SStefano Zampini       remotePointsBA[numLeavesBA] = leafdataB[l-minleaf];
203654729392SStefano Zampini       localPointsBA[numLeavesBA] = l;
203754729392SStefano Zampini       numLeavesBA++;
203854729392SStefano Zampini     }
2039d41018fbSJunchao Zhang     ierr = PetscFree(leafdataB);CHKERRQ(ierr);
2040d41018fbSJunchao Zhang   }
204154729392SStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA);CHKERRQ(ierr);
204220c24465SJunchao Zhang   ierr = PetscSFSetFromOptions(*sfBA);CHKERRQ(ierr);
204354729392SStefano Zampini   ierr = PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER);CHKERRQ(ierr);
2044a7b3aa13SAta Mesgarnejad   PetscFunctionReturn(0);
2045a7b3aa13SAta Mesgarnejad }
20461c6ba672SJunchao Zhang 
204704c0ada0SJunchao Zhang /*@
204854729392SStefano Zampini   PetscSFComposeInverse - Compose a new PetscSF by putting the inverse of the second SF under the first one
204904c0ada0SJunchao Zhang 
205004c0ada0SJunchao Zhang   Input Parameters:
205154729392SStefano Zampini + sfA - The first PetscSF
205254729392SStefano Zampini - sfB - The second PetscSF
205304c0ada0SJunchao Zhang 
205404c0ada0SJunchao Zhang   Output Parameters:
205554729392SStefano Zampini . sfBA - The composite SF.
205604c0ada0SJunchao Zhang 
205704c0ada0SJunchao Zhang   Level: developer
205804c0ada0SJunchao Zhang 
205954729392SStefano Zampini   Notes:
206054729392SStefano Zampini   Currently, the two SFs must be defined on congruent communicators and they must be true star
206154729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots. Even more, all roots of the
206254729392SStefano Zampini   second SF must have a degree of 1, i.e., no roots have more than one leaf connected.
206354729392SStefano Zampini 
206454729392SStefano Zampini   sfA's leaf space and sfB's leaf space might be partially overlapped. The composition builds
206554729392SStefano Zampini   a graph with sfA's roots and sfB's roots only when there is a path between them. Unconnected
206654729392SStefano Zampini   roots are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a Bcast on sfA, then
206754729392SStefano Zampini   a Reduce on sfB, on connected roots.
206854729392SStefano Zampini 
206954729392SStefano Zampini .seealso: PetscSF, PetscSFCompose(), PetscSFGetGraph(), PetscSFSetGraph(), PetscSFCreateInverseSF()
207004c0ada0SJunchao Zhang @*/
207104c0ada0SJunchao Zhang PetscErrorCode PetscSFComposeInverse(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
207204c0ada0SJunchao Zhang {
207304c0ada0SJunchao Zhang   PetscErrorCode    ierr;
207404c0ada0SJunchao Zhang   const PetscSFNode *remotePointsA,*remotePointsB;
207504c0ada0SJunchao Zhang   PetscSFNode       *remotePointsBA;
207604c0ada0SJunchao Zhang   const PetscInt    *localPointsA,*localPointsB;
207754729392SStefano Zampini   PetscSFNode       *reorderedRemotePointsA = NULL;
207854729392SStefano Zampini   PetscInt          i,numRootsA,numLeavesA,numLeavesBA,numRootsB,numLeavesB,minleaf,maxleaf,*localPointsBA;
20795b0d146aSStefano Zampini   MPI_Op            op;
20805b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
20815b0d146aSStefano Zampini   PetscBool         iswin;
20825b0d146aSStefano Zampini #endif
208304c0ada0SJunchao Zhang 
208404c0ada0SJunchao Zhang   PetscFunctionBegin;
208504c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1);
208604c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfA,1);
208704c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2);
208804c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfB,2);
208954729392SStefano Zampini   PetscCheckSameComm(sfA,1,sfB,2);
209004c0ada0SJunchao Zhang   PetscValidPointer(sfBA,3);
209154729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfA);CHKERRQ(ierr);
209254729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfB);CHKERRQ(ierr);
209354729392SStefano Zampini 
209404c0ada0SJunchao Zhang   ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr);
209504c0ada0SJunchao Zhang   ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr);
20965b0d146aSStefano Zampini 
20975b0d146aSStefano Zampini   /* TODO: Check roots of sfB have degree of 1 */
20985b0d146aSStefano Zampini   /* Once we implement it, we can replace the MPI_MAXLOC
20995b0d146aSStefano Zampini      with MPIU_REPLACE. In that case, MPI_MAXLOC and MPIU_REPLACE have the same effect.
21005b0d146aSStefano Zampini      We use MPI_MAXLOC only to have a deterministic output from this routine if
21015b0d146aSStefano Zampini      the root condition is not meet.
21025b0d146aSStefano Zampini    */
21035b0d146aSStefano Zampini   op = MPI_MAXLOC;
21045b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
21055b0d146aSStefano Zampini   /* we accept a non-deterministic output (if any) with PETSCSFWINDOW, since MPI_MAXLOC cannot operate on MPIU_2INT with MPI_Accumulate */
21065b0d146aSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)sfB,PETSCSFWINDOW,&iswin);CHKERRQ(ierr);
21075b0d146aSStefano Zampini   if (iswin) op = MPIU_REPLACE;
21085b0d146aSStefano Zampini #endif
21095b0d146aSStefano Zampini 
211054729392SStefano Zampini   ierr = PetscSFGetLeafRange(sfB, &minleaf, &maxleaf);CHKERRQ(ierr);
211154729392SStefano Zampini   ierr = PetscMalloc1(maxleaf - minleaf + 1,&reorderedRemotePointsA);CHKERRQ(ierr);
211254729392SStefano Zampini   for (i=0; i<maxleaf - minleaf + 1; i++) {
211354729392SStefano Zampini     reorderedRemotePointsA[i].rank = -1;
211454729392SStefano Zampini     reorderedRemotePointsA[i].index = -1;
211554729392SStefano Zampini   }
211654729392SStefano Zampini   if (localPointsA) {
211754729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
211854729392SStefano Zampini       if (localPointsA[i] > maxleaf || localPointsA[i] < minleaf) continue;
211954729392SStefano Zampini       reorderedRemotePointsA[localPointsA[i] - minleaf] = remotePointsA[i];
212054729392SStefano Zampini     }
212154729392SStefano Zampini   } else {
212254729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
212354729392SStefano Zampini       if (i > maxleaf || i < minleaf) continue;
212454729392SStefano Zampini       reorderedRemotePointsA[i - minleaf] = remotePointsA[i];
212554729392SStefano Zampini     }
212654729392SStefano Zampini   }
212754729392SStefano Zampini 
212854729392SStefano Zampini   ierr = PetscMalloc1(numRootsB,&localPointsBA);CHKERRQ(ierr);
212904c0ada0SJunchao Zhang   ierr = PetscMalloc1(numRootsB,&remotePointsBA);CHKERRQ(ierr);
213054729392SStefano Zampini   for (i=0; i<numRootsB; i++) {
213154729392SStefano Zampini     remotePointsBA[i].rank = -1;
213254729392SStefano Zampini     remotePointsBA[i].index = -1;
213354729392SStefano Zampini   }
213454729392SStefano Zampini 
21355b0d146aSStefano Zampini   ierr = PetscSFReduceBegin(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,op);CHKERRQ(ierr);
21365b0d146aSStefano Zampini   ierr = PetscSFReduceEnd(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,op);CHKERRQ(ierr);
213754729392SStefano Zampini   ierr = PetscFree(reorderedRemotePointsA);CHKERRQ(ierr);
213854729392SStefano Zampini   for (i=0,numLeavesBA=0; i<numRootsB; i++) {
213954729392SStefano Zampini     if (remotePointsBA[i].rank == -1) continue;
214054729392SStefano Zampini     remotePointsBA[numLeavesBA].rank = remotePointsBA[i].rank;
214154729392SStefano Zampini     remotePointsBA[numLeavesBA].index = remotePointsBA[i].index;
214254729392SStefano Zampini     localPointsBA[numLeavesBA] = i;
214354729392SStefano Zampini     numLeavesBA++;
214454729392SStefano Zampini   }
214554729392SStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA);CHKERRQ(ierr);
214620c24465SJunchao Zhang   ierr = PetscSFSetFromOptions(*sfBA);CHKERRQ(ierr);
214754729392SStefano Zampini   ierr = PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER);CHKERRQ(ierr);
214804c0ada0SJunchao Zhang   PetscFunctionReturn(0);
214904c0ada0SJunchao Zhang }
215004c0ada0SJunchao Zhang 
21511c6ba672SJunchao Zhang /*
21521c6ba672SJunchao Zhang   PetscSFCreateLocalSF_Private - Creates a local PetscSF that only has intra-process edges of the global PetscSF
21531c6ba672SJunchao Zhang 
21541c6ba672SJunchao Zhang   Input Parameters:
21551c6ba672SJunchao Zhang . sf - The global PetscSF
21561c6ba672SJunchao Zhang 
21571c6ba672SJunchao Zhang   Output Parameters:
21581c6ba672SJunchao Zhang . out - The local PetscSF
21591c6ba672SJunchao Zhang  */
21601c6ba672SJunchao Zhang PetscErrorCode PetscSFCreateLocalSF_Private(PetscSF sf,PetscSF *out)
21611c6ba672SJunchao Zhang {
21621c6ba672SJunchao Zhang   MPI_Comm           comm;
21631c6ba672SJunchao Zhang   PetscMPIInt        myrank;
21641c6ba672SJunchao Zhang   const PetscInt     *ilocal;
21651c6ba672SJunchao Zhang   const PetscSFNode  *iremote;
21661c6ba672SJunchao Zhang   PetscInt           i,j,nroots,nleaves,lnleaves,*lilocal;
21671c6ba672SJunchao Zhang   PetscSFNode        *liremote;
21681c6ba672SJunchao Zhang   PetscSF            lsf;
21691c6ba672SJunchao Zhang   PetscErrorCode     ierr;
21701c6ba672SJunchao Zhang 
21711c6ba672SJunchao Zhang   PetscFunctionBegin;
21721c6ba672SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
21731c6ba672SJunchao Zhang   if (sf->ops->CreateLocalSF) {
21741c6ba672SJunchao Zhang     ierr = (*sf->ops->CreateLocalSF)(sf,out);CHKERRQ(ierr);
21751c6ba672SJunchao Zhang   } else {
21761c6ba672SJunchao Zhang     /* Could use PetscSFCreateEmbeddedLeafSF, but since we know the comm is PETSC_COMM_SELF, we can make it fast */
21771c6ba672SJunchao Zhang     ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
21781c6ba672SJunchao Zhang     ierr = MPI_Comm_rank(comm,&myrank);CHKERRQ(ierr);
21791c6ba672SJunchao Zhang 
21801c6ba672SJunchao Zhang     /* Find out local edges and build a local SF */
21811c6ba672SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
21821c6ba672SJunchao Zhang     for (i=lnleaves=0; i<nleaves; i++) {if (iremote[i].rank == (PetscInt)myrank) lnleaves++;}
21831c6ba672SJunchao Zhang     ierr = PetscMalloc1(lnleaves,&lilocal);CHKERRQ(ierr);
21841c6ba672SJunchao Zhang     ierr = PetscMalloc1(lnleaves,&liremote);CHKERRQ(ierr);
21851c6ba672SJunchao Zhang 
21861c6ba672SJunchao Zhang     for (i=j=0; i<nleaves; i++) {
21871c6ba672SJunchao Zhang       if (iremote[i].rank == (PetscInt)myrank) {
21881c6ba672SJunchao Zhang         lilocal[j]        = ilocal? ilocal[i] : i; /* ilocal=NULL for contiguous storage */
21891c6ba672SJunchao Zhang         liremote[j].rank  = 0; /* rank in PETSC_COMM_SELF */
21901c6ba672SJunchao Zhang         liremote[j].index = iremote[i].index;
21911c6ba672SJunchao Zhang         j++;
21921c6ba672SJunchao Zhang       }
21931c6ba672SJunchao Zhang     }
21941c6ba672SJunchao Zhang     ierr = PetscSFCreate(PETSC_COMM_SELF,&lsf);CHKERRQ(ierr);
219520c24465SJunchao Zhang     ierr = PetscSFSetFromOptions(lsf);CHKERRQ(ierr);
21961c6ba672SJunchao Zhang     ierr = PetscSFSetGraph(lsf,nroots,lnleaves,lilocal,PETSC_OWN_POINTER,liremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
21971c6ba672SJunchao Zhang     ierr = PetscSFSetUp(lsf);CHKERRQ(ierr);
21981c6ba672SJunchao Zhang     *out = lsf;
21991c6ba672SJunchao Zhang   }
22001c6ba672SJunchao Zhang   PetscFunctionReturn(0);
22011c6ba672SJunchao Zhang }
2202dd5b3ca6SJunchao Zhang 
2203dd5b3ca6SJunchao Zhang /* Similar to PetscSFBcast, but only Bcast to leaves on rank 0 */
2204dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFBcastToZero_Private(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
2205dd5b3ca6SJunchao Zhang {
2206dd5b3ca6SJunchao Zhang   PetscErrorCode     ierr;
2207eb02082bSJunchao Zhang   PetscMemType       rootmtype,leafmtype;
2208dd5b3ca6SJunchao Zhang 
2209dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
2210dd5b3ca6SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
2211dd5b3ca6SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
2212dd5b3ca6SJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
2213eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
2214eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
2215dd5b3ca6SJunchao Zhang   if (sf->ops->BcastToZero) {
2216eb02082bSJunchao Zhang     ierr = (*sf->ops->BcastToZero)(sf,unit,rootmtype,rootdata,leafmtype,leafdata);CHKERRQ(ierr);
2217dd5b3ca6SJunchao Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFBcastToZero_Private is not supported on this SF type");
2218dd5b3ca6SJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
2219dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
2220dd5b3ca6SJunchao Zhang }
2221dd5b3ca6SJunchao Zhang 
2222