xref: /petsc/src/vec/is/sf/interface/sf.c (revision 546078ac55d0f4e36dc638bf7d35f3433941479a)
1af0996ceSBarry Smith #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/
2c4e6a40aSLawrence Mitchell #include <petsc/private/hashseti.h>
353dd6d7dSJunchao Zhang #include <petsc/private/viewerimpl.h>
495fce210SBarry Smith #include <petscctable.h>
595fce210SBarry Smith 
67fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
77fd2d3dbSJunchao Zhang   #include <cuda_runtime.h>
87fd2d3dbSJunchao Zhang #endif
97fd2d3dbSJunchao Zhang 
107fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_HIP)
117fd2d3dbSJunchao Zhang   #include <hip/hip_runtime.h>
127fd2d3dbSJunchao Zhang #endif
137fd2d3dbSJunchao Zhang 
142abc8c78SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
152abc8c78SJacob Faibussowitsch void PetscSFCheckGraphSet(PetscSF,int);
162abc8c78SJacob Faibussowitsch #else
1795fce210SBarry Smith #if defined(PETSC_USE_DEBUG)
1895fce210SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) do {                          \
1995fce210SBarry Smith     if (PetscUnlikely(!(sf)->graphset))                              \
202abc8c78SJacob Faibussowitsch       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() or PetscSFSetGraphWithPattern() on argument %d \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME); \
2195fce210SBarry Smith   } while (0)
2295fce210SBarry Smith #else
2395fce210SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) do {} while (0)
2495fce210SBarry Smith #endif
252abc8c78SJacob Faibussowitsch #endif
2695fce210SBarry Smith 
274c8fdceaSLisandro Dalcin const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",NULL};
2895fce210SBarry Smith 
297fd2d3dbSJunchao Zhang PETSC_STATIC_INLINE PetscErrorCode PetscGetMemType(const void *data,PetscMemType *type)
307fd2d3dbSJunchao Zhang {
317fd2d3dbSJunchao Zhang   PetscFunctionBegin;
327fd2d3dbSJunchao Zhang   PetscValidPointer(type,2);
337fd2d3dbSJunchao Zhang   *type = PETSC_MEMTYPE_HOST;
347fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
35a4af0ceeSJacob Faibussowitsch   if (PetscDeviceInitialized(PETSC_DEVICE_CUDA) && data) {
367fd2d3dbSJunchao Zhang     cudaError_t                  cerr;
377fd2d3dbSJunchao Zhang     struct cudaPointerAttributes attr;
387fd2d3dbSJunchao Zhang     enum cudaMemoryType          mtype;
397fd2d3dbSJunchao Zhang     cerr = cudaPointerGetAttributes(&attr,data); /* Do not check error since before CUDA 11.0, passing a host pointer returns cudaErrorInvalidValue */
407fd2d3dbSJunchao Zhang     cudaGetLastError(); /* Reset the last error */
417fd2d3dbSJunchao Zhang     #if (CUDART_VERSION < 10000)
427fd2d3dbSJunchao Zhang       mtype = attr.memoryType;
437fd2d3dbSJunchao Zhang     #else
447fd2d3dbSJunchao Zhang       mtype = attr.type;
457fd2d3dbSJunchao Zhang     #endif
467fd2d3dbSJunchao Zhang     if (cerr == cudaSuccess && mtype == cudaMemoryTypeDevice) *type = PETSC_MEMTYPE_DEVICE;
477fd2d3dbSJunchao Zhang   }
487fd2d3dbSJunchao Zhang #endif
497fd2d3dbSJunchao Zhang 
507fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_HIP)
51a4af0ceeSJacob Faibussowitsch   if (PetscDeviceInitialized(PETSC_DEVICE_HIP) && data) {
527fd2d3dbSJunchao Zhang     hipError_t                   cerr;
537fd2d3dbSJunchao Zhang     struct hipPointerAttribute_t attr;
547fd2d3dbSJunchao Zhang     enum hipMemoryType           mtype;
5559af0bd3SScott Kruger     cerr = hipPointerGetAttributes(&attr,data);
567fd2d3dbSJunchao Zhang     hipGetLastError(); /* Reset the last error */
577fd2d3dbSJunchao Zhang     mtype = attr.memoryType;
587fd2d3dbSJunchao Zhang     if (cerr == hipSuccess && mtype == hipMemoryTypeDevice) *type = PETSC_MEMTYPE_DEVICE;
597fd2d3dbSJunchao Zhang   }
607fd2d3dbSJunchao Zhang #endif
617fd2d3dbSJunchao Zhang   PetscFunctionReturn(0);
627fd2d3dbSJunchao Zhang }
637fd2d3dbSJunchao Zhang 
648af6ec1cSBarry Smith /*@
6595fce210SBarry Smith    PetscSFCreate - create a star forest communication context
6695fce210SBarry Smith 
67d083f849SBarry Smith    Collective
6895fce210SBarry Smith 
694165533cSJose E. Roman    Input Parameter:
7095fce210SBarry Smith .  comm - communicator on which the star forest will operate
7195fce210SBarry Smith 
724165533cSJose E. Roman    Output Parameter:
7395fce210SBarry Smith .  sf - new star forest context
7495fce210SBarry Smith 
75dd5b3ca6SJunchao Zhang    Options Database Keys:
76dd5b3ca6SJunchao Zhang +  -sf_type basic     -Use MPI persistent Isend/Irecv for communication (Default)
77dd5b3ca6SJunchao Zhang .  -sf_type window    -Use MPI-3 one-sided window for communication
78dd5b3ca6SJunchao Zhang -  -sf_type neighbor  -Use MPI-3 neighborhood collectives for communication
79dd5b3ca6SJunchao Zhang 
8095fce210SBarry Smith    Level: intermediate
8195fce210SBarry Smith 
82dd5b3ca6SJunchao Zhang    Notes:
83dd5b3ca6SJunchao Zhang    When one knows the communication graph is one of the predefined graph, such as MPI_Alltoall, MPI_Allgatherv,
84dd5b3ca6SJunchao Zhang    MPI_Gatherv, one can create a PetscSF and then set its graph with PetscSFSetGraphWithPattern(). These special
85dd5b3ca6SJunchao Zhang    SFs are optimized and they have better performance than general SFs.
86dd5b3ca6SJunchao Zhang 
87dd5b3ca6SJunchao Zhang .seealso: PetscSFSetGraph(), PetscSFSetGraphWithPattern(), PetscSFDestroy()
8895fce210SBarry Smith @*/
8995fce210SBarry Smith PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf)
9095fce210SBarry Smith {
9195fce210SBarry Smith   PetscErrorCode ierr;
9295fce210SBarry Smith   PetscSF        b;
9395fce210SBarry Smith 
9495fce210SBarry Smith   PetscFunctionBegin;
9595fce210SBarry Smith   PetscValidPointer(sf,2);
96607a6623SBarry Smith   ierr = PetscSFInitializePackage();CHKERRQ(ierr);
9795fce210SBarry Smith 
9873107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(b,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView);CHKERRQ(ierr);
9995fce210SBarry Smith 
10095fce210SBarry Smith   b->nroots    = -1;
10195fce210SBarry Smith   b->nleaves   = -1;
10229046d53SLisandro Dalcin   b->minleaf   = PETSC_MAX_INT;
10329046d53SLisandro Dalcin   b->maxleaf   = PETSC_MIN_INT;
10495fce210SBarry Smith   b->nranks    = -1;
10595fce210SBarry Smith   b->rankorder = PETSC_TRUE;
10695fce210SBarry Smith   b->ingroup   = MPI_GROUP_NULL;
10795fce210SBarry Smith   b->outgroup  = MPI_GROUP_NULL;
10895fce210SBarry Smith   b->graphset  = PETSC_FALSE;
10920c24465SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
11020c24465SJunchao Zhang   b->use_gpu_aware_mpi    = use_gpu_aware_mpi;
11120c24465SJunchao Zhang   b->use_stream_aware_mpi = PETSC_FALSE;
11271438e86SJunchao Zhang   b->unknown_input_stream= PETSC_FALSE;
11327f636e8SJunchao Zhang   #if defined(PETSC_HAVE_KOKKOS) /* Prefer kokkos over cuda*/
11420c24465SJunchao Zhang     b->backend = PETSCSF_BACKEND_KOKKOS;
11527f636e8SJunchao Zhang   #elif defined(PETSC_HAVE_CUDA)
11627f636e8SJunchao Zhang     b->backend = PETSCSF_BACKEND_CUDA;
11759af0bd3SScott Kruger   #elif defined(PETSC_HAVE_HIP)
11859af0bd3SScott Kruger     b->backend = PETSCSF_BACKEND_HIP;
11920c24465SJunchao Zhang   #endif
12071438e86SJunchao Zhang 
12171438e86SJunchao Zhang   #if defined(PETSC_HAVE_NVSHMEM)
12271438e86SJunchao Zhang     b->use_nvshmem     = PETSC_FALSE; /* Default is not to try NVSHMEM */
12371438e86SJunchao Zhang     b->use_nvshmem_get = PETSC_FALSE; /* Default is to use nvshmem_put based protocol */
12471438e86SJunchao Zhang     ierr = PetscOptionsGetBool(NULL,NULL,"-use_nvshmem",&b->use_nvshmem,NULL);CHKERRQ(ierr);
12571438e86SJunchao Zhang     ierr = PetscOptionsGetBool(NULL,NULL,"-use_nvshmem_get",&b->use_nvshmem_get,NULL);CHKERRQ(ierr);
12671438e86SJunchao Zhang   #endif
12720c24465SJunchao Zhang #endif
12860c22052SBarry Smith   b->vscat.from_n = -1;
12960c22052SBarry Smith   b->vscat.to_n   = -1;
13060c22052SBarry Smith   b->vscat.unit   = MPIU_SCALAR;
13195fce210SBarry Smith  *sf = b;
13295fce210SBarry Smith   PetscFunctionReturn(0);
13395fce210SBarry Smith }
13495fce210SBarry Smith 
13529046d53SLisandro Dalcin /*@
13695fce210SBarry Smith    PetscSFReset - Reset a star forest so that different sizes or neighbors can be used
13795fce210SBarry Smith 
13895fce210SBarry Smith    Collective
13995fce210SBarry Smith 
1404165533cSJose E. Roman    Input Parameter:
14195fce210SBarry Smith .  sf - star forest
14295fce210SBarry Smith 
14395fce210SBarry Smith    Level: advanced
14495fce210SBarry Smith 
14595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph(), PetscSFDestroy()
14695fce210SBarry Smith @*/
14795fce210SBarry Smith PetscErrorCode PetscSFReset(PetscSF sf)
14895fce210SBarry Smith {
14995fce210SBarry Smith   PetscErrorCode ierr;
15095fce210SBarry Smith 
15195fce210SBarry Smith   PetscFunctionBegin;
15295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
15379715d56SJed Brown   if (sf->ops->Reset) {ierr = (*sf->ops->Reset)(sf);CHKERRQ(ierr);}
15429046d53SLisandro Dalcin   sf->nroots   = -1;
15529046d53SLisandro Dalcin   sf->nleaves  = -1;
15629046d53SLisandro Dalcin   sf->minleaf  = PETSC_MAX_INT;
15729046d53SLisandro Dalcin   sf->maxleaf  = PETSC_MIN_INT;
15895fce210SBarry Smith   sf->mine     = NULL;
15995fce210SBarry Smith   sf->remote   = NULL;
16029046d53SLisandro Dalcin   sf->graphset = PETSC_FALSE;
16129046d53SLisandro Dalcin   ierr = PetscFree(sf->mine_alloc);CHKERRQ(ierr);
16295fce210SBarry Smith   ierr = PetscFree(sf->remote_alloc);CHKERRQ(ierr);
16321c688dcSJed Brown   sf->nranks = -1;
16429046d53SLisandro Dalcin   ierr = PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote);CHKERRQ(ierr);
16529046d53SLisandro Dalcin   sf->degreeknown = PETSC_FALSE;
16695fce210SBarry Smith   ierr = PetscFree(sf->degree);CHKERRQ(ierr);
167ffc4695bSBarry Smith   if (sf->ingroup  != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->ingroup);CHKERRMPI(ierr);}
168ffc4695bSBarry Smith   if (sf->outgroup != MPI_GROUP_NULL) {ierr = MPI_Group_free(&sf->outgroup);CHKERRMPI(ierr);}
169013b3241SStefano Zampini   if (sf->multi) sf->multi->multi = NULL;
17095fce210SBarry Smith   ierr = PetscSFDestroy(&sf->multi);CHKERRQ(ierr);
171dd5b3ca6SJunchao Zhang   ierr = PetscLayoutDestroy(&sf->map);CHKERRQ(ierr);
17271438e86SJunchao Zhang 
17371438e86SJunchao Zhang  #if defined(PETSC_HAVE_DEVICE)
17471438e86SJunchao Zhang   for (PetscInt i=0; i<2; i++) {ierr = PetscSFFree(sf,PETSC_MEMTYPE_DEVICE,sf->rmine_d[i]);CHKERRQ(ierr);}
17571438e86SJunchao Zhang  #endif
17671438e86SJunchao Zhang 
17795fce210SBarry Smith   sf->setupcalled = PETSC_FALSE;
17895fce210SBarry Smith   PetscFunctionReturn(0);
17995fce210SBarry Smith }
18095fce210SBarry Smith 
18195fce210SBarry Smith /*@C
18229046d53SLisandro Dalcin    PetscSFSetType - Set the PetscSF communication implementation
18395fce210SBarry Smith 
18495fce210SBarry Smith    Collective on PetscSF
18595fce210SBarry Smith 
18695fce210SBarry Smith    Input Parameters:
18795fce210SBarry Smith +  sf - the PetscSF context
18895fce210SBarry Smith -  type - a known method
18995fce210SBarry Smith 
19095fce210SBarry Smith    Options Database Key:
19195fce210SBarry Smith .  -sf_type <type> - Sets the method; use -help for a list
19270616304SStefano Zampini    of available methods (for instance, window, basic, neighbor)
19395fce210SBarry Smith 
19495fce210SBarry Smith    Notes:
19595fce210SBarry Smith    See "include/petscsf.h" for available methods (for instance)
19695fce210SBarry Smith +    PETSCSFWINDOW - MPI-2/3 one-sided
19795fce210SBarry Smith -    PETSCSFBASIC - basic implementation using MPI-1 two-sided
19895fce210SBarry Smith 
19995fce210SBarry Smith   Level: intermediate
20095fce210SBarry Smith 
20195fce210SBarry Smith .seealso: PetscSFType, PetscSFCreate()
20295fce210SBarry Smith @*/
20395fce210SBarry Smith PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type)
20495fce210SBarry Smith {
20595fce210SBarry Smith   PetscErrorCode ierr,(*r)(PetscSF);
20695fce210SBarry Smith   PetscBool      match;
20795fce210SBarry Smith 
20895fce210SBarry Smith   PetscFunctionBegin;
20995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
21095fce210SBarry Smith   PetscValidCharPointer(type,2);
21195fce210SBarry Smith 
21295fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)sf,type,&match);CHKERRQ(ierr);
21395fce210SBarry Smith   if (match) PetscFunctionReturn(0);
21495fce210SBarry Smith 
215adc40e5bSBarry Smith   ierr = PetscFunctionListFind(PetscSFList,type,&r);CHKERRQ(ierr);
21695fce210SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type);
21729046d53SLisandro Dalcin   /* Destroy the previous PetscSF implementation context */
21829046d53SLisandro Dalcin   if (sf->ops->Destroy) {ierr = (*(sf)->ops->Destroy)(sf);CHKERRQ(ierr);}
21995fce210SBarry Smith   ierr = PetscMemzero(sf->ops,sizeof(*sf->ops));CHKERRQ(ierr);
22095fce210SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)sf,type);CHKERRQ(ierr);
22195fce210SBarry Smith   ierr = (*r)(sf);CHKERRQ(ierr);
22295fce210SBarry Smith   PetscFunctionReturn(0);
22395fce210SBarry Smith }
22495fce210SBarry Smith 
22529046d53SLisandro Dalcin /*@C
22629046d53SLisandro Dalcin   PetscSFGetType - Get the PetscSF communication implementation
22729046d53SLisandro Dalcin 
22829046d53SLisandro Dalcin   Not Collective
22929046d53SLisandro Dalcin 
23029046d53SLisandro Dalcin   Input Parameter:
23129046d53SLisandro Dalcin . sf  - the PetscSF context
23229046d53SLisandro Dalcin 
23329046d53SLisandro Dalcin   Output Parameter:
23429046d53SLisandro Dalcin . type - the PetscSF type name
23529046d53SLisandro Dalcin 
23629046d53SLisandro Dalcin   Level: intermediate
23729046d53SLisandro Dalcin 
23829046d53SLisandro Dalcin .seealso: PetscSFSetType(), PetscSFCreate()
23929046d53SLisandro Dalcin @*/
24029046d53SLisandro Dalcin PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type)
24129046d53SLisandro Dalcin {
24229046d53SLisandro Dalcin   PetscFunctionBegin;
24329046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID,1);
24429046d53SLisandro Dalcin   PetscValidPointer(type,2);
24529046d53SLisandro Dalcin   *type = ((PetscObject)sf)->type_name;
24629046d53SLisandro Dalcin   PetscFunctionReturn(0);
24729046d53SLisandro Dalcin }
24829046d53SLisandro Dalcin 
2491fb7b255SJunchao Zhang /*@C
25095fce210SBarry Smith    PetscSFDestroy - destroy star forest
25195fce210SBarry Smith 
25295fce210SBarry Smith    Collective
25395fce210SBarry Smith 
2544165533cSJose E. Roman    Input Parameter:
25595fce210SBarry Smith .  sf - address of star forest
25695fce210SBarry Smith 
25795fce210SBarry Smith    Level: intermediate
25895fce210SBarry Smith 
25995fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFReset()
26095fce210SBarry Smith @*/
26195fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf)
26295fce210SBarry Smith {
26395fce210SBarry Smith   PetscErrorCode ierr;
26495fce210SBarry Smith 
26595fce210SBarry Smith   PetscFunctionBegin;
26695fce210SBarry Smith   if (!*sf) PetscFunctionReturn(0);
26795fce210SBarry Smith   PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1);
26829046d53SLisandro Dalcin   if (--((PetscObject)(*sf))->refct > 0) {*sf = NULL; PetscFunctionReturn(0);}
26995fce210SBarry Smith   ierr = PetscSFReset(*sf);CHKERRQ(ierr);
27095fce210SBarry Smith   if ((*sf)->ops->Destroy) {ierr = (*(*sf)->ops->Destroy)(*sf);CHKERRQ(ierr);}
27197929ea7SJunchao Zhang   ierr = PetscSFDestroy(&(*sf)->vscat.lsf);CHKERRQ(ierr);
27255b25c41SPierre Jolivet   if ((*sf)->vscat.bs > 1) {ierr = MPI_Type_free(&(*sf)->vscat.unit);CHKERRMPI(ierr);}
27395fce210SBarry Smith   ierr = PetscHeaderDestroy(sf);CHKERRQ(ierr);
27495fce210SBarry Smith   PetscFunctionReturn(0);
27595fce210SBarry Smith }
27695fce210SBarry Smith 
277c4e6a40aSLawrence Mitchell static PetscErrorCode PetscSFCheckGraphValid_Private(PetscSF sf)
278c4e6a40aSLawrence Mitchell {
279c4e6a40aSLawrence Mitchell   PetscInt           i, nleaves;
280c4e6a40aSLawrence Mitchell   PetscMPIInt        size;
281c4e6a40aSLawrence Mitchell   const PetscInt    *ilocal;
282c4e6a40aSLawrence Mitchell   const PetscSFNode *iremote;
283c4e6a40aSLawrence Mitchell   PetscErrorCode     ierr;
284c4e6a40aSLawrence Mitchell 
285c4e6a40aSLawrence Mitchell   PetscFunctionBegin;
28676bd3646SJed Brown   if (!sf->graphset || !PetscDefined(USE_DEBUG)) PetscFunctionReturn(0);
287c4e6a40aSLawrence Mitchell   ierr = PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
288ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRMPI(ierr);
289c4e6a40aSLawrence Mitchell   for (i = 0; i < nleaves; i++) {
290c4e6a40aSLawrence Mitchell     const PetscInt rank = iremote[i].rank;
291c4e6a40aSLawrence Mitchell     const PetscInt remote = iremote[i].index;
292c4e6a40aSLawrence Mitchell     const PetscInt leaf = ilocal ? ilocal[i] : i;
2932abc8c78SJacob Faibussowitsch     if (rank < 0 || rank >= size) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided rank (%" PetscInt_FMT ") for remote %" PetscInt_FMT " is invalid, should be in [0, %d)",rank,i,size);
2942abc8c78SJacob Faibussowitsch     if (remote < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided index (%" PetscInt_FMT ") for remote %" PetscInt_FMT " is invalid, should be >= 0",remote,i);
2952abc8c78SJacob Faibussowitsch     if (leaf < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided location (%" PetscInt_FMT ") for leaf %" PetscInt_FMT " is invalid, should be >= 0",leaf,i);
296c4e6a40aSLawrence Mitchell   }
297c4e6a40aSLawrence Mitchell   PetscFunctionReturn(0);
298c4e6a40aSLawrence Mitchell }
299c4e6a40aSLawrence Mitchell 
30095fce210SBarry Smith /*@
30195fce210SBarry Smith    PetscSFSetUp - set up communication structures
30295fce210SBarry Smith 
30395fce210SBarry Smith    Collective
30495fce210SBarry Smith 
3054165533cSJose E. Roman    Input Parameter:
30695fce210SBarry Smith .  sf - star forest communication object
30795fce210SBarry Smith 
30895fce210SBarry Smith    Level: beginner
30995fce210SBarry Smith 
31095fce210SBarry Smith .seealso: PetscSFSetFromOptions(), PetscSFSetType()
31195fce210SBarry Smith @*/
31295fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf)
31395fce210SBarry Smith {
31495fce210SBarry Smith   PetscErrorCode ierr;
31595fce210SBarry Smith 
31695fce210SBarry Smith   PetscFunctionBegin;
31729046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
31829046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
31995fce210SBarry Smith   if (sf->setupcalled) PetscFunctionReturn(0);
32029046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr);
32120c24465SJunchao Zhang   ierr = PetscSFCheckGraphValid_Private(sf);CHKERRQ(ierr);
32220c24465SJunchao Zhang   if (!((PetscObject)sf)->type_name) {ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);} /* Zero all sf->ops */
32395fce210SBarry Smith   if (sf->ops->SetUp) {ierr = (*sf->ops->SetUp)(sf);CHKERRQ(ierr);}
32420c24465SJunchao Zhang #if defined(PETSC_HAVE_CUDA)
32520c24465SJunchao Zhang   if (sf->backend == PETSCSF_BACKEND_CUDA) {
32671438e86SJunchao Zhang     sf->ops->Malloc = PetscSFMalloc_CUDA;
32771438e86SJunchao Zhang     sf->ops->Free   = PetscSFFree_CUDA;
32820c24465SJunchao Zhang   }
32920c24465SJunchao Zhang #endif
33059af0bd3SScott Kruger #if defined(PETSC_HAVE_HIP)
33159af0bd3SScott Kruger   if (sf->backend == PETSCSF_BACKEND_HIP) {
33259af0bd3SScott Kruger     sf->ops->Malloc = PetscSFMalloc_HIP;
33359af0bd3SScott Kruger     sf->ops->Free   = PetscSFFree_HIP;
33459af0bd3SScott Kruger   }
33559af0bd3SScott Kruger #endif
33620c24465SJunchao Zhang 
33759af0bd3SScott Kruger #
33820c24465SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
33920c24465SJunchao Zhang   if (sf->backend == PETSCSF_BACKEND_KOKKOS) {
34020c24465SJunchao Zhang     sf->ops->Malloc = PetscSFMalloc_Kokkos;
34120c24465SJunchao Zhang     sf->ops->Free   = PetscSFFree_Kokkos;
34220c24465SJunchao Zhang   }
34320c24465SJunchao Zhang #endif
34429046d53SLisandro Dalcin   ierr = PetscLogEventEnd(PETSCSF_SetUp,sf,0,0,0);CHKERRQ(ierr);
34595fce210SBarry Smith   sf->setupcalled = PETSC_TRUE;
34695fce210SBarry Smith   PetscFunctionReturn(0);
34795fce210SBarry Smith }
34895fce210SBarry Smith 
3498af6ec1cSBarry Smith /*@
35095fce210SBarry Smith    PetscSFSetFromOptions - set PetscSF options using the options database
35195fce210SBarry Smith 
35295fce210SBarry Smith    Logically Collective
35395fce210SBarry Smith 
3544165533cSJose E. Roman    Input Parameter:
35595fce210SBarry Smith .  sf - star forest
35695fce210SBarry Smith 
35795fce210SBarry Smith    Options Database Keys:
35860263706SJed Brown +  -sf_type               - implementation type, see PetscSFSetType()
35951ccb202SJunchao Zhang .  -sf_rank_order         - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
360b85e67b7SJunchao Zhang .  -sf_use_default_stream - Assume callers of SF computed the input root/leafdata with the default cuda stream. SF will also
361c2a741eeSJunchao Zhang                             use the default stream to process data. Therefore, no stream synchronization is needed between SF and its caller (default: true).
362c06a8e02SRichard Tran Mills                             If true, this option only works with -use_gpu_aware_mpi 1.
36320c24465SJunchao 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).
364c06a8e02SRichard Tran Mills                                If true, this option only works with -use_gpu_aware_mpi 1.
36595fce210SBarry Smith 
36659af0bd3SScott Kruger -  -sf_backend cuda | hip | kokkos -Select the device backend SF uses. Currently SF has these backends: cuda, hip and Kokkos.
36759af0bd3SScott Kruger                               On CUDA (HIP) devices, one can choose cuda (hip) or kokkos with the default being kokkos. On other devices,
36820c24465SJunchao Zhang                               the only available is kokkos.
36920c24465SJunchao Zhang 
37095fce210SBarry Smith    Level: intermediate
37195fce210SBarry Smith @*/
37295fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
37395fce210SBarry Smith {
37495fce210SBarry Smith   PetscSFType    deft;
37595fce210SBarry Smith   char           type[256];
37695fce210SBarry Smith   PetscErrorCode ierr;
37795fce210SBarry Smith   PetscBool      flg;
37895fce210SBarry Smith 
37995fce210SBarry Smith   PetscFunctionBegin;
38095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
38195fce210SBarry Smith   ierr = PetscObjectOptionsBegin((PetscObject)sf);CHKERRQ(ierr);
38295fce210SBarry Smith   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
38329046d53SLisandro Dalcin   ierr = PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,sizeof(type),&flg);CHKERRQ(ierr);
38495fce210SBarry Smith   ierr = PetscSFSetType(sf,flg ? type : deft);CHKERRQ(ierr);
38595fce210SBarry 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);
3867fd2d3dbSJunchao Zhang  #if defined(PETSC_HAVE_DEVICE)
38720c24465SJunchao Zhang   {
38820c24465SJunchao Zhang     char        backendstr[32] = {0};
38959af0bd3SScott Kruger     PetscBool   isCuda = PETSC_FALSE,isHip = PETSC_FALSE,isKokkos = PETSC_FALSE,set;
39020c24465SJunchao Zhang     /* Change the defaults set in PetscSFCreate() with command line options */
39171438e86SJunchao Zhang     ierr = PetscOptionsBool("-sf_unknown_input_stream","SF root/leafdata is computed on arbitary streams unknown to SF","PetscSFSetFromOptions",sf->unknown_input_stream,&sf->unknown_input_stream,NULL);CHKERRQ(ierr);
392b85e67b7SJunchao 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);
39320c24465SJunchao Zhang     ierr = PetscOptionsString("-sf_backend","Select the device backend SF uses","PetscSFSetFromOptions",NULL,backendstr,sizeof(backendstr),&set);CHKERRQ(ierr);
39420c24465SJunchao Zhang     ierr = PetscStrcasecmp("cuda",backendstr,&isCuda);CHKERRQ(ierr);
39520c24465SJunchao Zhang     ierr = PetscStrcasecmp("kokkos",backendstr,&isKokkos);CHKERRQ(ierr);
39659af0bd3SScott Kruger     ierr = PetscStrcasecmp("hip",backendstr,&isHip);CHKERRQ(ierr);
39759af0bd3SScott Kruger   #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
39820c24465SJunchao Zhang     if (isCuda) sf->backend = PETSCSF_BACKEND_CUDA;
39920c24465SJunchao Zhang     else if (isKokkos) sf->backend = PETSCSF_BACKEND_KOKKOS;
40059af0bd3SScott Kruger     else if (isHip) sf->backend = PETSCSF_BACKEND_HIP;
40159af0bd3SScott 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);
40220c24465SJunchao Zhang   #elif defined(PETSC_HAVE_KOKKOS)
40320c24465SJunchao Zhang     if (set && !isKokkos) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"-sf_backend %s is not supported. You can only choose kokkos", backendstr);
40420c24465SJunchao Zhang    #endif
40520c24465SJunchao Zhang   }
406c2a741eeSJunchao Zhang  #endif
407e55864a3SBarry Smith   if (sf->ops->SetFromOptions) {ierr = (*sf->ops->SetFromOptions)(PetscOptionsObject,sf);CHKERRQ(ierr);}
40895fce210SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
40995fce210SBarry Smith   PetscFunctionReturn(0);
41095fce210SBarry Smith }
41195fce210SBarry Smith 
41229046d53SLisandro Dalcin /*@
41395fce210SBarry Smith    PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order
41495fce210SBarry Smith 
41595fce210SBarry Smith    Logically Collective
41695fce210SBarry Smith 
4174165533cSJose E. Roman    Input Parameters:
41895fce210SBarry Smith +  sf - star forest
41995fce210SBarry Smith -  flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic)
42095fce210SBarry Smith 
42195fce210SBarry Smith    Level: advanced
42295fce210SBarry Smith 
42395fce210SBarry Smith .seealso: PetscSFGatherBegin(), PetscSFScatterBegin()
42495fce210SBarry Smith @*/
42595fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg)
42695fce210SBarry Smith {
42795fce210SBarry Smith   PetscFunctionBegin;
42895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
42995fce210SBarry Smith   PetscValidLogicalCollectiveBool(sf,flg,2);
43095fce210SBarry Smith   if (sf->multi) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
43195fce210SBarry Smith   sf->rankorder = flg;
43295fce210SBarry Smith   PetscFunctionReturn(0);
43395fce210SBarry Smith }
43495fce210SBarry Smith 
4358dbb0df6SBarry Smith /*@C
43695fce210SBarry Smith    PetscSFSetGraph - Set a parallel star forest
43795fce210SBarry Smith 
43895fce210SBarry Smith    Collective
43995fce210SBarry Smith 
4404165533cSJose E. Roman    Input Parameters:
44195fce210SBarry Smith +  sf - star forest
44295fce210SBarry Smith .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
44395fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
444c4e6a40aSLawrence Mitchell .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage (locations must be >= 0, enforced
445c4e6a40aSLawrence Mitchell during setup in debug mode)
44695fce210SBarry Smith .  localmode - copy mode for ilocal
447c4e6a40aSLawrence Mitchell .  iremote - remote locations of root vertices for each leaf on the current process (locations must be >= 0, enforced
448c4e6a40aSLawrence Mitchell during setup in debug mode)
44995fce210SBarry Smith -  remotemode - copy mode for iremote
45095fce210SBarry Smith 
45195fce210SBarry Smith    Level: intermediate
45295fce210SBarry Smith 
45395452b02SPatrick Sanan    Notes:
45495452b02SPatrick Sanan     In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode
45538ab3f8aSBarry Smith 
4562ad1e87fSLisandro Dalcin    Developers Note: Local indices which are the identity permutation in the range [0,nleaves) are discarded as they
4572ad1e87fSLisandro Dalcin    encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not
4582ad1e87fSLisandro Dalcin    needed
4592ad1e87fSLisandro Dalcin 
460c4e6a40aSLawrence Mitchell    Developers Note: This object does not necessarily encode a true star forest in the graph theoretic sense, since leaf
461c4e6a40aSLawrence Mitchell    indices are not required to be unique. Some functions, however, rely on unique leaf indices (checked in debug mode).
462c4e6a40aSLawrence Mitchell 
46395fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFGetGraph()
46495fce210SBarry Smith @*/
46595fce210SBarry Smith PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscSFNode *iremote,PetscCopyMode remotemode)
46695fce210SBarry Smith {
46795fce210SBarry Smith   PetscErrorCode ierr;
46895fce210SBarry Smith 
46995fce210SBarry Smith   PetscFunctionBegin;
47095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
47129046d53SLisandro Dalcin   if (nleaves > 0 && ilocal) PetscValidIntPointer(ilocal,4);
47229046d53SLisandro Dalcin   if (nleaves > 0) PetscValidPointer(iremote,6);
4732abc8c78SJacob Faibussowitsch   if (nroots  < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nroots %" PetscInt_FMT ", cannot be negative",nroots);
4742abc8c78SJacob Faibussowitsch   if (nleaves < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %" PetscInt_FMT ", cannot be negative",nleaves);
47529046d53SLisandro Dalcin 
4762a67d2daSStefano Zampini   if (sf->nroots >= 0) { /* Reset only if graph already set */
47795fce210SBarry Smith     ierr = PetscSFReset(sf);CHKERRQ(ierr);
4782a67d2daSStefano Zampini   }
4792a67d2daSStefano Zampini 
48029046d53SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
48129046d53SLisandro Dalcin 
48295fce210SBarry Smith   sf->nroots  = nroots;
48395fce210SBarry Smith   sf->nleaves = nleaves;
48429046d53SLisandro Dalcin 
48529046d53SLisandro Dalcin   if (nleaves && ilocal) {
48621c688dcSJed Brown     PetscInt i;
48729046d53SLisandro Dalcin     PetscInt minleaf = PETSC_MAX_INT;
48829046d53SLisandro Dalcin     PetscInt maxleaf = PETSC_MIN_INT;
4892ad1e87fSLisandro Dalcin     int      contiguous = 1;
49029046d53SLisandro Dalcin     for (i=0; i<nleaves; i++) {
49129046d53SLisandro Dalcin       minleaf = PetscMin(minleaf,ilocal[i]);
49229046d53SLisandro Dalcin       maxleaf = PetscMax(maxleaf,ilocal[i]);
4932ad1e87fSLisandro Dalcin       contiguous &= (ilocal[i] == i);
49429046d53SLisandro Dalcin     }
49529046d53SLisandro Dalcin     sf->minleaf = minleaf;
49629046d53SLisandro Dalcin     sf->maxleaf = maxleaf;
4972ad1e87fSLisandro Dalcin     if (contiguous) {
4982ad1e87fSLisandro Dalcin       if (localmode == PETSC_OWN_POINTER) {
4992ad1e87fSLisandro Dalcin         ierr = PetscFree(ilocal);CHKERRQ(ierr);
5002ad1e87fSLisandro Dalcin       }
5012ad1e87fSLisandro Dalcin       ilocal = NULL;
5022ad1e87fSLisandro Dalcin     }
50329046d53SLisandro Dalcin   } else {
50429046d53SLisandro Dalcin     sf->minleaf = 0;
50529046d53SLisandro Dalcin     sf->maxleaf = nleaves - 1;
50629046d53SLisandro Dalcin   }
50729046d53SLisandro Dalcin 
50829046d53SLisandro Dalcin   if (ilocal) {
50995fce210SBarry Smith     switch (localmode) {
51095fce210SBarry Smith     case PETSC_COPY_VALUES:
511785e854fSJed Brown       ierr = PetscMalloc1(nleaves,&sf->mine_alloc);CHKERRQ(ierr);
512580bdb30SBarry Smith       ierr = PetscArraycpy(sf->mine_alloc,ilocal,nleaves);CHKERRQ(ierr);
51395fce210SBarry Smith       sf->mine = sf->mine_alloc;
51495fce210SBarry Smith       break;
51595fce210SBarry Smith     case PETSC_OWN_POINTER:
51695fce210SBarry Smith       sf->mine_alloc = (PetscInt*)ilocal;
51795fce210SBarry Smith       sf->mine       = sf->mine_alloc;
51895fce210SBarry Smith       break;
51995fce210SBarry Smith     case PETSC_USE_POINTER:
52029046d53SLisandro Dalcin       sf->mine_alloc = NULL;
52195fce210SBarry Smith       sf->mine       = (PetscInt*)ilocal;
52295fce210SBarry Smith       break;
52395fce210SBarry Smith     default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown localmode");
52495fce210SBarry Smith     }
52595fce210SBarry Smith   }
52629046d53SLisandro Dalcin 
52795fce210SBarry Smith   switch (remotemode) {
52895fce210SBarry Smith   case PETSC_COPY_VALUES:
529785e854fSJed Brown     ierr = PetscMalloc1(nleaves,&sf->remote_alloc);CHKERRQ(ierr);
530580bdb30SBarry Smith     ierr = PetscArraycpy(sf->remote_alloc,iremote,nleaves);CHKERRQ(ierr);
53195fce210SBarry Smith     sf->remote = sf->remote_alloc;
53295fce210SBarry Smith     break;
53395fce210SBarry Smith   case PETSC_OWN_POINTER:
53495fce210SBarry Smith     sf->remote_alloc = (PetscSFNode*)iremote;
53595fce210SBarry Smith     sf->remote       = sf->remote_alloc;
53695fce210SBarry Smith     break;
53795fce210SBarry Smith   case PETSC_USE_POINTER:
53829046d53SLisandro Dalcin     sf->remote_alloc = NULL;
53995fce210SBarry Smith     sf->remote       = (PetscSFNode*)iremote;
54095fce210SBarry Smith     break;
54195fce210SBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_OUTOFRANGE,"Unknown remotemode");
54295fce210SBarry Smith   }
54395fce210SBarry Smith 
544563ffbabSMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0);CHKERRQ(ierr);
54529046d53SLisandro Dalcin   sf->graphset = PETSC_TRUE;
54695fce210SBarry Smith   PetscFunctionReturn(0);
54795fce210SBarry Smith }
54895fce210SBarry Smith 
54929046d53SLisandro Dalcin /*@
550dd5b3ca6SJunchao Zhang   PetscSFSetGraphWithPattern - Sets the graph of an SF with a specific pattern
551dd5b3ca6SJunchao Zhang 
552dd5b3ca6SJunchao Zhang   Collective
553dd5b3ca6SJunchao Zhang 
554dd5b3ca6SJunchao Zhang   Input Parameters:
555dd5b3ca6SJunchao Zhang + sf      - The PetscSF
556dd5b3ca6SJunchao Zhang . map     - Layout of roots over all processes (insignificant when pattern is PETSCSF_PATTERN_ALLTOALL)
557dd5b3ca6SJunchao Zhang - pattern - One of PETSCSF_PATTERN_ALLGATHER, PETSCSF_PATTERN_GATHER, PETSCSF_PATTERN_ALLTOALL
558dd5b3ca6SJunchao Zhang 
559dd5b3ca6SJunchao Zhang   Notes:
560dd5b3ca6SJunchao Zhang   It is easier to explain PetscSFPattern using vectors. Suppose we have an MPI vector x and its layout is map.
561dd5b3ca6SJunchao Zhang   n and N are local and global sizes of x respectively.
562dd5b3ca6SJunchao Zhang 
563dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLGATHER, the routine creates a graph that if one does Bcast on it, it will copy x to
564dd5b3ca6SJunchao Zhang   sequential vectors y on all ranks.
565dd5b3ca6SJunchao Zhang 
566dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_GATHER, the routine creates a graph that if one does Bcast on it, it will copy x to a
567dd5b3ca6SJunchao Zhang   sequential vector y on rank 0.
568dd5b3ca6SJunchao Zhang 
569dd5b3ca6SJunchao Zhang   In above cases, entries of x are roots and entries of y are leaves.
570dd5b3ca6SJunchao Zhang 
571dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLTOALL, map is insignificant. Suppose NP is size of sf's communicator. The routine
572dd5b3ca6SJunchao Zhang   creates a graph that every rank has NP leaves and NP roots. On rank i, its leaf j is connected to root i
573dd5b3ca6SJunchao 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
574dd5b3ca6SJunchao Zhang   not mean one can not send multiple items. One just needs to create a new MPI datatype for the mulptiple data
575dd5b3ca6SJunchao Zhang   items with MPI_Type_contiguous() and use that as the <unit> argument in SF routines.
576dd5b3ca6SJunchao Zhang 
577dd5b3ca6SJunchao Zhang   In this case, roots and leaves are symmetric.
578dd5b3ca6SJunchao Zhang 
579dd5b3ca6SJunchao Zhang   Level: intermediate
580dd5b3ca6SJunchao Zhang  @*/
581dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFSetGraphWithPattern(PetscSF sf,PetscLayout map,PetscSFPattern pattern)
582dd5b3ca6SJunchao Zhang {
583dd5b3ca6SJunchao Zhang   MPI_Comm       comm;
584dd5b3ca6SJunchao Zhang   PetscInt       n,N,res[2];
585dd5b3ca6SJunchao Zhang   PetscMPIInt    rank,size;
586dd5b3ca6SJunchao Zhang   PetscSFType    type;
587dd5b3ca6SJunchao Zhang   PetscErrorCode ierr;
588dd5b3ca6SJunchao Zhang 
589dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
5902abc8c78SJacob Faibussowitsch   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
5912abc8c78SJacob Faibussowitsch   if (pattern != PETSCSF_PATTERN_ALLTOALL) PetscValidPointer(map,2);
592dd5b3ca6SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf, &comm);CHKERRQ(ierr);
593*546078acSJacob Faibussowitsch   if (PetscUnlikely((pattern < PETSCSF_PATTERN_ALLGATHER) || (pattern > PETSCSF_PATTERN_ALLTOALL))) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Unsupported PetscSFPattern %d",pattern);
594ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
595ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
596dd5b3ca6SJunchao Zhang 
597dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLTOALL) {
598dd5b3ca6SJunchao Zhang     type = PETSCSFALLTOALL;
599dd5b3ca6SJunchao Zhang     ierr = PetscLayoutCreate(comm,&sf->map);CHKERRQ(ierr);
600dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetLocalSize(sf->map,size);CHKERRQ(ierr);
601dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetSize(sf->map,((PetscInt)size)*size);CHKERRQ(ierr);
602dd5b3ca6SJunchao Zhang     ierr = PetscLayoutSetUp(sf->map);CHKERRQ(ierr);
603dd5b3ca6SJunchao Zhang   } else {
604dd5b3ca6SJunchao Zhang     ierr   = PetscLayoutGetLocalSize(map,&n);CHKERRQ(ierr);
605dd5b3ca6SJunchao Zhang     ierr   = PetscLayoutGetSize(map,&N);CHKERRQ(ierr);
606dd5b3ca6SJunchao Zhang     res[0] = n;
607dd5b3ca6SJunchao Zhang     res[1] = -n;
608dd5b3ca6SJunchao Zhang     /* Check if n are same over all ranks so that we can optimize it */
609820f2d46SBarry Smith     ierr   = MPIU_Allreduce(MPI_IN_PLACE,res,2,MPIU_INT,MPI_MAX,comm);CHKERRMPI(ierr);
610dd5b3ca6SJunchao Zhang     if (res[0] == -res[1]) { /* same n */
611dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHER  : PETSCSFGATHER;
612dd5b3ca6SJunchao Zhang     } else {
613dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHERV : PETSCSFGATHERV;
614dd5b3ca6SJunchao Zhang     }
615dd5b3ca6SJunchao Zhang     ierr = PetscLayoutReference(map,&sf->map);CHKERRQ(ierr);
616dd5b3ca6SJunchao Zhang   }
617dd5b3ca6SJunchao Zhang   ierr = PetscSFSetType(sf,type);CHKERRQ(ierr);
618dd5b3ca6SJunchao Zhang 
619dd5b3ca6SJunchao Zhang   sf->pattern = pattern;
620dd5b3ca6SJunchao Zhang   sf->mine    = NULL; /* Contiguous */
621dd5b3ca6SJunchao Zhang 
622dd5b3ca6SJunchao Zhang   /* Set nleaves, nroots here in case user calls PetscSFGetGraph, which is legal to call even before PetscSFSetUp is called.
623dd5b3ca6SJunchao Zhang      Also set other easy stuff.
624dd5b3ca6SJunchao Zhang    */
625dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLGATHER) {
626dd5b3ca6SJunchao Zhang     sf->nleaves      = N;
627dd5b3ca6SJunchao Zhang     sf->nroots       = n;
628dd5b3ca6SJunchao Zhang     sf->nranks       = size;
629dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
630dd5b3ca6SJunchao Zhang     sf->maxleaf      = N - 1;
631dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_GATHER) {
632dd5b3ca6SJunchao Zhang     sf->nleaves      = rank ? 0 : N;
633dd5b3ca6SJunchao Zhang     sf->nroots       = n;
634dd5b3ca6SJunchao Zhang     sf->nranks       = rank ? 0 : size;
635dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
636dd5b3ca6SJunchao Zhang     sf->maxleaf      = rank ? -1 : N - 1;
637dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_ALLTOALL) {
638dd5b3ca6SJunchao Zhang     sf->nleaves      = size;
639dd5b3ca6SJunchao Zhang     sf->nroots       = size;
640dd5b3ca6SJunchao Zhang     sf->nranks       = size;
641dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
642dd5b3ca6SJunchao Zhang     sf->maxleaf      = size - 1;
643dd5b3ca6SJunchao Zhang   }
644dd5b3ca6SJunchao Zhang   sf->ndranks  = 0; /* We do not need to separate out distinguished ranks for patterned graphs to improve communication performance */
645dd5b3ca6SJunchao Zhang   sf->graphset = PETSC_TRUE;
646dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
647dd5b3ca6SJunchao Zhang }
648dd5b3ca6SJunchao Zhang 
649dd5b3ca6SJunchao Zhang /*@
65095fce210SBarry Smith    PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map
65195fce210SBarry Smith 
65295fce210SBarry Smith    Collective
65395fce210SBarry Smith 
6544165533cSJose E. Roman    Input Parameter:
65595fce210SBarry Smith .  sf - star forest to invert
65695fce210SBarry Smith 
6574165533cSJose E. Roman    Output Parameter:
65895fce210SBarry Smith .  isf - inverse of sf
6594165533cSJose E. Roman 
66095fce210SBarry Smith    Level: advanced
66195fce210SBarry Smith 
66295fce210SBarry Smith    Notes:
66395fce210SBarry Smith    All roots must have degree 1.
66495fce210SBarry Smith 
66595fce210SBarry Smith    The local space may be a permutation, but cannot be sparse.
66695fce210SBarry Smith 
66795fce210SBarry Smith .seealso: PetscSFSetGraph()
66895fce210SBarry Smith @*/
66995fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
67095fce210SBarry Smith {
67195fce210SBarry Smith   PetscErrorCode ierr;
67295fce210SBarry Smith   PetscMPIInt    rank;
67395fce210SBarry Smith   PetscInt       i,nroots,nleaves,maxlocal,count,*newilocal;
67495fce210SBarry Smith   const PetscInt *ilocal;
67595fce210SBarry Smith   PetscSFNode    *roots,*leaves;
67695fce210SBarry Smith 
67795fce210SBarry Smith   PetscFunctionBegin;
67829046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
67929046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
68029046d53SLisandro Dalcin   PetscValidPointer(isf,2);
68129046d53SLisandro Dalcin 
68295fce210SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
68329046d53SLisandro Dalcin   maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
68429046d53SLisandro Dalcin 
685ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRMPI(ierr);
686ae9aee6dSMatthew G. Knepley   ierr = PetscMalloc2(nroots,&roots,maxlocal,&leaves);CHKERRQ(ierr);
687ae9aee6dSMatthew G. Knepley   for (i=0; i<maxlocal; i++) {
68895fce210SBarry Smith     leaves[i].rank  = rank;
68995fce210SBarry Smith     leaves[i].index = i;
69095fce210SBarry Smith   }
69195fce210SBarry Smith   for (i=0; i <nroots; i++) {
69295fce210SBarry Smith     roots[i].rank  = -1;
69395fce210SBarry Smith     roots[i].index = -1;
69495fce210SBarry Smith   }
69583df288dSJunchao Zhang   ierr = PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPI_REPLACE);CHKERRQ(ierr);
69683df288dSJunchao Zhang   ierr = PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPI_REPLACE);CHKERRQ(ierr);
69795fce210SBarry Smith 
69895fce210SBarry Smith   /* Check whether our leaves are sparse */
69995fce210SBarry Smith   for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
70095fce210SBarry Smith   if (count == nroots) newilocal = NULL;
70195fce210SBarry Smith   else {                        /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
702785e854fSJed Brown     ierr = PetscMalloc1(count,&newilocal);CHKERRQ(ierr);
70395fce210SBarry Smith     for (i=0,count=0; i<nroots; i++) {
70495fce210SBarry Smith       if (roots[i].rank >= 0) {
70595fce210SBarry Smith         newilocal[count]   = i;
70695fce210SBarry Smith         roots[count].rank  = roots[i].rank;
70795fce210SBarry Smith         roots[count].index = roots[i].index;
70895fce210SBarry Smith         count++;
70995fce210SBarry Smith       }
71095fce210SBarry Smith     }
71195fce210SBarry Smith   }
71295fce210SBarry Smith 
71395fce210SBarry Smith   ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf);CHKERRQ(ierr);
71495fce210SBarry Smith   ierr = PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES);CHKERRQ(ierr);
71595fce210SBarry Smith   ierr = PetscFree2(roots,leaves);CHKERRQ(ierr);
71695fce210SBarry Smith   PetscFunctionReturn(0);
71795fce210SBarry Smith }
71895fce210SBarry Smith 
71995fce210SBarry Smith /*@
72095fce210SBarry Smith    PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph
72195fce210SBarry Smith 
72295fce210SBarry Smith    Collective
72395fce210SBarry Smith 
7244165533cSJose E. Roman    Input Parameters:
72595fce210SBarry Smith +  sf - communication object to duplicate
72695fce210SBarry Smith -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)
72795fce210SBarry Smith 
7284165533cSJose E. Roman    Output Parameter:
72995fce210SBarry Smith .  newsf - new communication object
73095fce210SBarry Smith 
73195fce210SBarry Smith    Level: beginner
73295fce210SBarry Smith 
73395fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetType(), PetscSFSetGraph()
73495fce210SBarry Smith @*/
73595fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
73695fce210SBarry Smith {
73729046d53SLisandro Dalcin   PetscSFType    type;
73895fce210SBarry Smith   PetscErrorCode ierr;
73997929ea7SJunchao Zhang   MPI_Datatype   dtype=MPIU_SCALAR;
74095fce210SBarry Smith 
74195fce210SBarry Smith   PetscFunctionBegin;
74229046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
74329046d53SLisandro Dalcin   PetscValidLogicalCollectiveEnum(sf,opt,2);
74429046d53SLisandro Dalcin   PetscValidPointer(newsf,3);
74595fce210SBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf);CHKERRQ(ierr);
74629046d53SLisandro Dalcin   ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr);
74729046d53SLisandro Dalcin   if (type) {ierr = PetscSFSetType(*newsf,type);CHKERRQ(ierr);}
74895fce210SBarry Smith   if (opt == PETSCSF_DUPLICATE_GRAPH) {
749dd5b3ca6SJunchao Zhang     PetscSFCheckGraphSet(sf,1);
750dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
75195fce210SBarry Smith       PetscInt          nroots,nleaves;
75295fce210SBarry Smith       const PetscInt    *ilocal;
75395fce210SBarry Smith       const PetscSFNode *iremote;
75495fce210SBarry Smith       ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
75595fce210SBarry Smith       ierr = PetscSFSetGraph(*newsf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr);
756dd5b3ca6SJunchao Zhang     } else {
757dd5b3ca6SJunchao Zhang       ierr = PetscSFSetGraphWithPattern(*newsf,sf->map,sf->pattern);CHKERRQ(ierr);
758dd5b3ca6SJunchao Zhang     }
75995fce210SBarry Smith   }
76097929ea7SJunchao Zhang   /* Since oldtype is committed, so is newtype, according to MPI */
76155b25c41SPierre Jolivet   if (sf->vscat.bs > 1) {ierr = MPI_Type_dup(sf->vscat.unit,&dtype);CHKERRMPI(ierr);}
76297929ea7SJunchao Zhang   (*newsf)->vscat.bs     = sf->vscat.bs;
76397929ea7SJunchao Zhang   (*newsf)->vscat.unit   = dtype;
76497929ea7SJunchao Zhang   (*newsf)->vscat.to_n   = sf->vscat.to_n;
76597929ea7SJunchao Zhang   (*newsf)->vscat.from_n = sf->vscat.from_n;
76697929ea7SJunchao Zhang   /* Do not copy lsf. Build it on demand since it is rarely used */
76797929ea7SJunchao Zhang 
76820c24465SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
76920c24465SJunchao Zhang   (*newsf)->backend              = sf->backend;
77071438e86SJunchao Zhang   (*newsf)->unknown_input_stream= sf->unknown_input_stream;
77120c24465SJunchao Zhang   (*newsf)->use_gpu_aware_mpi    = sf->use_gpu_aware_mpi;
77220c24465SJunchao Zhang   (*newsf)->use_stream_aware_mpi = sf->use_stream_aware_mpi;
77320c24465SJunchao Zhang #endif
77429046d53SLisandro Dalcin   if (sf->ops->Duplicate) {ierr = (*sf->ops->Duplicate)(sf,opt,*newsf);CHKERRQ(ierr);}
77520c24465SJunchao Zhang   /* Don't do PetscSFSetUp() since the new sf's graph might have not been set. */
77695fce210SBarry Smith   PetscFunctionReturn(0);
77795fce210SBarry Smith }
77895fce210SBarry Smith 
77995fce210SBarry Smith /*@C
78095fce210SBarry Smith    PetscSFGetGraph - Get the graph specifying a parallel star forest
78195fce210SBarry Smith 
78295fce210SBarry Smith    Not Collective
78395fce210SBarry Smith 
7844165533cSJose E. Roman    Input Parameter:
78595fce210SBarry Smith .  sf - star forest
78695fce210SBarry Smith 
7874165533cSJose E. Roman    Output Parameters:
78895fce210SBarry Smith +  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
78995fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
790bc6585dcSJunchao Zhang .  ilocal - locations of leaves in leafdata buffers (if returned value is NULL, it means leaves are in contiguous storage)
79195fce210SBarry Smith -  iremote - remote locations of root vertices for each leaf on the current process
79295fce210SBarry Smith 
793373e0d91SLisandro Dalcin    Notes:
794373e0d91SLisandro Dalcin      We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set yet
795373e0d91SLisandro Dalcin 
7968dbb0df6SBarry Smith    Fortran Notes:
7978dbb0df6SBarry Smith      The returned iremote array is a copy and must be deallocated after use. Consequently, if you
7988dbb0df6SBarry Smith      want to update the graph, you must call PetscSFSetGraph() after modifying the iremote array.
7998dbb0df6SBarry Smith 
8008dbb0df6SBarry Smith      To check for a NULL ilocal use
8018dbb0df6SBarry Smith $      if (loc(ilocal) == loc(PETSC_NULL_INTEGER)) then
802ca797d7aSLawrence Mitchell 
80395fce210SBarry Smith    Level: intermediate
80495fce210SBarry Smith 
80595fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph()
80695fce210SBarry Smith @*/
80795fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
80895fce210SBarry Smith {
809b8dee149SJunchao Zhang   PetscErrorCode ierr;
81095fce210SBarry Smith 
81195fce210SBarry Smith   PetscFunctionBegin;
81295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
813b8dee149SJunchao Zhang   if (sf->ops->GetGraph) {
814b8dee149SJunchao Zhang     ierr = (sf->ops->GetGraph)(sf,nroots,nleaves,ilocal,iremote);CHKERRQ(ierr);
815b8dee149SJunchao Zhang   } else {
81695fce210SBarry Smith     if (nroots) *nroots = sf->nroots;
81795fce210SBarry Smith     if (nleaves) *nleaves = sf->nleaves;
81895fce210SBarry Smith     if (ilocal) *ilocal = sf->mine;
81995fce210SBarry Smith     if (iremote) *iremote = sf->remote;
820b8dee149SJunchao Zhang   }
82195fce210SBarry Smith   PetscFunctionReturn(0);
82295fce210SBarry Smith }
82395fce210SBarry Smith 
82429046d53SLisandro Dalcin /*@
82595fce210SBarry Smith    PetscSFGetLeafRange - Get the active leaf ranges
82695fce210SBarry Smith 
82795fce210SBarry Smith    Not Collective
82895fce210SBarry Smith 
8294165533cSJose E. Roman    Input Parameter:
83095fce210SBarry Smith .  sf - star forest
83195fce210SBarry Smith 
8324165533cSJose E. Roman    Output Parameters:
833dd5b3ca6SJunchao Zhang +  minleaf - minimum active leaf on this process. Return 0 if there are no leaves.
834dd5b3ca6SJunchao Zhang -  maxleaf - maximum active leaf on this process. Return -1 if there are no leaves.
83595fce210SBarry Smith 
83695fce210SBarry Smith    Level: developer
83795fce210SBarry Smith 
83895fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
83995fce210SBarry Smith @*/
84095fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
84195fce210SBarry Smith {
84295fce210SBarry Smith   PetscFunctionBegin;
84395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
84429046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
84595fce210SBarry Smith   if (minleaf) *minleaf = sf->minleaf;
84695fce210SBarry Smith   if (maxleaf) *maxleaf = sf->maxleaf;
84795fce210SBarry Smith   PetscFunctionReturn(0);
84895fce210SBarry Smith }
84995fce210SBarry Smith 
85095fce210SBarry Smith /*@C
851fe2efc57SMark    PetscSFViewFromOptions - View from Options
852fe2efc57SMark 
853fe2efc57SMark    Collective on PetscSF
854fe2efc57SMark 
855fe2efc57SMark    Input Parameters:
856fe2efc57SMark +  A - the star forest
857736c3998SJose E. Roman .  obj - Optional object
858736c3998SJose E. Roman -  name - command line option
859fe2efc57SMark 
860fe2efc57SMark    Level: intermediate
861fe2efc57SMark .seealso:  PetscSF, PetscSFView, PetscObjectViewFromOptions(), PetscSFCreate()
862fe2efc57SMark @*/
863fe2efc57SMark PetscErrorCode  PetscSFViewFromOptions(PetscSF A,PetscObject obj,const char name[])
864fe2efc57SMark {
865fe2efc57SMark   PetscErrorCode ierr;
866fe2efc57SMark 
867fe2efc57SMark   PetscFunctionBegin;
868fe2efc57SMark   PetscValidHeaderSpecific(A,PETSCSF_CLASSID,1);
869fe2efc57SMark   ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr);
870fe2efc57SMark   PetscFunctionReturn(0);
871fe2efc57SMark }
872fe2efc57SMark 
873fe2efc57SMark /*@C
87495fce210SBarry Smith    PetscSFView - view a star forest
87595fce210SBarry Smith 
87695fce210SBarry Smith    Collective
87795fce210SBarry Smith 
8784165533cSJose E. Roman    Input Parameters:
87995fce210SBarry Smith +  sf - star forest
88095fce210SBarry Smith -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD
88195fce210SBarry Smith 
88295fce210SBarry Smith    Level: beginner
88395fce210SBarry Smith 
88495fce210SBarry Smith .seealso: PetscSFCreate(), PetscSFSetGraph()
88595fce210SBarry Smith @*/
88695fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
88795fce210SBarry Smith {
88895fce210SBarry Smith   PetscErrorCode    ierr;
88995fce210SBarry Smith   PetscBool         iascii;
89095fce210SBarry Smith   PetscViewerFormat format;
89195fce210SBarry Smith 
89295fce210SBarry Smith   PetscFunctionBegin;
89395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
89495fce210SBarry Smith   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer);CHKERRQ(ierr);}
89595fce210SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
89695fce210SBarry Smith   PetscCheckSameComm(sf,1,viewer,2);
89780153354SVaclav Hapla   if (sf->graphset) {ierr = PetscSFSetUp(sf);CHKERRQ(ierr);}
89895fce210SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
89953dd6d7dSJunchao Zhang   if (iascii && viewer->format != PETSC_VIEWER_ASCII_MATLAB) {
90095fce210SBarry Smith     PetscMPIInt rank;
90181bfa7aaSJed Brown     PetscInt    ii,i,j;
90295fce210SBarry Smith 
903dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer);CHKERRQ(ierr);
90495fce210SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
905dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
90680153354SVaclav Hapla       if (!sf->graphset) {
90780153354SVaclav Hapla         ierr = PetscViewerASCIIPrintf(viewer,"PetscSFSetGraph() has not been called yet\n");CHKERRQ(ierr);
90880153354SVaclav Hapla         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
90980153354SVaclav Hapla         PetscFunctionReturn(0);
91080153354SVaclav Hapla       }
911ffc4695bSBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRMPI(ierr);
9121575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
9132abc8c78SJacob Faibussowitsch       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%" PetscInt_FMT ", leaves=%" PetscInt_FMT ", remote ranks=%" PetscInt_FMT "\n",rank,sf->nroots,sf->nleaves,sf->nranks);CHKERRQ(ierr);
91495fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
9152abc8c78SJacob Faibussowitsch         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %" PetscInt_FMT " <- (%" PetscInt_FMT ",%" PetscInt_FMT ")\n",rank,sf->mine ? sf->mine[i] : i,sf->remote[i].rank,sf->remote[i].index);CHKERRQ(ierr);
91695fce210SBarry Smith       }
91795fce210SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
91895fce210SBarry Smith       ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
91995fce210SBarry Smith       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
92081bfa7aaSJed Brown         PetscMPIInt *tmpranks,*perm;
92181bfa7aaSJed Brown         ierr = PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm);CHKERRQ(ierr);
922580bdb30SBarry Smith         ierr = PetscArraycpy(tmpranks,sf->ranks,sf->nranks);CHKERRQ(ierr);
92381bfa7aaSJed Brown         for (i=0; i<sf->nranks; i++) perm[i] = i;
92481bfa7aaSJed Brown         ierr = PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm);CHKERRQ(ierr);
92595fce210SBarry Smith         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank);CHKERRQ(ierr);
92681bfa7aaSJed Brown         for (ii=0; ii<sf->nranks; ii++) {
92781bfa7aaSJed Brown           i = perm[ii];
9282abc8c78SJacob Faibussowitsch           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %" PetscInt_FMT " edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]);CHKERRQ(ierr);
92995fce210SBarry Smith           for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
9302abc8c78SJacob Faibussowitsch             ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d]    %" PetscInt_FMT " <- %" PetscInt_FMT "\n",rank,sf->rmine[j],sf->rremote[j]);CHKERRQ(ierr);
93195fce210SBarry Smith           }
93295fce210SBarry Smith         }
93381bfa7aaSJed Brown         ierr = PetscFree2(tmpranks,perm);CHKERRQ(ierr);
93495fce210SBarry Smith       }
93595fce210SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
9361575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
937dd5b3ca6SJunchao Zhang     }
93895fce210SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
93995fce210SBarry Smith   }
94062152dedSBarry Smith   if (sf->ops->View) {ierr = (*sf->ops->View)(sf,viewer);CHKERRQ(ierr);}
94195fce210SBarry Smith   PetscFunctionReturn(0);
94295fce210SBarry Smith }
94395fce210SBarry Smith 
94495fce210SBarry Smith /*@C
945dec1416fSJunchao Zhang    PetscSFGetRootRanks - Get root ranks and number of vertices referenced by leaves on this process
94695fce210SBarry Smith 
94795fce210SBarry Smith    Not Collective
94895fce210SBarry Smith 
9494165533cSJose E. Roman    Input Parameter:
95095fce210SBarry Smith .  sf - star forest
95195fce210SBarry Smith 
9524165533cSJose E. Roman    Output Parameters:
95395fce210SBarry Smith +  nranks - number of ranks referenced by local part
95495fce210SBarry Smith .  ranks - array of ranks
95595fce210SBarry Smith .  roffset - offset in rmine/rremote for each rank (length nranks+1)
95695fce210SBarry Smith .  rmine - concatenated array holding local indices referencing each remote rank
95795fce210SBarry Smith -  rremote - concatenated array holding remote indices referenced for each remote rank
95895fce210SBarry Smith 
95995fce210SBarry Smith    Level: developer
96095fce210SBarry Smith 
961dec1416fSJunchao Zhang .seealso: PetscSFGetLeafRanks()
96295fce210SBarry Smith @*/
963dec1416fSJunchao Zhang PetscErrorCode PetscSFGetRootRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
96495fce210SBarry Smith {
965dec1416fSJunchao Zhang   PetscErrorCode ierr;
96695fce210SBarry Smith 
96795fce210SBarry Smith   PetscFunctionBegin;
96895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
96929046d53SLisandro Dalcin   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
970dec1416fSJunchao Zhang   if (sf->ops->GetRootRanks) {
971dec1416fSJunchao Zhang     ierr = (sf->ops->GetRootRanks)(sf,nranks,ranks,roffset,rmine,rremote);CHKERRQ(ierr);
972dec1416fSJunchao Zhang   } else {
973dec1416fSJunchao Zhang     /* The generic implementation */
97495fce210SBarry Smith     if (nranks)  *nranks  = sf->nranks;
97595fce210SBarry Smith     if (ranks)   *ranks   = sf->ranks;
97695fce210SBarry Smith     if (roffset) *roffset = sf->roffset;
97795fce210SBarry Smith     if (rmine)   *rmine   = sf->rmine;
97895fce210SBarry Smith     if (rremote) *rremote = sf->rremote;
979dec1416fSJunchao Zhang   }
98095fce210SBarry Smith   PetscFunctionReturn(0);
98195fce210SBarry Smith }
98295fce210SBarry Smith 
9838750ddebSJunchao Zhang /*@C
9848750ddebSJunchao Zhang    PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process
9858750ddebSJunchao Zhang 
9868750ddebSJunchao Zhang    Not Collective
9878750ddebSJunchao Zhang 
9884165533cSJose E. Roman    Input Parameter:
9898750ddebSJunchao Zhang .  sf - star forest
9908750ddebSJunchao Zhang 
9914165533cSJose E. Roman    Output Parameters:
9928750ddebSJunchao Zhang +  niranks - number of leaf ranks referencing roots on this process
9938750ddebSJunchao Zhang .  iranks - array of ranks
9948750ddebSJunchao Zhang .  ioffset - offset in irootloc for each rank (length niranks+1)
9958750ddebSJunchao Zhang -  irootloc - concatenated array holding local indices of roots referenced by each leaf rank
9968750ddebSJunchao Zhang 
9978750ddebSJunchao Zhang    Level: developer
9988750ddebSJunchao Zhang 
999dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks()
10008750ddebSJunchao Zhang @*/
10018750ddebSJunchao Zhang PetscErrorCode PetscSFGetLeafRanks(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
10028750ddebSJunchao Zhang {
10038750ddebSJunchao Zhang   PetscErrorCode ierr;
10048750ddebSJunchao Zhang 
10058750ddebSJunchao Zhang   PetscFunctionBegin;
10068750ddebSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
10078750ddebSJunchao Zhang   if (!sf->setupcalled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
10088750ddebSJunchao Zhang   if (sf->ops->GetLeafRanks) {
10098750ddebSJunchao Zhang     ierr = (sf->ops->GetLeafRanks)(sf,niranks,iranks,ioffset,irootloc);CHKERRQ(ierr);
10108750ddebSJunchao Zhang   } else {
10118750ddebSJunchao Zhang     PetscSFType type;
10128750ddebSJunchao Zhang     ierr = PetscSFGetType(sf,&type);CHKERRQ(ierr);
10138750ddebSJunchao Zhang     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type);
10148750ddebSJunchao Zhang   }
10158750ddebSJunchao Zhang   PetscFunctionReturn(0);
10168750ddebSJunchao Zhang }
10178750ddebSJunchao Zhang 
1018b5a8e515SJed Brown static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) {
1019b5a8e515SJed Brown   PetscInt i;
1020b5a8e515SJed Brown   for (i=0; i<n; i++) {
1021b5a8e515SJed Brown     if (needle == list[i]) return PETSC_TRUE;
1022b5a8e515SJed Brown   }
1023b5a8e515SJed Brown   return PETSC_FALSE;
1024b5a8e515SJed Brown }
1025b5a8e515SJed Brown 
102695fce210SBarry Smith /*@C
102721c688dcSJed Brown    PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations.
102821c688dcSJed Brown 
102921c688dcSJed Brown    Collective
103021c688dcSJed Brown 
10314165533cSJose E. Roman    Input Parameters:
1032b5a8e515SJed Brown +  sf - PetscSF to set up; PetscSFSetGraph() must have been called
1033b5a8e515SJed Brown -  dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange)
103421c688dcSJed Brown 
103521c688dcSJed Brown    Level: developer
103621c688dcSJed Brown 
1037dec1416fSJunchao Zhang .seealso: PetscSFGetRootRanks()
103821c688dcSJed Brown @*/
1039b5a8e515SJed Brown PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup)
104021c688dcSJed Brown {
104121c688dcSJed Brown   PetscErrorCode     ierr;
104221c688dcSJed Brown   PetscTable         table;
104321c688dcSJed Brown   PetscTablePosition pos;
1044b5a8e515SJed Brown   PetscMPIInt        size,groupsize,*groupranks;
1045247e8311SStefano Zampini   PetscInt           *rcount,*ranks;
1046247e8311SStefano Zampini   PetscInt           i, irank = -1,orank = -1;
104721c688dcSJed Brown 
104821c688dcSJed Brown   PetscFunctionBegin;
104921c688dcSJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
105029046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1051ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size);CHKERRMPI(ierr);
105221c688dcSJed Brown   ierr = PetscTableCreate(10,size,&table);CHKERRQ(ierr);
105321c688dcSJed Brown   for (i=0; i<sf->nleaves; i++) {
105421c688dcSJed Brown     /* Log 1-based rank */
105521c688dcSJed Brown     ierr = PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES);CHKERRQ(ierr);
105621c688dcSJed Brown   }
105721c688dcSJed Brown   ierr = PetscTableGetCount(table,&sf->nranks);CHKERRQ(ierr);
105821c688dcSJed Brown   ierr = PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote);CHKERRQ(ierr);
105921c688dcSJed Brown   ierr = PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks);CHKERRQ(ierr);
106021c688dcSJed Brown   ierr = PetscTableGetHeadPosition(table,&pos);CHKERRQ(ierr);
106121c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
106221c688dcSJed Brown     ierr = PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]);CHKERRQ(ierr);
106321c688dcSJed Brown     ranks[i]--;             /* Convert back to 0-based */
106421c688dcSJed Brown   }
106521c688dcSJed Brown   ierr = PetscTableDestroy(&table);CHKERRQ(ierr);
1066b5a8e515SJed Brown 
1067b5a8e515SJed Brown   /* We expect that dgroup is reliably "small" while nranks could be large */
1068b5a8e515SJed Brown   {
10697fb8a5e4SKarl Rupp     MPI_Group group = MPI_GROUP_NULL;
1070b5a8e515SJed Brown     PetscMPIInt *dgroupranks;
1071ffc4695bSBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRMPI(ierr);
1072ffc4695bSBarry Smith     ierr = MPI_Group_size(dgroup,&groupsize);CHKERRMPI(ierr);
1073b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&dgroupranks);CHKERRQ(ierr);
1074b5a8e515SJed Brown     ierr = PetscMalloc1(groupsize,&groupranks);CHKERRQ(ierr);
1075b5a8e515SJed Brown     for (i=0; i<groupsize; i++) dgroupranks[i] = i;
1076ffc4695bSBarry Smith     if (groupsize) {ierr = MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks);CHKERRMPI(ierr);}
1077ffc4695bSBarry Smith     ierr = MPI_Group_free(&group);CHKERRMPI(ierr);
1078b5a8e515SJed Brown     ierr = PetscFree(dgroupranks);CHKERRQ(ierr);
1079b5a8e515SJed Brown   }
1080b5a8e515SJed Brown 
1081b5a8e515SJed Brown   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
1082b5a8e515SJed Brown   for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i;) {
1083b5a8e515SJed Brown     for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */
1084b5a8e515SJed Brown       if (InList(ranks[i],groupsize,groupranks)) break;
1085b5a8e515SJed Brown     }
1086b5a8e515SJed Brown     for (; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
1087b5a8e515SJed Brown       if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break;
1088b5a8e515SJed Brown     }
1089b5a8e515SJed Brown     if (sf->ndranks < i) {                         /* Swap ranks[sf->ndranks] with ranks[i] */
1090b5a8e515SJed Brown       PetscInt tmprank,tmpcount;
1091247e8311SStefano Zampini 
1092b5a8e515SJed Brown       tmprank             = ranks[i];
1093b5a8e515SJed Brown       tmpcount            = rcount[i];
1094b5a8e515SJed Brown       ranks[i]            = ranks[sf->ndranks];
1095b5a8e515SJed Brown       rcount[i]           = rcount[sf->ndranks];
1096b5a8e515SJed Brown       ranks[sf->ndranks]  = tmprank;
1097b5a8e515SJed Brown       rcount[sf->ndranks] = tmpcount;
1098b5a8e515SJed Brown       sf->ndranks++;
1099b5a8e515SJed Brown     }
1100b5a8e515SJed Brown   }
1101b5a8e515SJed Brown   ierr = PetscFree(groupranks);CHKERRQ(ierr);
1102b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->ndranks,ranks,rcount);CHKERRQ(ierr);
1103b5a8e515SJed Brown   ierr = PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks);CHKERRQ(ierr);
110421c688dcSJed Brown   sf->roffset[0] = 0;
110521c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
110621c688dcSJed Brown     ierr = PetscMPIIntCast(ranks[i],sf->ranks+i);CHKERRQ(ierr);
110721c688dcSJed Brown     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
110821c688dcSJed Brown     rcount[i]        = 0;
110921c688dcSJed Brown   }
1110247e8311SStefano Zampini   for (i=0, irank = -1, orank = -1; i<sf->nleaves; i++) {
1111247e8311SStefano Zampini     /* short circuit */
1112247e8311SStefano Zampini     if (orank != sf->remote[i].rank) {
111321c688dcSJed Brown       /* Search for index of iremote[i].rank in sf->ranks */
1114b5a8e515SJed Brown       ierr = PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank);CHKERRQ(ierr);
1115b5a8e515SJed Brown       if (irank < 0) {
1116b5a8e515SJed Brown         ierr = PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank);CHKERRQ(ierr);
1117b5a8e515SJed Brown         if (irank >= 0) irank += sf->ndranks;
111821c688dcSJed Brown       }
1119247e8311SStefano Zampini       orank = sf->remote[i].rank;
1120247e8311SStefano Zampini     }
11212abc8c78SJacob Faibussowitsch     if (irank < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %" PetscInt_FMT " in array",sf->remote[i].rank);
112221c688dcSJed Brown     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
112321c688dcSJed Brown     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
112421c688dcSJed Brown     rcount[irank]++;
112521c688dcSJed Brown   }
112621c688dcSJed Brown   ierr = PetscFree2(rcount,ranks);CHKERRQ(ierr);
112721c688dcSJed Brown   PetscFunctionReturn(0);
112821c688dcSJed Brown }
112921c688dcSJed Brown 
113021c688dcSJed Brown /*@C
113195fce210SBarry Smith    PetscSFGetGroups - gets incoming and outgoing process groups
113295fce210SBarry Smith 
113395fce210SBarry Smith    Collective
113495fce210SBarry Smith 
11354165533cSJose E. Roman    Input Parameter:
113695fce210SBarry Smith .  sf - star forest
113795fce210SBarry Smith 
11384165533cSJose E. Roman    Output Parameters:
113995fce210SBarry Smith +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
114095fce210SBarry Smith -  outgoing - group of destination processes for outgoing edges (roots that I reference)
114195fce210SBarry Smith 
114295fce210SBarry Smith    Level: developer
114395fce210SBarry Smith 
114495fce210SBarry Smith .seealso: PetscSFGetWindow(), PetscSFRestoreWindow()
114595fce210SBarry Smith @*/
114695fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
114795fce210SBarry Smith {
114895fce210SBarry Smith   PetscErrorCode ierr;
11497fb8a5e4SKarl Rupp   MPI_Group      group = MPI_GROUP_NULL;
115095fce210SBarry Smith 
115195fce210SBarry Smith   PetscFunctionBegin;
115244ee17edSStefano Zampini   if (sf->nranks < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUpRanks() before obtaining groups");
115395fce210SBarry Smith   if (sf->ingroup == MPI_GROUP_NULL) {
115495fce210SBarry Smith     PetscInt       i;
115595fce210SBarry Smith     const PetscInt *indegree;
115695fce210SBarry Smith     PetscMPIInt    rank,*outranks,*inranks;
115795fce210SBarry Smith     PetscSFNode    *remote;
115895fce210SBarry Smith     PetscSF        bgcount;
115995fce210SBarry Smith 
116095fce210SBarry Smith     /* Compute the number of incoming ranks */
1161785e854fSJed Brown     ierr = PetscMalloc1(sf->nranks,&remote);CHKERRQ(ierr);
116295fce210SBarry Smith     for (i=0; i<sf->nranks; i++) {
116395fce210SBarry Smith       remote[i].rank  = sf->ranks[i];
116495fce210SBarry Smith       remote[i].index = 0;
116595fce210SBarry Smith     }
116695fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount);CHKERRQ(ierr);
116795fce210SBarry Smith     ierr = PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
116895fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(bgcount,&indegree);CHKERRQ(ierr);
116995fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(bgcount,&indegree);CHKERRQ(ierr);
117095fce210SBarry Smith     /* Enumerate the incoming ranks */
1171dcca6d9dSJed Brown     ierr = PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks);CHKERRQ(ierr);
1172ffc4695bSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRMPI(ierr);
117395fce210SBarry Smith     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
117495fce210SBarry Smith     ierr = PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
117595fce210SBarry Smith     ierr = PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks);CHKERRQ(ierr);
1176ffc4695bSBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRMPI(ierr);
1177ffc4695bSBarry Smith     ierr = MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup);CHKERRMPI(ierr);
1178ffc4695bSBarry Smith     ierr = MPI_Group_free(&group);CHKERRMPI(ierr);
117995fce210SBarry Smith     ierr = PetscFree2(inranks,outranks);CHKERRQ(ierr);
118095fce210SBarry Smith     ierr = PetscSFDestroy(&bgcount);CHKERRQ(ierr);
118195fce210SBarry Smith   }
118295fce210SBarry Smith   *incoming = sf->ingroup;
118395fce210SBarry Smith 
118495fce210SBarry Smith   if (sf->outgroup == MPI_GROUP_NULL) {
1185ffc4695bSBarry Smith     ierr = MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group);CHKERRMPI(ierr);
1186ffc4695bSBarry Smith     ierr = MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup);CHKERRMPI(ierr);
1187ffc4695bSBarry Smith     ierr = MPI_Group_free(&group);CHKERRMPI(ierr);
118895fce210SBarry Smith   }
118995fce210SBarry Smith   *outgoing = sf->outgroup;
119095fce210SBarry Smith   PetscFunctionReturn(0);
119195fce210SBarry Smith }
119295fce210SBarry Smith 
119329046d53SLisandro Dalcin /*@
11943b8d980fSPierre Jolivet    PetscSFGetMultiSF - gets the inner SF implementing gathers and scatters
119595fce210SBarry Smith 
119695fce210SBarry Smith    Collective
119795fce210SBarry Smith 
11984165533cSJose E. Roman    Input Parameter:
119995fce210SBarry Smith .  sf - star forest that may contain roots with 0 or with more than 1 vertex
120095fce210SBarry Smith 
12014165533cSJose E. Roman    Output Parameter:
120295fce210SBarry Smith .  multi - star forest with split roots, such that each root has degree exactly 1
120395fce210SBarry Smith 
120495fce210SBarry Smith    Level: developer
120595fce210SBarry Smith 
120695fce210SBarry Smith    Notes:
120795fce210SBarry Smith 
120895fce210SBarry Smith    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
120995fce210SBarry Smith    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
121095fce210SBarry Smith    edge, it is a candidate for future optimization that might involve its removal.
121195fce210SBarry Smith 
1212673100f5SVaclav Hapla .seealso: PetscSFSetGraph(), PetscSFGatherBegin(), PetscSFScatterBegin(), PetscSFComputeMultiRootOriginalNumbering()
121395fce210SBarry Smith @*/
121495fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
121595fce210SBarry Smith {
121695fce210SBarry Smith   PetscErrorCode ierr;
121795fce210SBarry Smith 
121895fce210SBarry Smith   PetscFunctionBegin;
121995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
122095fce210SBarry Smith   PetscValidPointer(multi,2);
122195fce210SBarry Smith   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
122295fce210SBarry Smith     ierr   = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
122395fce210SBarry Smith     *multi = sf->multi;
1224013b3241SStefano Zampini     sf->multi->multi = sf->multi;
122595fce210SBarry Smith     PetscFunctionReturn(0);
122695fce210SBarry Smith   }
122795fce210SBarry Smith   if (!sf->multi) {
122895fce210SBarry Smith     const PetscInt *indegree;
12299837ea96SMatthew G. Knepley     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
123095fce210SBarry Smith     PetscSFNode    *remote;
123129046d53SLisandro Dalcin     maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
123295fce210SBarry Smith     ierr = PetscSFComputeDegreeBegin(sf,&indegree);CHKERRQ(ierr);
123395fce210SBarry Smith     ierr = PetscSFComputeDegreeEnd(sf,&indegree);CHKERRQ(ierr);
12349837ea96SMatthew G. Knepley     ierr = PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset);CHKERRQ(ierr);
123595fce210SBarry Smith     inoffset[0] = 0;
123695fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
12379837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) outones[i] = 1;
1238dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
1239dbd2ff41SBarry Smith     ierr = PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM);CHKERRQ(ierr);
124095fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
124176bd3646SJed Brown     if (PetscDefined(USE_DEBUG)) { /* Check that the expected number of increments occurred */
124295fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
124395fce210SBarry Smith         if (inoffset[i] + indegree[i] != inoffset[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
124495fce210SBarry Smith       }
124576bd3646SJed Brown     }
1246785e854fSJed Brown     ierr = PetscMalloc1(sf->nleaves,&remote);CHKERRQ(ierr);
124795fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
124895fce210SBarry Smith       remote[i].rank  = sf->remote[i].rank;
124938e7336fSToby Isaac       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
125095fce210SBarry Smith     }
125195fce210SBarry Smith     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi);CHKERRQ(ierr);
1252013b3241SStefano Zampini     sf->multi->multi = sf->multi;
125301365b40SToby Isaac     ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
125495fce210SBarry Smith     if (sf->rankorder) {        /* Sort the ranks */
125595fce210SBarry Smith       PetscMPIInt rank;
125695fce210SBarry Smith       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
125795fce210SBarry Smith       PetscSFNode *newremote;
1258ffc4695bSBarry Smith       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);CHKERRMPI(ierr);
125995fce210SBarry Smith       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
12609837ea96SMatthew G. Knepley       ierr = PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset);CHKERRQ(ierr);
12619837ea96SMatthew G. Knepley       for (i=0; i<maxlocal; i++) outranks[i] = rank;
126283df288dSJunchao Zhang       ierr = PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPI_REPLACE);CHKERRQ(ierr);
126383df288dSJunchao Zhang       ierr = PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPI_REPLACE);CHKERRQ(ierr);
126495fce210SBarry Smith       /* Sort the incoming ranks at each vertex, build the inverse map */
126595fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
126695fce210SBarry Smith         PetscInt j;
126795fce210SBarry Smith         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
126895fce210SBarry Smith         ierr = PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset);CHKERRQ(ierr);
126995fce210SBarry Smith         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
127095fce210SBarry Smith       }
1271ad227feaSJunchao Zhang       ierr = PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset,MPI_REPLACE);CHKERRQ(ierr);
1272ad227feaSJunchao Zhang       ierr = PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset,MPI_REPLACE);CHKERRQ(ierr);
1273785e854fSJed Brown       ierr = PetscMalloc1(sf->nleaves,&newremote);CHKERRQ(ierr);
127495fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
127595fce210SBarry Smith         newremote[i].rank  = sf->remote[i].rank;
127601365b40SToby Isaac         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
127795fce210SBarry Smith       }
127801365b40SToby Isaac       ierr = PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
127995fce210SBarry Smith       ierr = PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset);CHKERRQ(ierr);
128095fce210SBarry Smith     }
128195fce210SBarry Smith     ierr = PetscFree3(inoffset,outones,outoffset);CHKERRQ(ierr);
128295fce210SBarry Smith   }
128395fce210SBarry Smith   *multi = sf->multi;
128495fce210SBarry Smith   PetscFunctionReturn(0);
128595fce210SBarry Smith }
128695fce210SBarry Smith 
128795fce210SBarry Smith /*@C
128872502a1fSJunchao Zhang    PetscSFCreateEmbeddedRootSF - removes edges from all but the selected roots, does not remap indices
128995fce210SBarry Smith 
129095fce210SBarry Smith    Collective
129195fce210SBarry Smith 
12924165533cSJose E. Roman    Input Parameters:
129395fce210SBarry Smith +  sf - original star forest
1294ba2a7774SJunchao Zhang .  nselected  - number of selected roots on this process
1295ba2a7774SJunchao Zhang -  selected   - indices of the selected roots on this process
129695fce210SBarry Smith 
12974165533cSJose E. Roman    Output Parameter:
1298cd620004SJunchao Zhang .  esf - new star forest
129995fce210SBarry Smith 
130095fce210SBarry Smith    Level: advanced
130195fce210SBarry Smith 
130295fce210SBarry Smith    Note:
130395fce210SBarry Smith    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
130495fce210SBarry Smith    be done by calling PetscSFGetGraph().
130595fce210SBarry Smith 
130695fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFGetGraph()
130795fce210SBarry Smith @*/
130872502a1fSJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedRootSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *esf)
130995fce210SBarry Smith {
1310cd620004SJunchao Zhang   PetscInt          i,j,n,nroots,nleaves,esf_nleaves,*new_ilocal,minleaf,maxleaf,maxlocal;
1311cd620004SJunchao Zhang   const PetscInt    *ilocal;
1312cd620004SJunchao Zhang   signed char       *rootdata,*leafdata,*leafmem;
1313ba2a7774SJunchao Zhang   const PetscSFNode *iremote;
1314f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1315f659e5c7SJunchao Zhang   MPI_Comm          comm;
13160511a646SMatthew G. Knepley   PetscErrorCode    ierr;
131795fce210SBarry Smith 
131895fce210SBarry Smith   PetscFunctionBegin;
131995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
132029046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1321ba2a7774SJunchao Zhang   if (nselected) PetscValidPointer(selected,3);
1322cd620004SJunchao Zhang   PetscValidPointer(esf,4);
13230511a646SMatthew G. Knepley 
1324f659e5c7SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1325140a1472SStefano Zampini   ierr = PetscLogEventBegin(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr);
1326f659e5c7SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1327cd620004SJunchao Zhang   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
1328cd620004SJunchao Zhang 
132976bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {  /* Error out if selected[] has dups or  out of range indices */
1330cd620004SJunchao Zhang     PetscBool dups;
1331cd620004SJunchao Zhang     ierr = PetscCheckDupsInt(nselected,selected,&dups);CHKERRQ(ierr);
1332cd620004SJunchao Zhang     if (dups) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"selected[] has dups");
1333cd620004SJunchao Zhang     for (i=0; i<nselected; i++)
13342abc8c78SJacob Faibussowitsch       if (selected[i] < 0 || selected[i] >= nroots) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"selected root indice %" PetscInt_FMT " is out of [0,%" PetscInt_FMT ")",selected[i],nroots);
1335cd620004SJunchao Zhang   }
1336f659e5c7SJunchao Zhang 
133772502a1fSJunchao Zhang   if (sf->ops->CreateEmbeddedRootSF) {
133872502a1fSJunchao Zhang     ierr = (*sf->ops->CreateEmbeddedRootSF)(sf,nselected,selected,esf);CHKERRQ(ierr);
1339f659e5c7SJunchao Zhang   } else {
1340cd620004SJunchao Zhang     /* A generic version of creating embedded sf */
1341cd620004SJunchao Zhang     ierr = PetscSFGetLeafRange(sf,&minleaf,&maxleaf);CHKERRQ(ierr);
1342cd620004SJunchao Zhang     maxlocal = maxleaf - minleaf + 1;
1343cd620004SJunchao Zhang     ierr = PetscCalloc2(nroots,&rootdata,maxlocal,&leafmem);CHKERRQ(ierr);
1344cd620004SJunchao Zhang     leafdata = leafmem - minleaf;
1345cd620004SJunchao Zhang     /* Tag selected roots and bcast to leaves */
1346cd620004SJunchao Zhang     for (i=0; i<nselected; i++) rootdata[selected[i]] = 1;
1347ad227feaSJunchao Zhang     ierr = PetscSFBcastBegin(sf,MPI_SIGNED_CHAR,rootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr);
1348ad227feaSJunchao Zhang     ierr = PetscSFBcastEnd(sf,MPI_SIGNED_CHAR,rootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr);
1349ba2a7774SJunchao Zhang 
1350cd620004SJunchao Zhang     /* Build esf with leaves that are still connected */
1351cd620004SJunchao Zhang     esf_nleaves = 0;
1352cd620004SJunchao Zhang     for (i=0; i<nleaves; i++) {
1353cd620004SJunchao Zhang       j = ilocal ? ilocal[i] : i;
1354cd620004SJunchao Zhang       /* esf_nleaves += leafdata[j] should work in theory, but failed with SFWindow bugs
1355cd620004SJunchao Zhang          with PetscSFBcast. See https://gitlab.com/petsc/petsc/issues/555
1356cd620004SJunchao Zhang       */
1357cd620004SJunchao Zhang       esf_nleaves += (leafdata[j] ? 1 : 0);
1358cd620004SJunchao Zhang     }
1359cd620004SJunchao Zhang     ierr = PetscMalloc1(esf_nleaves,&new_ilocal);CHKERRQ(ierr);
1360cd620004SJunchao Zhang     ierr = PetscMalloc1(esf_nleaves,&new_iremote);CHKERRQ(ierr);
1361cd620004SJunchao Zhang     for (i=n=0; i<nleaves; i++) {
1362cd620004SJunchao Zhang       j = ilocal ? ilocal[i] : i;
1363cd620004SJunchao Zhang       if (leafdata[j]) {
1364cd620004SJunchao Zhang         new_ilocal[n]        = j;
1365cd620004SJunchao Zhang         new_iremote[n].rank  = iremote[i].rank;
1366cd620004SJunchao Zhang         new_iremote[n].index = iremote[i].index;
1367fc1ede2bSMatthew G. Knepley         ++n;
136895fce210SBarry Smith       }
136995fce210SBarry Smith     }
1370cd620004SJunchao Zhang     ierr = PetscSFCreate(comm,esf);CHKERRQ(ierr);
1371cd620004SJunchao Zhang     ierr = PetscSFSetFromOptions(*esf);CHKERRQ(ierr);
1372cd620004SJunchao Zhang     ierr = PetscSFSetGraph(*esf,nroots,esf_nleaves,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1373cd620004SJunchao Zhang     ierr = PetscFree2(rootdata,leafmem);CHKERRQ(ierr);
1374f659e5c7SJunchao Zhang   }
1375140a1472SStefano Zampini   ierr = PetscLogEventEnd(PETSCSF_EmbedSF,sf,0,0,0);CHKERRQ(ierr);
137695fce210SBarry Smith   PetscFunctionReturn(0);
137795fce210SBarry Smith }
137895fce210SBarry Smith 
13792f5fb4c2SMatthew G. Knepley /*@C
13802f5fb4c2SMatthew G. Knepley   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
13812f5fb4c2SMatthew G. Knepley 
13822f5fb4c2SMatthew G. Knepley   Collective
13832f5fb4c2SMatthew G. Knepley 
13844165533cSJose E. Roman   Input Parameters:
13852f5fb4c2SMatthew G. Knepley + sf - original star forest
1386f659e5c7SJunchao Zhang . nselected  - number of selected leaves on this process
1387f659e5c7SJunchao Zhang - selected   - indices of the selected leaves on this process
13882f5fb4c2SMatthew G. Knepley 
13894165533cSJose E. Roman   Output Parameter:
13902f5fb4c2SMatthew G. Knepley .  newsf - new star forest
13912f5fb4c2SMatthew G. Knepley 
13922f5fb4c2SMatthew G. Knepley   Level: advanced
13932f5fb4c2SMatthew G. Knepley 
139472502a1fSJunchao Zhang .seealso: PetscSFCreateEmbeddedRootSF(), PetscSFSetGraph(), PetscSFGetGraph()
13952f5fb4c2SMatthew G. Knepley @*/
1396f659e5c7SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
13972f5fb4c2SMatthew G. Knepley {
1398f659e5c7SJunchao Zhang   const PetscSFNode *iremote;
1399f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1400f659e5c7SJunchao Zhang   const PetscInt    *ilocal;
1401f659e5c7SJunchao Zhang   PetscInt          i,nroots,*leaves,*new_ilocal;
1402f659e5c7SJunchao Zhang   MPI_Comm          comm;
14032f5fb4c2SMatthew G. Knepley   PetscErrorCode    ierr;
14042f5fb4c2SMatthew G. Knepley 
14052f5fb4c2SMatthew G. Knepley   PetscFunctionBegin;
14062f5fb4c2SMatthew G. Knepley   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
140729046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1408f659e5c7SJunchao Zhang   if (nselected) PetscValidPointer(selected,3);
14092f5fb4c2SMatthew G. Knepley   PetscValidPointer(newsf,4);
14102f5fb4c2SMatthew G. Knepley 
1411f659e5c7SJunchao Zhang   /* Uniq selected[] and put results in leaves[] */
1412f659e5c7SJunchao Zhang   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1413f659e5c7SJunchao Zhang   ierr = PetscMalloc1(nselected,&leaves);CHKERRQ(ierr);
1414dd5b3ca6SJunchao Zhang   ierr = PetscArraycpy(leaves,selected,nselected);CHKERRQ(ierr);
1415f659e5c7SJunchao Zhang   ierr = PetscSortedRemoveDupsInt(&nselected,leaves);CHKERRQ(ierr);
14162abc8c78SJacob Faibussowitsch   if (nselected && (leaves[0] < 0 || leaves[nselected-1] >= sf->nleaves)) SETERRQ3(comm,PETSC_ERR_ARG_OUTOFRANGE,"Min/Max leaf indices %" PetscInt_FMT "/%" PetscInt_FMT " are not in [0,%" PetscInt_FMT ")",leaves[0],leaves[nselected-1],sf->nleaves);
1417f659e5c7SJunchao Zhang 
1418f659e5c7SJunchao Zhang   /* Optimize the routine only when sf is setup and hence we can reuse sf's communication pattern */
1419f659e5c7SJunchao Zhang   if (sf->setupcalled && sf->ops->CreateEmbeddedLeafSF) {
1420f659e5c7SJunchao Zhang     ierr = (*sf->ops->CreateEmbeddedLeafSF)(sf,nselected,leaves,newsf);CHKERRQ(ierr);
1421f659e5c7SJunchao Zhang   } else {
1422f659e5c7SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,NULL,&ilocal,&iremote);CHKERRQ(ierr);
1423f659e5c7SJunchao Zhang     ierr = PetscMalloc1(nselected,&new_ilocal);CHKERRQ(ierr);
1424f659e5c7SJunchao Zhang     ierr = PetscMalloc1(nselected,&new_iremote);CHKERRQ(ierr);
1425f659e5c7SJunchao Zhang     for (i=0; i<nselected; ++i) {
1426f659e5c7SJunchao Zhang       const PetscInt l     = leaves[i];
1427f659e5c7SJunchao Zhang       new_ilocal[i]        = ilocal ? ilocal[l] : l;
1428f659e5c7SJunchao Zhang       new_iremote[i].rank  = iremote[l].rank;
1429f659e5c7SJunchao Zhang       new_iremote[i].index = iremote[l].index;
14302f5fb4c2SMatthew G. Knepley     }
14314cc19a0cSStefano Zampini     ierr = PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,newsf);CHKERRQ(ierr);
1432f659e5c7SJunchao Zhang     ierr = PetscSFSetGraph(*newsf,nroots,nselected,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1433f659e5c7SJunchao Zhang   }
1434f659e5c7SJunchao Zhang   ierr = PetscFree(leaves);CHKERRQ(ierr);
14352f5fb4c2SMatthew G. Knepley   PetscFunctionReturn(0);
14362f5fb4c2SMatthew G. Knepley }
14372f5fb4c2SMatthew G. Knepley 
143895fce210SBarry Smith /*@C
1439ad227feaSJunchao Zhang    PetscSFBcastBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastEnd()
14403482bfa8SJunchao Zhang 
14413482bfa8SJunchao Zhang    Collective on PetscSF
14423482bfa8SJunchao Zhang 
14434165533cSJose E. Roman    Input Parameters:
14443482bfa8SJunchao Zhang +  sf - star forest on which to communicate
14453482bfa8SJunchao Zhang .  unit - data type associated with each node
14463482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
14473482bfa8SJunchao Zhang -  op - operation to use for reduction
14483482bfa8SJunchao Zhang 
14494165533cSJose E. Roman    Output Parameter:
14503482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
14513482bfa8SJunchao Zhang 
14523482bfa8SJunchao Zhang    Level: intermediate
14533482bfa8SJunchao Zhang 
1454d0295fc0SJunchao Zhang    Notes:
1455d0295fc0SJunchao Zhang     When petsc is configured with device support, it will use its own mechanism to figure out whether the given data pointers
1456d0295fc0SJunchao Zhang     are host pointers or device pointers, which may incur a noticable cost. If you already knew the info, you should
1457ad227feaSJunchao Zhang     use PetscSFBcastWithMemTypeBegin() instead.
1458ad227feaSJunchao Zhang .seealso: PetscSFBcastEnd(), PetscSFBcastWithMemTypeBegin()
14593482bfa8SJunchao Zhang @*/
1460ad227feaSJunchao Zhang PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
14613482bfa8SJunchao Zhang {
14623482bfa8SJunchao Zhang   PetscErrorCode ierr;
1463eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
14643482bfa8SJunchao Zhang 
14653482bfa8SJunchao Zhang   PetscFunctionBegin;
14663482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
14673482bfa8SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1468ad227feaSJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);}
1469eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1470eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1471ad227feaSJunchao Zhang   ierr = (*sf->ops->BcastBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op);CHKERRQ(ierr);
1472ad227feaSJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);}
14733482bfa8SJunchao Zhang   PetscFunctionReturn(0);
14743482bfa8SJunchao Zhang }
14753482bfa8SJunchao Zhang 
14763482bfa8SJunchao Zhang /*@C
1477ad227feaSJunchao Zhang    PetscSFBcastWithMemTypeBegin - begin pointwise broadcast with root value being reduced to leaf value with explicit memory types, to be concluded with call to PetscSFBcastEnd()
1478d0295fc0SJunchao Zhang 
1479d0295fc0SJunchao Zhang    Collective on PetscSF
1480d0295fc0SJunchao Zhang 
14814165533cSJose E. Roman    Input Parameters:
1482d0295fc0SJunchao Zhang +  sf - star forest on which to communicate
1483d0295fc0SJunchao Zhang .  unit - data type associated with each node
1484d0295fc0SJunchao Zhang .  rootmtype - memory type of rootdata
1485d0295fc0SJunchao Zhang .  rootdata - buffer to broadcast
1486d0295fc0SJunchao Zhang .  leafmtype - memory type of leafdata
1487d0295fc0SJunchao Zhang -  op - operation to use for reduction
1488d0295fc0SJunchao Zhang 
14894165533cSJose E. Roman    Output Parameter:
1490d0295fc0SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
1491d0295fc0SJunchao Zhang 
1492d0295fc0SJunchao Zhang    Level: intermediate
1493d0295fc0SJunchao Zhang 
1494ad227feaSJunchao Zhang .seealso: PetscSFBcastEnd(), PetscSFBcastBegin()
1495d0295fc0SJunchao Zhang @*/
1496ad227feaSJunchao Zhang PetscErrorCode PetscSFBcastWithMemTypeBegin(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,const void *rootdata,PetscMemType leafmtype,void *leafdata,MPI_Op op)
1497d0295fc0SJunchao Zhang {
1498d0295fc0SJunchao Zhang   PetscErrorCode ierr;
1499d0295fc0SJunchao Zhang 
1500d0295fc0SJunchao Zhang   PetscFunctionBegin;
1501d0295fc0SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1502d0295fc0SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1503ad227feaSJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);}
1504ad227feaSJunchao Zhang   ierr = (*sf->ops->BcastBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,op);CHKERRQ(ierr);
1505ad227feaSJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);}
1506d0295fc0SJunchao Zhang   PetscFunctionReturn(0);
1507d0295fc0SJunchao Zhang }
1508d0295fc0SJunchao Zhang 
1509d0295fc0SJunchao Zhang /*@C
1510ad227feaSJunchao Zhang    PetscSFBcastEnd - end a broadcast & reduce operation started with PetscSFBcastBegin()
15113482bfa8SJunchao Zhang 
15123482bfa8SJunchao Zhang    Collective
15133482bfa8SJunchao Zhang 
15144165533cSJose E. Roman    Input Parameters:
15153482bfa8SJunchao Zhang +  sf - star forest
15163482bfa8SJunchao Zhang .  unit - data type
15173482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
15183482bfa8SJunchao Zhang -  op - operation to use for reduction
15193482bfa8SJunchao Zhang 
15204165533cSJose E. Roman    Output Parameter:
15213482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
15223482bfa8SJunchao Zhang 
15233482bfa8SJunchao Zhang    Level: intermediate
15243482bfa8SJunchao Zhang 
15253482bfa8SJunchao Zhang .seealso: PetscSFSetGraph(), PetscSFReduceEnd()
15263482bfa8SJunchao Zhang @*/
1527ad227feaSJunchao Zhang PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
15283482bfa8SJunchao Zhang {
15293482bfa8SJunchao Zhang   PetscErrorCode ierr;
15303482bfa8SJunchao Zhang 
15313482bfa8SJunchao Zhang   PetscFunctionBegin;
15323482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1533ad227feaSJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);}
1534ad227feaSJunchao Zhang   ierr = (*sf->ops->BcastEnd)(sf,unit,rootdata,leafdata,op);CHKERRQ(ierr);
1535ad227feaSJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0);CHKERRQ(ierr);}
15363482bfa8SJunchao Zhang   PetscFunctionReturn(0);
15373482bfa8SJunchao Zhang }
15383482bfa8SJunchao Zhang 
15393482bfa8SJunchao Zhang /*@C
154095fce210SBarry Smith    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
154195fce210SBarry Smith 
154295fce210SBarry Smith    Collective
154395fce210SBarry Smith 
15444165533cSJose E. Roman    Input Parameters:
154595fce210SBarry Smith +  sf - star forest
154695fce210SBarry Smith .  unit - data type
154795fce210SBarry Smith .  leafdata - values to reduce
154895fce210SBarry Smith -  op - reduction operation
154995fce210SBarry Smith 
15504165533cSJose E. Roman    Output Parameter:
155195fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
155295fce210SBarry Smith 
155395fce210SBarry Smith    Level: intermediate
155495fce210SBarry Smith 
1555d0295fc0SJunchao Zhang    Notes:
1556d0295fc0SJunchao Zhang     When petsc is configured with device support, it will use its own mechanism to figure out whether the given data pointers
1557d0295fc0SJunchao Zhang     are host pointers or device pointers, which may incur a noticable cost. If you already knew the info, you should
1558d0295fc0SJunchao Zhang     use PetscSFReduceWithMemTypeBegin() instead.
1559d0295fc0SJunchao Zhang 
1560d0295fc0SJunchao Zhang .seealso: PetscSFBcastBegin(), PetscSFReduceWithMemTypeBegin()
156195fce210SBarry Smith @*/
156295fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
156395fce210SBarry Smith {
156495fce210SBarry Smith   PetscErrorCode ierr;
1565eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
156695fce210SBarry Smith 
156795fce210SBarry Smith   PetscFunctionBegin;
156895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
156995fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
157097929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);}
1571eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1572eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1573eb02082bSJunchao Zhang   ierr = (sf->ops->ReduceBegin)(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op);CHKERRQ(ierr);
157497929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);}
157595fce210SBarry Smith   PetscFunctionReturn(0);
157695fce210SBarry Smith }
157795fce210SBarry Smith 
157895fce210SBarry Smith /*@C
1579d0295fc0SJunchao Zhang    PetscSFReduceWithMemTypeBegin - begin reduction of leafdata into rootdata with explicit memory types, to be completed with call to PetscSFReduceEnd()
1580d0295fc0SJunchao Zhang 
1581d0295fc0SJunchao Zhang    Collective
1582d0295fc0SJunchao Zhang 
15834165533cSJose E. Roman    Input Parameters:
1584d0295fc0SJunchao Zhang +  sf - star forest
1585d0295fc0SJunchao Zhang .  unit - data type
1586d0295fc0SJunchao Zhang .  leafmtype - memory type of leafdata
1587d0295fc0SJunchao Zhang .  leafdata - values to reduce
1588d0295fc0SJunchao Zhang .  rootmtype - memory type of rootdata
1589d0295fc0SJunchao Zhang -  op - reduction operation
1590d0295fc0SJunchao Zhang 
15914165533cSJose E. Roman    Output Parameter:
1592d0295fc0SJunchao Zhang .  rootdata - result of reduction of values from all leaves of each root
1593d0295fc0SJunchao Zhang 
1594d0295fc0SJunchao Zhang    Level: intermediate
1595d0295fc0SJunchao Zhang 
1596d0295fc0SJunchao Zhang .seealso: PetscSFBcastBegin(), PetscSFReduceBegin()
1597d0295fc0SJunchao Zhang @*/
1598d0295fc0SJunchao Zhang PetscErrorCode PetscSFReduceWithMemTypeBegin(PetscSF sf,MPI_Datatype unit,PetscMemType leafmtype,const void *leafdata,PetscMemType rootmtype,void *rootdata,MPI_Op op)
1599d0295fc0SJunchao Zhang {
1600d0295fc0SJunchao Zhang   PetscErrorCode ierr;
1601d0295fc0SJunchao Zhang 
1602d0295fc0SJunchao Zhang   PetscFunctionBegin;
1603d0295fc0SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1604d0295fc0SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
160597929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);}
1606d0295fc0SJunchao Zhang   ierr = (sf->ops->ReduceBegin)(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op);CHKERRQ(ierr);
160797929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0);CHKERRQ(ierr);}
1608d0295fc0SJunchao Zhang   PetscFunctionReturn(0);
1609d0295fc0SJunchao Zhang }
1610d0295fc0SJunchao Zhang 
1611d0295fc0SJunchao Zhang /*@C
161295fce210SBarry Smith    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
161395fce210SBarry Smith 
161495fce210SBarry Smith    Collective
161595fce210SBarry Smith 
16164165533cSJose E. Roman    Input Parameters:
161795fce210SBarry Smith +  sf - star forest
161895fce210SBarry Smith .  unit - data type
161995fce210SBarry Smith .  leafdata - values to reduce
162095fce210SBarry Smith -  op - reduction operation
162195fce210SBarry Smith 
16224165533cSJose E. Roman    Output Parameter:
162395fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
162495fce210SBarry Smith 
162595fce210SBarry Smith    Level: intermediate
162695fce210SBarry Smith 
162795fce210SBarry Smith .seealso: PetscSFSetGraph(), PetscSFBcastEnd()
162895fce210SBarry Smith @*/
162995fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
163095fce210SBarry Smith {
163195fce210SBarry Smith   PetscErrorCode ierr;
163295fce210SBarry Smith 
163395fce210SBarry Smith   PetscFunctionBegin;
163495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
163597929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);}
163600816365SJunchao Zhang   ierr = (*sf->ops->ReduceEnd)(sf,unit,leafdata,rootdata,op);CHKERRQ(ierr);
163797929ea7SJunchao Zhang   if (!sf->vscat.logging) {ierr = PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0);CHKERRQ(ierr);}
163895fce210SBarry Smith   PetscFunctionReturn(0);
163995fce210SBarry Smith }
164095fce210SBarry Smith 
164195fce210SBarry Smith /*@C
1642a1729e3fSJunchao Zhang    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
1643a1729e3fSJunchao Zhang 
1644a1729e3fSJunchao Zhang    Collective
1645a1729e3fSJunchao Zhang 
16464165533cSJose E. Roman    Input Parameters:
1647a1729e3fSJunchao Zhang +  sf - star forest
1648a1729e3fSJunchao Zhang .  unit - data type
1649a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1650a1729e3fSJunchao Zhang -  op - operation to use for reduction
1651a1729e3fSJunchao Zhang 
16524165533cSJose E. Roman    Output Parameters:
1653a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1654a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1655a1729e3fSJunchao Zhang 
1656a1729e3fSJunchao Zhang    Level: advanced
1657a1729e3fSJunchao Zhang 
1658a1729e3fSJunchao Zhang    Note:
1659a1729e3fSJunchao Zhang    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
1660a1729e3fSJunchao Zhang    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
1661a1729e3fSJunchao Zhang    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
1662a1729e3fSJunchao Zhang    integers.
1663a1729e3fSJunchao Zhang 
1664a1729e3fSJunchao Zhang .seealso: PetscSFComputeDegreeBegin(), PetscSFReduceBegin(), PetscSFSetGraph()
1665a1729e3fSJunchao Zhang @*/
1666a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1667a1729e3fSJunchao Zhang {
1668a1729e3fSJunchao Zhang   PetscErrorCode ierr;
1669eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype,leafupdatemtype;
1670a1729e3fSJunchao Zhang 
1671a1729e3fSJunchao Zhang   PetscFunctionBegin;
1672a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1673a1729e3fSJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1674a1729e3fSJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1675eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
1676eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
1677eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafupdate,&leafupdatemtype);CHKERRQ(ierr);
1678eb02082bSJunchao Zhang   if (leafmtype != leafupdatemtype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for leafdata and leafupdate in different memory types");
1679eb02082bSJunchao Zhang   ierr = (*sf->ops->FetchAndOpBegin)(sf,unit,rootmtype,rootdata,leafmtype,leafdata,leafupdate,op);CHKERRQ(ierr);
1680a1729e3fSJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0);CHKERRQ(ierr);
1681a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1682a1729e3fSJunchao Zhang }
1683a1729e3fSJunchao Zhang 
1684a1729e3fSJunchao Zhang /*@C
1685a1729e3fSJunchao 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
1686a1729e3fSJunchao Zhang 
1687a1729e3fSJunchao Zhang    Collective
1688a1729e3fSJunchao Zhang 
16894165533cSJose E. Roman    Input Parameters:
1690a1729e3fSJunchao Zhang +  sf - star forest
1691a1729e3fSJunchao Zhang .  unit - data type
1692a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1693a1729e3fSJunchao Zhang -  op - operation to use for reduction
1694a1729e3fSJunchao Zhang 
16954165533cSJose E. Roman    Output Parameters:
1696a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1697a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1698a1729e3fSJunchao Zhang 
1699a1729e3fSJunchao Zhang    Level: advanced
1700a1729e3fSJunchao Zhang 
1701a1729e3fSJunchao Zhang .seealso: PetscSFComputeDegreeEnd(), PetscSFReduceEnd(), PetscSFSetGraph()
1702a1729e3fSJunchao Zhang @*/
1703a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1704a1729e3fSJunchao Zhang {
1705a1729e3fSJunchao Zhang   PetscErrorCode ierr;
1706a1729e3fSJunchao Zhang 
1707a1729e3fSJunchao Zhang   PetscFunctionBegin;
1708a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1709a1729e3fSJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
171000816365SJunchao Zhang   ierr = (*sf->ops->FetchAndOpEnd)(sf,unit,rootdata,leafdata,leafupdate,op);CHKERRQ(ierr);
1711a1729e3fSJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0);CHKERRQ(ierr);
1712a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1713a1729e3fSJunchao Zhang }
1714a1729e3fSJunchao Zhang 
1715a1729e3fSJunchao Zhang /*@C
171695fce210SBarry Smith    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
171795fce210SBarry Smith 
171895fce210SBarry Smith    Collective
171995fce210SBarry Smith 
17204165533cSJose E. Roman    Input Parameter:
172195fce210SBarry Smith .  sf - star forest
172295fce210SBarry Smith 
17234165533cSJose E. Roman    Output Parameter:
172495fce210SBarry Smith .  degree - degree of each root vertex
172595fce210SBarry Smith 
172695fce210SBarry Smith    Level: advanced
172795fce210SBarry Smith 
1728ffe67aa5SVáclav Hapla    Notes:
1729ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1730ffe67aa5SVáclav Hapla 
173195fce210SBarry Smith .seealso: PetscSFGatherBegin()
173295fce210SBarry Smith @*/
173395fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
173495fce210SBarry Smith {
173595fce210SBarry Smith   PetscErrorCode ierr;
173695fce210SBarry Smith 
173795fce210SBarry Smith   PetscFunctionBegin;
173895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
173995fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
174095fce210SBarry Smith   PetscValidPointer(degree,2);
1741803bd9e8SMatthew G. Knepley   if (!sf->degreeknown) {
17425b0d146aSStefano Zampini     PetscInt i, nroots = sf->nroots, maxlocal;
1743803bd9e8SMatthew G. Knepley     if (sf->degree) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
17445b0d146aSStefano Zampini     maxlocal = sf->maxleaf-sf->minleaf+1;
174529046d53SLisandro Dalcin     ierr = PetscMalloc1(nroots,&sf->degree);CHKERRQ(ierr);
174629046d53SLisandro Dalcin     ierr = PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp);CHKERRQ(ierr); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
174729046d53SLisandro Dalcin     for (i=0; i<nroots; i++) sf->degree[i] = 0;
17489837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
17495b0d146aSStefano Zampini     ierr = PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp-sf->minleaf,sf->degree,MPI_SUM);CHKERRQ(ierr);
175095fce210SBarry Smith   }
175195fce210SBarry Smith   *degree = NULL;
175295fce210SBarry Smith   PetscFunctionReturn(0);
175395fce210SBarry Smith }
175495fce210SBarry Smith 
175595fce210SBarry Smith /*@C
175695fce210SBarry Smith    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
175795fce210SBarry Smith 
175895fce210SBarry Smith    Collective
175995fce210SBarry Smith 
17604165533cSJose E. Roman    Input Parameter:
176195fce210SBarry Smith .  sf - star forest
176295fce210SBarry Smith 
17634165533cSJose E. Roman    Output Parameter:
176495fce210SBarry Smith .  degree - degree of each root vertex
176595fce210SBarry Smith 
176695fce210SBarry Smith    Level: developer
176795fce210SBarry Smith 
1768ffe67aa5SVáclav Hapla    Notes:
1769ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1770ffe67aa5SVáclav Hapla 
177195fce210SBarry Smith .seealso:
177295fce210SBarry Smith @*/
177395fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
177495fce210SBarry Smith {
177595fce210SBarry Smith   PetscErrorCode ierr;
177695fce210SBarry Smith 
177795fce210SBarry Smith   PetscFunctionBegin;
177895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
177995fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
178029046d53SLisandro Dalcin   PetscValidPointer(degree,2);
178195fce210SBarry Smith   if (!sf->degreeknown) {
178229046d53SLisandro Dalcin     if (!sf->degreetmp) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
17835b0d146aSStefano Zampini     ierr = PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp-sf->minleaf,sf->degree,MPI_SUM);CHKERRQ(ierr);
178495fce210SBarry Smith     ierr = PetscFree(sf->degreetmp);CHKERRQ(ierr);
178595fce210SBarry Smith     sf->degreeknown = PETSC_TRUE;
178695fce210SBarry Smith   }
178795fce210SBarry Smith   *degree = sf->degree;
178895fce210SBarry Smith   PetscFunctionReturn(0);
178995fce210SBarry Smith }
179095fce210SBarry Smith 
1791673100f5SVaclav Hapla /*@C
179266dfcd1aSVaclav Hapla    PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()).
179366dfcd1aSVaclav Hapla    Each multi-root is assigned index of the corresponding original root.
1794673100f5SVaclav Hapla 
1795673100f5SVaclav Hapla    Collective
1796673100f5SVaclav Hapla 
17974165533cSJose E. Roman    Input Parameters:
1798673100f5SVaclav Hapla +  sf - star forest
1799673100f5SVaclav Hapla -  degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd()
1800673100f5SVaclav Hapla 
18014165533cSJose E. Roman    Output Parameters:
180266dfcd1aSVaclav Hapla +  nMultiRoots - (optional) number of multi-roots (roots of multi-SF)
180366dfcd1aSVaclav Hapla -  multiRootsOrigNumbering - original indices of multi-roots; length of this array is nMultiRoots
1804673100f5SVaclav Hapla 
1805673100f5SVaclav Hapla    Level: developer
1806673100f5SVaclav Hapla 
1807ffe67aa5SVáclav Hapla    Notes:
1808ffe67aa5SVáclav Hapla    The returned array multiRootsOrigNumbering is newly allocated and should be destroyed with PetscFree() when no longer needed.
1809ffe67aa5SVáclav Hapla 
1810673100f5SVaclav Hapla .seealso: PetscSFComputeDegreeBegin(), PetscSFComputeDegreeEnd(), PetscSFGetMultiSF()
1811673100f5SVaclav Hapla @*/
181266dfcd1aSVaclav Hapla PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[])
1813673100f5SVaclav Hapla {
1814673100f5SVaclav Hapla   PetscSF             msf;
1815673100f5SVaclav Hapla   PetscInt            i, j, k, nroots, nmroots;
1816673100f5SVaclav Hapla   PetscErrorCode      ierr;
1817673100f5SVaclav Hapla 
1818673100f5SVaclav Hapla   PetscFunctionBegin;
1819673100f5SVaclav Hapla   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
1820673100f5SVaclav Hapla   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
182129328920SVaclav Hapla   if (nroots) PetscValidIntPointer(degree,2);
182266dfcd1aSVaclav Hapla   if (nMultiRoots) PetscValidIntPointer(nMultiRoots,3);
182366dfcd1aSVaclav Hapla   PetscValidPointer(multiRootsOrigNumbering,4);
1824673100f5SVaclav Hapla   ierr = PetscSFGetMultiSF(sf,&msf);CHKERRQ(ierr);
1825673100f5SVaclav Hapla   ierr = PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL);CHKERRQ(ierr);
182666dfcd1aSVaclav Hapla   ierr = PetscMalloc1(nmroots, multiRootsOrigNumbering);CHKERRQ(ierr);
1827673100f5SVaclav Hapla   for (i=0,j=0,k=0; i<nroots; i++) {
1828673100f5SVaclav Hapla     if (!degree[i]) continue;
1829673100f5SVaclav Hapla     for (j=0; j<degree[i]; j++,k++) {
183066dfcd1aSVaclav Hapla       (*multiRootsOrigNumbering)[k] = i;
1831673100f5SVaclav Hapla     }
1832673100f5SVaclav Hapla   }
1833673100f5SVaclav Hapla   if (PetscUnlikely(k != nmroots)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail");
183466dfcd1aSVaclav Hapla   if (nMultiRoots) *nMultiRoots = nmroots;
1835673100f5SVaclav Hapla   PetscFunctionReturn(0);
1836673100f5SVaclav Hapla }
1837673100f5SVaclav Hapla 
183895fce210SBarry Smith /*@C
183995fce210SBarry Smith    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
184095fce210SBarry Smith 
184195fce210SBarry Smith    Collective
184295fce210SBarry Smith 
18434165533cSJose E. Roman    Input Parameters:
184495fce210SBarry Smith +  sf - star forest
184595fce210SBarry Smith .  unit - data type
184695fce210SBarry Smith -  leafdata - leaf data to gather to roots
184795fce210SBarry Smith 
18484165533cSJose E. Roman    Output Parameter:
184995fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
185095fce210SBarry Smith 
185195fce210SBarry Smith    Level: intermediate
185295fce210SBarry Smith 
185395fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
185495fce210SBarry Smith @*/
185595fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
185695fce210SBarry Smith {
185795fce210SBarry Smith   PetscErrorCode ierr;
1858a5526d50SJunchao Zhang   PetscSF        multi = NULL;
185995fce210SBarry Smith 
186095fce210SBarry Smith   PetscFunctionBegin;
186195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
186229046d53SLisandro Dalcin   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
186395fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
186483df288dSJunchao Zhang   ierr = PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPI_REPLACE);CHKERRQ(ierr);
186595fce210SBarry Smith   PetscFunctionReturn(0);
186695fce210SBarry Smith }
186795fce210SBarry Smith 
186895fce210SBarry Smith /*@C
186995fce210SBarry Smith    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
187095fce210SBarry Smith 
187195fce210SBarry Smith    Collective
187295fce210SBarry Smith 
18734165533cSJose E. Roman    Input Parameters:
187495fce210SBarry Smith +  sf - star forest
187595fce210SBarry Smith .  unit - data type
187695fce210SBarry Smith -  leafdata - leaf data to gather to roots
187795fce210SBarry Smith 
18784165533cSJose E. Roman    Output Parameter:
187995fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
188095fce210SBarry Smith 
188195fce210SBarry Smith    Level: intermediate
188295fce210SBarry Smith 
188395fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
188495fce210SBarry Smith @*/
188595fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
188695fce210SBarry Smith {
188795fce210SBarry Smith   PetscErrorCode ierr;
1888a5526d50SJunchao Zhang   PetscSF        multi = NULL;
188995fce210SBarry Smith 
189095fce210SBarry Smith   PetscFunctionBegin;
189195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
189295fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
189383df288dSJunchao Zhang   ierr = PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPI_REPLACE);CHKERRQ(ierr);
189495fce210SBarry Smith   PetscFunctionReturn(0);
189595fce210SBarry Smith }
189695fce210SBarry Smith 
189795fce210SBarry Smith /*@C
189895fce210SBarry Smith    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
189995fce210SBarry Smith 
190095fce210SBarry Smith    Collective
190195fce210SBarry Smith 
19024165533cSJose E. Roman    Input Parameters:
190395fce210SBarry Smith +  sf - star forest
190495fce210SBarry Smith .  unit - data type
190595fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
190695fce210SBarry Smith 
19074165533cSJose E. Roman    Output Parameter:
190895fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
190995fce210SBarry Smith 
191095fce210SBarry Smith    Level: intermediate
191195fce210SBarry Smith 
191295fce210SBarry Smith .seealso: PetscSFComputeDegreeBegin(), PetscSFScatterBegin()
191395fce210SBarry Smith @*/
191495fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
191595fce210SBarry Smith {
191695fce210SBarry Smith   PetscErrorCode ierr;
1917a5526d50SJunchao Zhang   PetscSF        multi = NULL;
191895fce210SBarry Smith 
191995fce210SBarry Smith   PetscFunctionBegin;
192095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
192195fce210SBarry Smith   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
192295fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1923ad227feaSJunchao Zhang   ierr = PetscSFBcastBegin(multi,unit,multirootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr);
192495fce210SBarry Smith   PetscFunctionReturn(0);
192595fce210SBarry Smith }
192695fce210SBarry Smith 
192795fce210SBarry Smith /*@C
192895fce210SBarry Smith    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
192995fce210SBarry Smith 
193095fce210SBarry Smith    Collective
193195fce210SBarry Smith 
19324165533cSJose E. Roman    Input Parameters:
193395fce210SBarry Smith +  sf - star forest
193495fce210SBarry Smith .  unit - data type
193595fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
193695fce210SBarry Smith 
19374165533cSJose E. Roman    Output Parameter:
193895fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
193995fce210SBarry Smith 
194095fce210SBarry Smith    Level: intermediate
194195fce210SBarry Smith 
194295fce210SBarry Smith .seealso: PetscSFComputeDegreeEnd(), PetscSFScatterEnd()
194395fce210SBarry Smith @*/
194495fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
194595fce210SBarry Smith {
194695fce210SBarry Smith   PetscErrorCode ierr;
1947a5526d50SJunchao Zhang   PetscSF        multi = NULL;
194895fce210SBarry Smith 
194995fce210SBarry Smith   PetscFunctionBegin;
195095fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
195195fce210SBarry Smith   ierr = PetscSFGetMultiSF(sf,&multi);CHKERRQ(ierr);
1952ad227feaSJunchao Zhang   ierr = PetscSFBcastEnd(multi,unit,multirootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr);
195395fce210SBarry Smith   PetscFunctionReturn(0);
195495fce210SBarry Smith }
1955a7b3aa13SAta Mesgarnejad 
1956a072220fSLawrence Mitchell static PetscErrorCode PetscSFCheckLeavesUnique_Private(PetscSF sf)
1957a072220fSLawrence Mitchell {
1958a072220fSLawrence Mitchell   PetscInt        i, n, nleaves;
1959a072220fSLawrence Mitchell   const PetscInt *ilocal = NULL;
1960a072220fSLawrence Mitchell   PetscHSetI      seen;
1961a072220fSLawrence Mitchell   PetscErrorCode  ierr;
1962a072220fSLawrence Mitchell 
1963a072220fSLawrence Mitchell   PetscFunctionBegin;
1964b458e8f1SJose E. Roman   if (PetscDefined(USE_DEBUG)) {
1965a072220fSLawrence Mitchell     ierr = PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
1966a072220fSLawrence Mitchell     ierr = PetscHSetICreate(&seen);CHKERRQ(ierr);
1967a072220fSLawrence Mitchell     for (i = 0; i < nleaves; i++) {
1968a072220fSLawrence Mitchell       const PetscInt leaf = ilocal ? ilocal[i] : i;
1969a072220fSLawrence Mitchell       ierr = PetscHSetIAdd(seen,leaf);CHKERRQ(ierr);
1970a072220fSLawrence Mitchell     }
1971a072220fSLawrence Mitchell     ierr = PetscHSetIGetSize(seen,&n);CHKERRQ(ierr);
1972a072220fSLawrence Mitchell     if (n != nleaves) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided leaves have repeated values: all leaves must be unique");
1973a072220fSLawrence Mitchell     ierr = PetscHSetIDestroy(&seen);CHKERRQ(ierr);
1974b458e8f1SJose E. Roman   }
1975a072220fSLawrence Mitchell   PetscFunctionReturn(0);
1976a072220fSLawrence Mitchell }
197754729392SStefano Zampini 
1978a7b3aa13SAta Mesgarnejad /*@
197904c0ada0SJunchao Zhang   PetscSFCompose - Compose a new PetscSF by putting the second SF under the first one in a top (roots) down (leaves) view
1980a7b3aa13SAta Mesgarnejad 
1981a7b3aa13SAta Mesgarnejad   Input Parameters:
198254729392SStefano Zampini + sfA - The first PetscSF
198354729392SStefano Zampini - sfB - The second PetscSF
1984a7b3aa13SAta Mesgarnejad 
1985a7b3aa13SAta Mesgarnejad   Output Parameters:
198654729392SStefano Zampini . sfBA - The composite SF
1987a7b3aa13SAta Mesgarnejad 
1988a7b3aa13SAta Mesgarnejad   Level: developer
1989a7b3aa13SAta Mesgarnejad 
1990a072220fSLawrence Mitchell   Notes:
199154729392SStefano Zampini   Currently, the two SFs must be defined on congruent communicators and they must be true star
199254729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots.
199354729392SStefano Zampini 
199454729392SStefano Zampini   sfA's leaf space and sfB's root space might be partially overlapped. The composition builds
199554729392SStefano Zampini   a graph with sfA's roots and sfB's leaves only when there is a path between them. Unconnected
199654729392SStefano Zampini   nodes (roots or leaves) are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a
199754729392SStefano Zampini   Bcast on sfA, then a Bcast on sfB, on connected nodes.
1998a072220fSLawrence Mitchell 
199904c0ada0SJunchao Zhang .seealso: PetscSF, PetscSFComposeInverse(), PetscSFGetGraph(), PetscSFSetGraph()
2000a7b3aa13SAta Mesgarnejad @*/
2001a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
2002a7b3aa13SAta Mesgarnejad {
200304c0ada0SJunchao Zhang   PetscErrorCode    ierr;
2004a7b3aa13SAta Mesgarnejad   const PetscSFNode *remotePointsA,*remotePointsB;
2005d41018fbSJunchao Zhang   PetscSFNode       *remotePointsBA=NULL,*reorderedRemotePointsA = NULL,*leafdataB;
200654729392SStefano Zampini   const PetscInt    *localPointsA,*localPointsB;
200754729392SStefano Zampini   PetscInt          *localPointsBA;
200854729392SStefano Zampini   PetscInt          i,numRootsA,numLeavesA,numRootsB,numLeavesB,minleaf,maxleaf,numLeavesBA;
200954729392SStefano Zampini   PetscBool         denseB;
2010a7b3aa13SAta Mesgarnejad 
2011a7b3aa13SAta Mesgarnejad   PetscFunctionBegin;
2012a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1);
201329046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfA,1);
201429046d53SLisandro Dalcin   PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2);
201529046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfB,2);
201654729392SStefano Zampini   PetscCheckSameComm(sfA,1,sfB,2);
201729046d53SLisandro Dalcin   PetscValidPointer(sfBA,3);
201854729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfA);CHKERRQ(ierr);
201954729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfB);CHKERRQ(ierr);
202054729392SStefano Zampini 
2021a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfA,&numRootsA,&numLeavesA,&localPointsA,&remotePointsA);CHKERRQ(ierr);
2022a7b3aa13SAta Mesgarnejad   ierr = PetscSFGetGraph(sfB,&numRootsB,&numLeavesB,&localPointsB,&remotePointsB);CHKERRQ(ierr);
2023d41018fbSJunchao Zhang   if (localPointsA) {
202454729392SStefano Zampini     ierr = PetscMalloc1(numRootsB,&reorderedRemotePointsA);CHKERRQ(ierr);
202554729392SStefano Zampini     for (i=0; i<numRootsB; i++) {
202654729392SStefano Zampini       reorderedRemotePointsA[i].rank = -1;
202754729392SStefano Zampini       reorderedRemotePointsA[i].index = -1;
202854729392SStefano Zampini     }
202954729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
203054729392SStefano Zampini       if (localPointsA[i] >= numRootsB) continue;
203154729392SStefano Zampini       reorderedRemotePointsA[localPointsA[i]] = remotePointsA[i];
203254729392SStefano Zampini     }
2033d41018fbSJunchao Zhang     remotePointsA = reorderedRemotePointsA;
2034d41018fbSJunchao Zhang   }
2035d41018fbSJunchao Zhang   ierr = PetscSFGetLeafRange(sfB,&minleaf,&maxleaf);CHKERRQ(ierr);
2036d41018fbSJunchao Zhang   ierr = PetscMalloc1(maxleaf-minleaf+1,&leafdataB);CHKERRQ(ierr);
2037ad227feaSJunchao Zhang   ierr = PetscSFBcastBegin(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf,MPI_REPLACE);CHKERRQ(ierr);
2038ad227feaSJunchao Zhang   ierr = PetscSFBcastEnd(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf,MPI_REPLACE);CHKERRQ(ierr);
2039d41018fbSJunchao Zhang   ierr = PetscFree(reorderedRemotePointsA);CHKERRQ(ierr);
2040d41018fbSJunchao Zhang 
204154729392SStefano Zampini   denseB = (PetscBool)!localPointsB;
204254729392SStefano Zampini   for (i=0,numLeavesBA=0; i<numLeavesB; i++) {
204354729392SStefano Zampini     if (leafdataB[localPointsB ? localPointsB[i]-minleaf : i].rank == -1) denseB = PETSC_FALSE;
204454729392SStefano Zampini     else numLeavesBA++;
204554729392SStefano Zampini   }
204654729392SStefano Zampini   if (denseB) {
2047d41018fbSJunchao Zhang     localPointsBA  = NULL;
2048d41018fbSJunchao Zhang     remotePointsBA = leafdataB;
2049d41018fbSJunchao Zhang   } else {
205054729392SStefano Zampini     ierr = PetscMalloc1(numLeavesBA,&localPointsBA);CHKERRQ(ierr);
205154729392SStefano Zampini     ierr = PetscMalloc1(numLeavesBA,&remotePointsBA);CHKERRQ(ierr);
205254729392SStefano Zampini     for (i=0,numLeavesBA=0; i<numLeavesB; i++) {
205354729392SStefano Zampini       const PetscInt l = localPointsB ? localPointsB[i] : i;
205454729392SStefano Zampini 
205554729392SStefano Zampini       if (leafdataB[l-minleaf].rank == -1) continue;
205654729392SStefano Zampini       remotePointsBA[numLeavesBA] = leafdataB[l-minleaf];
205754729392SStefano Zampini       localPointsBA[numLeavesBA] = l;
205854729392SStefano Zampini       numLeavesBA++;
205954729392SStefano Zampini     }
2060d41018fbSJunchao Zhang     ierr = PetscFree(leafdataB);CHKERRQ(ierr);
2061d41018fbSJunchao Zhang   }
206254729392SStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA);CHKERRQ(ierr);
206320c24465SJunchao Zhang   ierr = PetscSFSetFromOptions(*sfBA);CHKERRQ(ierr);
206454729392SStefano Zampini   ierr = PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER);CHKERRQ(ierr);
2065a7b3aa13SAta Mesgarnejad   PetscFunctionReturn(0);
2066a7b3aa13SAta Mesgarnejad }
20671c6ba672SJunchao Zhang 
206804c0ada0SJunchao Zhang /*@
206954729392SStefano Zampini   PetscSFComposeInverse - Compose a new PetscSF by putting the inverse of the second SF under the first one
207004c0ada0SJunchao Zhang 
207104c0ada0SJunchao Zhang   Input Parameters:
207254729392SStefano Zampini + sfA - The first PetscSF
207354729392SStefano Zampini - sfB - The second PetscSF
207404c0ada0SJunchao Zhang 
207504c0ada0SJunchao Zhang   Output Parameters:
207654729392SStefano Zampini . sfBA - The composite SF.
207704c0ada0SJunchao Zhang 
207804c0ada0SJunchao Zhang   Level: developer
207904c0ada0SJunchao Zhang 
208054729392SStefano Zampini   Notes:
208154729392SStefano Zampini   Currently, the two SFs must be defined on congruent communicators and they must be true star
208254729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots. Even more, all roots of the
208354729392SStefano Zampini   second SF must have a degree of 1, i.e., no roots have more than one leaf connected.
208454729392SStefano Zampini 
208554729392SStefano Zampini   sfA's leaf space and sfB's leaf space might be partially overlapped. The composition builds
208654729392SStefano Zampini   a graph with sfA's roots and sfB's roots only when there is a path between them. Unconnected
208754729392SStefano Zampini   roots are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a Bcast on sfA, then
208854729392SStefano Zampini   a Reduce on sfB, on connected roots.
208954729392SStefano Zampini 
209054729392SStefano Zampini .seealso: PetscSF, PetscSFCompose(), PetscSFGetGraph(), PetscSFSetGraph(), PetscSFCreateInverseSF()
209104c0ada0SJunchao Zhang @*/
209204c0ada0SJunchao Zhang PetscErrorCode PetscSFComposeInverse(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
209304c0ada0SJunchao Zhang {
209404c0ada0SJunchao Zhang   PetscErrorCode    ierr;
209504c0ada0SJunchao Zhang   const PetscSFNode *remotePointsA,*remotePointsB;
209604c0ada0SJunchao Zhang   PetscSFNode       *remotePointsBA;
209704c0ada0SJunchao Zhang   const PetscInt    *localPointsA,*localPointsB;
209854729392SStefano Zampini   PetscSFNode       *reorderedRemotePointsA = NULL;
209954729392SStefano Zampini   PetscInt          i,numRootsA,numLeavesA,numLeavesBA,numRootsB,numLeavesB,minleaf,maxleaf,*localPointsBA;
21005b0d146aSStefano Zampini   MPI_Op            op;
21015b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
21025b0d146aSStefano Zampini   PetscBool         iswin;
21035b0d146aSStefano Zampini #endif
210404c0ada0SJunchao Zhang 
210504c0ada0SJunchao Zhang   PetscFunctionBegin;
210604c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1);
210704c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfA,1);
210804c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2);
210904c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfB,2);
211054729392SStefano Zampini   PetscCheckSameComm(sfA,1,sfB,2);
211104c0ada0SJunchao Zhang   PetscValidPointer(sfBA,3);
211254729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfA);CHKERRQ(ierr);
211354729392SStefano Zampini   ierr = PetscSFCheckLeavesUnique_Private(sfB);CHKERRQ(ierr);
211454729392SStefano Zampini 
211504c0ada0SJunchao Zhang   ierr = PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA);CHKERRQ(ierr);
211604c0ada0SJunchao Zhang   ierr = PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB);CHKERRQ(ierr);
21175b0d146aSStefano Zampini 
21185b0d146aSStefano Zampini   /* TODO: Check roots of sfB have degree of 1 */
21195b0d146aSStefano Zampini   /* Once we implement it, we can replace the MPI_MAXLOC
212083df288dSJunchao Zhang      with MPI_REPLACE. In that case, MPI_MAXLOC and MPI_REPLACE have the same effect.
21215b0d146aSStefano Zampini      We use MPI_MAXLOC only to have a deterministic output from this routine if
21225b0d146aSStefano Zampini      the root condition is not meet.
21235b0d146aSStefano Zampini    */
21245b0d146aSStefano Zampini   op = MPI_MAXLOC;
21255b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
21265b0d146aSStefano Zampini   /* we accept a non-deterministic output (if any) with PETSCSFWINDOW, since MPI_MAXLOC cannot operate on MPIU_2INT with MPI_Accumulate */
21275b0d146aSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)sfB,PETSCSFWINDOW,&iswin);CHKERRQ(ierr);
212883df288dSJunchao Zhang   if (iswin) op = MPI_REPLACE;
21295b0d146aSStefano Zampini #endif
21305b0d146aSStefano Zampini 
213154729392SStefano Zampini   ierr = PetscSFGetLeafRange(sfB, &minleaf, &maxleaf);CHKERRQ(ierr);
213254729392SStefano Zampini   ierr = PetscMalloc1(maxleaf - minleaf + 1,&reorderedRemotePointsA);CHKERRQ(ierr);
213354729392SStefano Zampini   for (i=0; i<maxleaf - minleaf + 1; i++) {
213454729392SStefano Zampini     reorderedRemotePointsA[i].rank = -1;
213554729392SStefano Zampini     reorderedRemotePointsA[i].index = -1;
213654729392SStefano Zampini   }
213754729392SStefano Zampini   if (localPointsA) {
213854729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
213954729392SStefano Zampini       if (localPointsA[i] > maxleaf || localPointsA[i] < minleaf) continue;
214054729392SStefano Zampini       reorderedRemotePointsA[localPointsA[i] - minleaf] = remotePointsA[i];
214154729392SStefano Zampini     }
214254729392SStefano Zampini   } else {
214354729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
214454729392SStefano Zampini       if (i > maxleaf || i < minleaf) continue;
214554729392SStefano Zampini       reorderedRemotePointsA[i - minleaf] = remotePointsA[i];
214654729392SStefano Zampini     }
214754729392SStefano Zampini   }
214854729392SStefano Zampini 
214954729392SStefano Zampini   ierr = PetscMalloc1(numRootsB,&localPointsBA);CHKERRQ(ierr);
215004c0ada0SJunchao Zhang   ierr = PetscMalloc1(numRootsB,&remotePointsBA);CHKERRQ(ierr);
215154729392SStefano Zampini   for (i=0; i<numRootsB; i++) {
215254729392SStefano Zampini     remotePointsBA[i].rank = -1;
215354729392SStefano Zampini     remotePointsBA[i].index = -1;
215454729392SStefano Zampini   }
215554729392SStefano Zampini 
21565b0d146aSStefano Zampini   ierr = PetscSFReduceBegin(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,op);CHKERRQ(ierr);
21575b0d146aSStefano Zampini   ierr = PetscSFReduceEnd(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,op);CHKERRQ(ierr);
215854729392SStefano Zampini   ierr = PetscFree(reorderedRemotePointsA);CHKERRQ(ierr);
215954729392SStefano Zampini   for (i=0,numLeavesBA=0; i<numRootsB; i++) {
216054729392SStefano Zampini     if (remotePointsBA[i].rank == -1) continue;
216154729392SStefano Zampini     remotePointsBA[numLeavesBA].rank = remotePointsBA[i].rank;
216254729392SStefano Zampini     remotePointsBA[numLeavesBA].index = remotePointsBA[i].index;
216354729392SStefano Zampini     localPointsBA[numLeavesBA] = i;
216454729392SStefano Zampini     numLeavesBA++;
216554729392SStefano Zampini   }
216654729392SStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA);CHKERRQ(ierr);
216720c24465SJunchao Zhang   ierr = PetscSFSetFromOptions(*sfBA);CHKERRQ(ierr);
216854729392SStefano Zampini   ierr = PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER);CHKERRQ(ierr);
216904c0ada0SJunchao Zhang   PetscFunctionReturn(0);
217004c0ada0SJunchao Zhang }
217104c0ada0SJunchao Zhang 
21721c6ba672SJunchao Zhang /*
21731c6ba672SJunchao Zhang   PetscSFCreateLocalSF_Private - Creates a local PetscSF that only has intra-process edges of the global PetscSF
21741c6ba672SJunchao Zhang 
21751c6ba672SJunchao Zhang   Input Parameters:
21761c6ba672SJunchao Zhang . sf - The global PetscSF
21771c6ba672SJunchao Zhang 
21781c6ba672SJunchao Zhang   Output Parameters:
21791c6ba672SJunchao Zhang . out - The local PetscSF
21801c6ba672SJunchao Zhang  */
21811c6ba672SJunchao Zhang PetscErrorCode PetscSFCreateLocalSF_Private(PetscSF sf,PetscSF *out)
21821c6ba672SJunchao Zhang {
21831c6ba672SJunchao Zhang   MPI_Comm           comm;
21841c6ba672SJunchao Zhang   PetscMPIInt        myrank;
21851c6ba672SJunchao Zhang   const PetscInt     *ilocal;
21861c6ba672SJunchao Zhang   const PetscSFNode  *iremote;
21871c6ba672SJunchao Zhang   PetscInt           i,j,nroots,nleaves,lnleaves,*lilocal;
21881c6ba672SJunchao Zhang   PetscSFNode        *liremote;
21891c6ba672SJunchao Zhang   PetscSF            lsf;
21901c6ba672SJunchao Zhang   PetscErrorCode     ierr;
21911c6ba672SJunchao Zhang 
21921c6ba672SJunchao Zhang   PetscFunctionBegin;
21931c6ba672SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
21941c6ba672SJunchao Zhang   if (sf->ops->CreateLocalSF) {
21951c6ba672SJunchao Zhang     ierr = (*sf->ops->CreateLocalSF)(sf,out);CHKERRQ(ierr);
21961c6ba672SJunchao Zhang   } else {
21971c6ba672SJunchao Zhang     /* Could use PetscSFCreateEmbeddedLeafSF, but since we know the comm is PETSC_COMM_SELF, we can make it fast */
21981c6ba672SJunchao Zhang     ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
2199ffc4695bSBarry Smith     ierr = MPI_Comm_rank(comm,&myrank);CHKERRMPI(ierr);
22001c6ba672SJunchao Zhang 
22011c6ba672SJunchao Zhang     /* Find out local edges and build a local SF */
22021c6ba672SJunchao Zhang     ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
22031c6ba672SJunchao Zhang     for (i=lnleaves=0; i<nleaves; i++) {if (iremote[i].rank == (PetscInt)myrank) lnleaves++;}
22041c6ba672SJunchao Zhang     ierr = PetscMalloc1(lnleaves,&lilocal);CHKERRQ(ierr);
22051c6ba672SJunchao Zhang     ierr = PetscMalloc1(lnleaves,&liremote);CHKERRQ(ierr);
22061c6ba672SJunchao Zhang 
22071c6ba672SJunchao Zhang     for (i=j=0; i<nleaves; i++) {
22081c6ba672SJunchao Zhang       if (iremote[i].rank == (PetscInt)myrank) {
22091c6ba672SJunchao Zhang         lilocal[j]        = ilocal? ilocal[i] : i; /* ilocal=NULL for contiguous storage */
22101c6ba672SJunchao Zhang         liremote[j].rank  = 0; /* rank in PETSC_COMM_SELF */
22111c6ba672SJunchao Zhang         liremote[j].index = iremote[i].index;
22121c6ba672SJunchao Zhang         j++;
22131c6ba672SJunchao Zhang       }
22141c6ba672SJunchao Zhang     }
22151c6ba672SJunchao Zhang     ierr = PetscSFCreate(PETSC_COMM_SELF,&lsf);CHKERRQ(ierr);
221620c24465SJunchao Zhang     ierr = PetscSFSetFromOptions(lsf);CHKERRQ(ierr);
22171c6ba672SJunchao Zhang     ierr = PetscSFSetGraph(lsf,nroots,lnleaves,lilocal,PETSC_OWN_POINTER,liremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
22181c6ba672SJunchao Zhang     ierr = PetscSFSetUp(lsf);CHKERRQ(ierr);
22191c6ba672SJunchao Zhang     *out = lsf;
22201c6ba672SJunchao Zhang   }
22211c6ba672SJunchao Zhang   PetscFunctionReturn(0);
22221c6ba672SJunchao Zhang }
2223dd5b3ca6SJunchao Zhang 
2224dd5b3ca6SJunchao Zhang /* Similar to PetscSFBcast, but only Bcast to leaves on rank 0 */
2225dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFBcastToZero_Private(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
2226dd5b3ca6SJunchao Zhang {
2227dd5b3ca6SJunchao Zhang   PetscErrorCode     ierr;
2228eb02082bSJunchao Zhang   PetscMemType       rootmtype,leafmtype;
2229dd5b3ca6SJunchao Zhang 
2230dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
2231dd5b3ca6SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
2232dd5b3ca6SJunchao Zhang   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
2233ad227feaSJunchao Zhang   ierr = PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
2234eb02082bSJunchao Zhang   ierr = PetscGetMemType(rootdata,&rootmtype);CHKERRQ(ierr);
2235eb02082bSJunchao Zhang   ierr = PetscGetMemType(leafdata,&leafmtype);CHKERRQ(ierr);
2236dd5b3ca6SJunchao Zhang   if (sf->ops->BcastToZero) {
2237eb02082bSJunchao Zhang     ierr = (*sf->ops->BcastToZero)(sf,unit,rootmtype,rootdata,leafmtype,leafdata);CHKERRQ(ierr);
2238dd5b3ca6SJunchao Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFBcastToZero_Private is not supported on this SF type");
2239ad227feaSJunchao Zhang   ierr = PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0);CHKERRQ(ierr);
2240dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
2241dd5b3ca6SJunchao Zhang }
2242dd5b3ca6SJunchao Zhang 
2243