xref: /petsc/src/vec/is/sf/interface/sf.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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)
187a46b595SBarry Smith     #define PetscSFCheckGraphSet(sf, arg) PetscCheck((sf)->graphset, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFSetGraph() or PetscSFSetGraphWithPattern() on argument %d \"%s\" before %s()", (arg), #sf, PETSC_FUNCTION_NAME);
1995fce210SBarry Smith   #else
209371c9d4SSatish Balay     #define PetscSFCheckGraphSet(sf, arg) \
219371c9d4SSatish Balay       do { \
229371c9d4SSatish Balay       } while (0)
2395fce210SBarry Smith   #endif
242abc8c78SJacob Faibussowitsch #endif
2595fce210SBarry Smith 
264c8fdceaSLisandro Dalcin const char *const PetscSFDuplicateOptions[] = {"CONFONLY", "RANKS", "GRAPH", "PetscSFDuplicateOption", "PETSCSF_DUPLICATE_", NULL};
2795fce210SBarry Smith 
288af6ec1cSBarry Smith /*@
2995fce210SBarry Smith    PetscSFCreate - create a star forest communication context
3095fce210SBarry Smith 
31d083f849SBarry Smith    Collective
3295fce210SBarry Smith 
334165533cSJose E. Roman    Input Parameter:
3495fce210SBarry Smith .  comm - communicator on which the star forest will operate
3595fce210SBarry Smith 
364165533cSJose E. Roman    Output Parameter:
3795fce210SBarry Smith .  sf - new star forest context
3895fce210SBarry Smith 
39dd5b3ca6SJunchao Zhang    Options Database Keys:
40dd5b3ca6SJunchao Zhang +  -sf_type basic     -Use MPI persistent Isend/Irecv for communication (Default)
41dd5b3ca6SJunchao Zhang .  -sf_type window    -Use MPI-3 one-sided window for communication
42dd5b3ca6SJunchao Zhang -  -sf_type neighbor  -Use MPI-3 neighborhood collectives for communication
43dd5b3ca6SJunchao Zhang 
4495fce210SBarry Smith    Level: intermediate
4595fce210SBarry Smith 
46dd5b3ca6SJunchao Zhang    Notes:
47dd5b3ca6SJunchao Zhang    When one knows the communication graph is one of the predefined graph, such as MPI_Alltoall, MPI_Allgatherv,
48dd5b3ca6SJunchao Zhang    MPI_Gatherv, one can create a PetscSF and then set its graph with PetscSFSetGraphWithPattern(). These special
49dd5b3ca6SJunchao Zhang    SFs are optimized and they have better performance than general SFs.
50dd5b3ca6SJunchao Zhang 
51db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`, `PetscSFSetGraphWithPattern()`, `PetscSFDestroy()`
5295fce210SBarry Smith @*/
53*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreate(MPI_Comm comm, PetscSF *sf)
54*d71ae5a4SJacob Faibussowitsch {
5595fce210SBarry Smith   PetscSF b;
5695fce210SBarry Smith 
5795fce210SBarry Smith   PetscFunctionBegin;
5895fce210SBarry Smith   PetscValidPointer(sf, 2);
599566063dSJacob Faibussowitsch   PetscCall(PetscSFInitializePackage());
6095fce210SBarry Smith 
619566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(b, PETSCSF_CLASSID, "PetscSF", "Star Forest", "PetscSF", comm, PetscSFDestroy, PetscSFView));
6295fce210SBarry Smith 
6395fce210SBarry Smith   b->nroots    = -1;
6495fce210SBarry Smith   b->nleaves   = -1;
6529046d53SLisandro Dalcin   b->minleaf   = PETSC_MAX_INT;
6629046d53SLisandro Dalcin   b->maxleaf   = PETSC_MIN_INT;
6795fce210SBarry Smith   b->nranks    = -1;
6895fce210SBarry Smith   b->rankorder = PETSC_TRUE;
6995fce210SBarry Smith   b->ingroup   = MPI_GROUP_NULL;
7095fce210SBarry Smith   b->outgroup  = MPI_GROUP_NULL;
7195fce210SBarry Smith   b->graphset  = PETSC_FALSE;
7220c24465SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
7320c24465SJunchao Zhang   b->use_gpu_aware_mpi    = use_gpu_aware_mpi;
7420c24465SJunchao Zhang   b->use_stream_aware_mpi = PETSC_FALSE;
7571438e86SJunchao Zhang   b->unknown_input_stream = PETSC_FALSE;
7627f636e8SJunchao Zhang   #if defined(PETSC_HAVE_KOKKOS) /* Prefer kokkos over cuda*/
7720c24465SJunchao Zhang   b->backend = PETSCSF_BACKEND_KOKKOS;
7827f636e8SJunchao Zhang   #elif defined(PETSC_HAVE_CUDA)
7927f636e8SJunchao Zhang   b->backend = PETSCSF_BACKEND_CUDA;
8059af0bd3SScott Kruger   #elif defined(PETSC_HAVE_HIP)
8159af0bd3SScott Kruger   b->backend = PETSCSF_BACKEND_HIP;
8220c24465SJunchao Zhang   #endif
8371438e86SJunchao Zhang 
8471438e86SJunchao Zhang   #if defined(PETSC_HAVE_NVSHMEM)
8571438e86SJunchao Zhang   b->use_nvshmem     = PETSC_FALSE; /* Default is not to try NVSHMEM */
8671438e86SJunchao Zhang   b->use_nvshmem_get = PETSC_FALSE; /* Default is to use nvshmem_put based protocol */
879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_nvshmem", &b->use_nvshmem, NULL));
889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_nvshmem_get", &b->use_nvshmem_get, NULL));
8971438e86SJunchao Zhang   #endif
9020c24465SJunchao Zhang #endif
9160c22052SBarry Smith   b->vscat.from_n = -1;
9260c22052SBarry Smith   b->vscat.to_n   = -1;
9360c22052SBarry Smith   b->vscat.unit   = MPIU_SCALAR;
9495fce210SBarry Smith   *sf             = b;
9595fce210SBarry Smith   PetscFunctionReturn(0);
9695fce210SBarry Smith }
9795fce210SBarry Smith 
9829046d53SLisandro Dalcin /*@
9995fce210SBarry Smith    PetscSFReset - Reset a star forest so that different sizes or neighbors can be used
10095fce210SBarry Smith 
10195fce210SBarry Smith    Collective
10295fce210SBarry Smith 
1034165533cSJose E. Roman    Input Parameter:
10495fce210SBarry Smith .  sf - star forest
10595fce210SBarry Smith 
10695fce210SBarry Smith    Level: advanced
10795fce210SBarry Smith 
108db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFSetGraph()`, `PetscSFDestroy()`
10995fce210SBarry Smith @*/
110*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFReset(PetscSF sf)
111*d71ae5a4SJacob Faibussowitsch {
11295fce210SBarry Smith   PetscFunctionBegin;
11395fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
114dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf, Reset);
11529046d53SLisandro Dalcin   sf->nroots   = -1;
11629046d53SLisandro Dalcin   sf->nleaves  = -1;
11729046d53SLisandro Dalcin   sf->minleaf  = PETSC_MAX_INT;
11829046d53SLisandro Dalcin   sf->maxleaf  = PETSC_MIN_INT;
11995fce210SBarry Smith   sf->mine     = NULL;
12095fce210SBarry Smith   sf->remote   = NULL;
12129046d53SLisandro Dalcin   sf->graphset = PETSC_FALSE;
1229566063dSJacob Faibussowitsch   PetscCall(PetscFree(sf->mine_alloc));
1239566063dSJacob Faibussowitsch   PetscCall(PetscFree(sf->remote_alloc));
12421c688dcSJed Brown   sf->nranks = -1;
1259566063dSJacob Faibussowitsch   PetscCall(PetscFree4(sf->ranks, sf->roffset, sf->rmine, sf->rremote));
12629046d53SLisandro Dalcin   sf->degreeknown = PETSC_FALSE;
1279566063dSJacob Faibussowitsch   PetscCall(PetscFree(sf->degree));
1289566063dSJacob Faibussowitsch   if (sf->ingroup != MPI_GROUP_NULL) PetscCallMPI(MPI_Group_free(&sf->ingroup));
1299566063dSJacob Faibussowitsch   if (sf->outgroup != MPI_GROUP_NULL) PetscCallMPI(MPI_Group_free(&sf->outgroup));
130013b3241SStefano Zampini   if (sf->multi) sf->multi->multi = NULL;
1319566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf->multi));
1329566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&sf->map));
13371438e86SJunchao Zhang 
13471438e86SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
1359566063dSJacob Faibussowitsch   for (PetscInt i = 0; i < 2; i++) PetscCall(PetscSFFree(sf, PETSC_MEMTYPE_DEVICE, sf->rmine_d[i]));
13671438e86SJunchao Zhang #endif
13771438e86SJunchao Zhang 
13895fce210SBarry Smith   sf->setupcalled = PETSC_FALSE;
13995fce210SBarry Smith   PetscFunctionReturn(0);
14095fce210SBarry Smith }
14195fce210SBarry Smith 
14295fce210SBarry Smith /*@C
14329046d53SLisandro Dalcin    PetscSFSetType - Set the PetscSF communication implementation
14495fce210SBarry Smith 
14595fce210SBarry Smith    Collective on PetscSF
14695fce210SBarry Smith 
14795fce210SBarry Smith    Input Parameters:
14895fce210SBarry Smith +  sf - the PetscSF context
14995fce210SBarry Smith -  type - a known method
15095fce210SBarry Smith 
15195fce210SBarry Smith    Options Database Key:
15295fce210SBarry Smith .  -sf_type <type> - Sets the method; use -help for a list
15370616304SStefano Zampini    of available methods (for instance, window, basic, neighbor)
15495fce210SBarry Smith 
15595fce210SBarry Smith    Notes:
15695fce210SBarry Smith    See "include/petscsf.h" for available methods (for instance)
15795fce210SBarry Smith +    PETSCSFWINDOW - MPI-2/3 one-sided
15895fce210SBarry Smith -    PETSCSFBASIC - basic implementation using MPI-1 two-sided
15995fce210SBarry Smith 
16095fce210SBarry Smith   Level: intermediate
16195fce210SBarry Smith 
162db781477SPatrick Sanan .seealso: `PetscSFType`, `PetscSFCreate()`
16395fce210SBarry Smith @*/
164*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetType(PetscSF sf, PetscSFType type)
165*d71ae5a4SJacob Faibussowitsch {
16695fce210SBarry Smith   PetscBool match;
1675f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(PetscSF);
16895fce210SBarry Smith 
16995fce210SBarry Smith   PetscFunctionBegin;
17095fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
17195fce210SBarry Smith   PetscValidCharPointer(type, 2);
17295fce210SBarry Smith 
1739566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)sf, type, &match));
17495fce210SBarry Smith   if (match) PetscFunctionReturn(0);
17595fce210SBarry Smith 
1769566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(PetscSFList, type, &r));
17728b400f6SJacob Faibussowitsch   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested PetscSF type %s", type);
17829046d53SLisandro Dalcin   /* Destroy the previous PetscSF implementation context */
179dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf, Destroy);
1809566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(sf->ops, sizeof(*sf->ops)));
1819566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)sf, type));
1829566063dSJacob Faibussowitsch   PetscCall((*r)(sf));
18395fce210SBarry Smith   PetscFunctionReturn(0);
18495fce210SBarry Smith }
18595fce210SBarry Smith 
18629046d53SLisandro Dalcin /*@C
18729046d53SLisandro Dalcin   PetscSFGetType - Get the PetscSF communication implementation
18829046d53SLisandro Dalcin 
18929046d53SLisandro Dalcin   Not Collective
19029046d53SLisandro Dalcin 
19129046d53SLisandro Dalcin   Input Parameter:
19229046d53SLisandro Dalcin . sf  - the PetscSF context
19329046d53SLisandro Dalcin 
19429046d53SLisandro Dalcin   Output Parameter:
19529046d53SLisandro Dalcin . type - the PetscSF type name
19629046d53SLisandro Dalcin 
19729046d53SLisandro Dalcin   Level: intermediate
19829046d53SLisandro Dalcin 
199db781477SPatrick Sanan .seealso: `PetscSFSetType()`, `PetscSFCreate()`
20029046d53SLisandro Dalcin @*/
201*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type)
202*d71ae5a4SJacob Faibussowitsch {
20329046d53SLisandro Dalcin   PetscFunctionBegin;
20429046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
20529046d53SLisandro Dalcin   PetscValidPointer(type, 2);
20629046d53SLisandro Dalcin   *type = ((PetscObject)sf)->type_name;
20729046d53SLisandro Dalcin   PetscFunctionReturn(0);
20829046d53SLisandro Dalcin }
20929046d53SLisandro Dalcin 
2101fb7b255SJunchao Zhang /*@C
21195fce210SBarry Smith    PetscSFDestroy - destroy star forest
21295fce210SBarry Smith 
21395fce210SBarry Smith    Collective
21495fce210SBarry Smith 
2154165533cSJose E. Roman    Input Parameter:
21695fce210SBarry Smith .  sf - address of star forest
21795fce210SBarry Smith 
21895fce210SBarry Smith    Level: intermediate
21995fce210SBarry Smith 
220db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFReset()`
22195fce210SBarry Smith @*/
222*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFDestroy(PetscSF *sf)
223*d71ae5a4SJacob Faibussowitsch {
22495fce210SBarry Smith   PetscFunctionBegin;
22595fce210SBarry Smith   if (!*sf) PetscFunctionReturn(0);
22695fce210SBarry Smith   PetscValidHeaderSpecific((*sf), PETSCSF_CLASSID, 1);
2279371c9d4SSatish Balay   if (--((PetscObject)(*sf))->refct > 0) {
2289371c9d4SSatish Balay     *sf = NULL;
2299371c9d4SSatish Balay     PetscFunctionReturn(0);
2309371c9d4SSatish Balay   }
2319566063dSJacob Faibussowitsch   PetscCall(PetscSFReset(*sf));
232dbbe0bcdSBarry Smith   PetscTryTypeMethod((*sf), Destroy);
2339566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&(*sf)->vscat.lsf));
2349566063dSJacob Faibussowitsch   if ((*sf)->vscat.bs > 1) PetscCallMPI(MPI_Type_free(&(*sf)->vscat.unit));
2359566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(sf));
23695fce210SBarry Smith   PetscFunctionReturn(0);
23795fce210SBarry Smith }
23895fce210SBarry Smith 
239*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFCheckGraphValid_Private(PetscSF sf)
240*d71ae5a4SJacob Faibussowitsch {
241c4e6a40aSLawrence Mitchell   PetscInt           i, nleaves;
242c4e6a40aSLawrence Mitchell   PetscMPIInt        size;
243c4e6a40aSLawrence Mitchell   const PetscInt    *ilocal;
244c4e6a40aSLawrence Mitchell   const PetscSFNode *iremote;
245c4e6a40aSLawrence Mitchell 
246c4e6a40aSLawrence Mitchell   PetscFunctionBegin;
24776bd3646SJed Brown   if (!sf->graphset || !PetscDefined(USE_DEBUG)) PetscFunctionReturn(0);
2489566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, &ilocal, &iremote));
2499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf), &size));
250c4e6a40aSLawrence Mitchell   for (i = 0; i < nleaves; i++) {
251c4e6a40aSLawrence Mitchell     const PetscInt rank   = iremote[i].rank;
252c4e6a40aSLawrence Mitchell     const PetscInt remote = iremote[i].index;
253c4e6a40aSLawrence Mitchell     const PetscInt leaf   = ilocal ? ilocal[i] : i;
254c9cc58a2SBarry Smith     PetscCheck(rank >= 0 && rank < size, 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);
25508401ef6SPierre Jolivet     PetscCheck(remote >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Provided index (%" PetscInt_FMT ") for remote %" PetscInt_FMT " is invalid, should be >= 0", remote, i);
25608401ef6SPierre Jolivet     PetscCheck(leaf >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Provided location (%" PetscInt_FMT ") for leaf %" PetscInt_FMT " is invalid, should be >= 0", leaf, i);
257c4e6a40aSLawrence Mitchell   }
258c4e6a40aSLawrence Mitchell   PetscFunctionReturn(0);
259c4e6a40aSLawrence Mitchell }
260c4e6a40aSLawrence Mitchell 
26195fce210SBarry Smith /*@
26295fce210SBarry Smith    PetscSFSetUp - set up communication structures
26395fce210SBarry Smith 
26495fce210SBarry Smith    Collective
26595fce210SBarry Smith 
2664165533cSJose E. Roman    Input Parameter:
26795fce210SBarry Smith .  sf - star forest communication object
26895fce210SBarry Smith 
26995fce210SBarry Smith    Level: beginner
27095fce210SBarry Smith 
271db781477SPatrick Sanan .seealso: `PetscSFSetFromOptions()`, `PetscSFSetType()`
27295fce210SBarry Smith @*/
273*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetUp(PetscSF sf)
274*d71ae5a4SJacob Faibussowitsch {
27595fce210SBarry Smith   PetscFunctionBegin;
27629046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
27729046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf, 1);
27895fce210SBarry Smith   if (sf->setupcalled) PetscFunctionReturn(0);
2799566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_SetUp, sf, 0, 0, 0));
2809566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckGraphValid_Private(sf));
2819566063dSJacob Faibussowitsch   if (!((PetscObject)sf)->type_name) PetscCall(PetscSFSetType(sf, PETSCSFBASIC)); /* Zero all sf->ops */
282dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf, SetUp);
28320c24465SJunchao Zhang #if defined(PETSC_HAVE_CUDA)
28420c24465SJunchao Zhang   if (sf->backend == PETSCSF_BACKEND_CUDA) {
28571438e86SJunchao Zhang     sf->ops->Malloc = PetscSFMalloc_CUDA;
28671438e86SJunchao Zhang     sf->ops->Free   = PetscSFFree_CUDA;
28720c24465SJunchao Zhang   }
28820c24465SJunchao Zhang #endif
28959af0bd3SScott Kruger #if defined(PETSC_HAVE_HIP)
29059af0bd3SScott Kruger   if (sf->backend == PETSCSF_BACKEND_HIP) {
29159af0bd3SScott Kruger     sf->ops->Malloc = PetscSFMalloc_HIP;
29259af0bd3SScott Kruger     sf->ops->Free   = PetscSFFree_HIP;
29359af0bd3SScott Kruger   }
29459af0bd3SScott Kruger #endif
29520c24465SJunchao Zhang 
29659af0bd3SScott Kruger #
29720c24465SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
29820c24465SJunchao Zhang   if (sf->backend == PETSCSF_BACKEND_KOKKOS) {
29920c24465SJunchao Zhang     sf->ops->Malloc = PetscSFMalloc_Kokkos;
30020c24465SJunchao Zhang     sf->ops->Free   = PetscSFFree_Kokkos;
30120c24465SJunchao Zhang   }
30220c24465SJunchao Zhang #endif
3039566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_SetUp, sf, 0, 0, 0));
30495fce210SBarry Smith   sf->setupcalled = PETSC_TRUE;
30595fce210SBarry Smith   PetscFunctionReturn(0);
30695fce210SBarry Smith }
30795fce210SBarry Smith 
3088af6ec1cSBarry Smith /*@
30995fce210SBarry Smith    PetscSFSetFromOptions - set PetscSF options using the options database
31095fce210SBarry Smith 
31195fce210SBarry Smith    Logically Collective
31295fce210SBarry Smith 
3134165533cSJose E. Roman    Input Parameter:
31495fce210SBarry Smith .  sf - star forest
31595fce210SBarry Smith 
31695fce210SBarry Smith    Options Database Keys:
31760263706SJed Brown +  -sf_type               - implementation type, see PetscSFSetType()
31851ccb202SJunchao Zhang .  -sf_rank_order         - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
319b85e67b7SJunchao Zhang .  -sf_use_default_stream - Assume callers of SF computed the input root/leafdata with the default cuda stream. SF will also
320c2a741eeSJunchao Zhang                             use the default stream to process data. Therefore, no stream synchronization is needed between SF and its caller (default: true).
321c06a8e02SRichard Tran Mills                             If true, this option only works with -use_gpu_aware_mpi 1.
32220c24465SJunchao 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).
323c06a8e02SRichard Tran Mills                                If true, this option only works with -use_gpu_aware_mpi 1.
32495fce210SBarry Smith 
32559af0bd3SScott Kruger -  -sf_backend cuda | hip | kokkos -Select the device backend SF uses. Currently SF has these backends: cuda, hip and Kokkos.
32659af0bd3SScott Kruger                               On CUDA (HIP) devices, one can choose cuda (hip) or kokkos with the default being kokkos. On other devices,
32720c24465SJunchao Zhang                               the only available is kokkos.
32820c24465SJunchao Zhang 
32995fce210SBarry Smith    Level: intermediate
33095fce210SBarry Smith @*/
331*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
332*d71ae5a4SJacob Faibussowitsch {
33395fce210SBarry Smith   PetscSFType deft;
33495fce210SBarry Smith   char        type[256];
33595fce210SBarry Smith   PetscBool   flg;
33695fce210SBarry Smith 
33795fce210SBarry Smith   PetscFunctionBegin;
33895fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
339d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)sf);
34095fce210SBarry Smith   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
3419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-sf_type", "PetscSF implementation type", "PetscSFSetType", PetscSFList, deft, type, sizeof(type), &flg));
3429566063dSJacob Faibussowitsch   PetscCall(PetscSFSetType(sf, flg ? type : deft));
3439566063dSJacob Faibussowitsch   PetscCall(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));
3447fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
34520c24465SJunchao Zhang   {
34620c24465SJunchao Zhang     char      backendstr[32] = {0};
34759af0bd3SScott Kruger     PetscBool isCuda = PETSC_FALSE, isHip = PETSC_FALSE, isKokkos = PETSC_FALSE, set;
34820c24465SJunchao Zhang     /* Change the defaults set in PetscSFCreate() with command line options */
3499566063dSJacob Faibussowitsch     PetscCall(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));
3509566063dSJacob Faibussowitsch     PetscCall(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));
3519566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-sf_backend", "Select the device backend SF uses", "PetscSFSetFromOptions", NULL, backendstr, sizeof(backendstr), &set));
3529566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp("cuda", backendstr, &isCuda));
3539566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp("kokkos", backendstr, &isKokkos));
3549566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp("hip", backendstr, &isHip));
35559af0bd3SScott Kruger   #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
35620c24465SJunchao Zhang     if (isCuda) sf->backend = PETSCSF_BACKEND_CUDA;
35720c24465SJunchao Zhang     else if (isKokkos) sf->backend = PETSCSF_BACKEND_KOKKOS;
35859af0bd3SScott Kruger     else if (isHip) sf->backend = PETSCSF_BACKEND_HIP;
35928b400f6SJacob Faibussowitsch     else PetscCheck(!set, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "-sf_backend %s is not supported. You may choose cuda, hip or kokkos (if installed)", backendstr);
36020c24465SJunchao Zhang   #elif defined(PETSC_HAVE_KOKKOS)
36108401ef6SPierre Jolivet     PetscCheck(!set || isKokkos, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "-sf_backend %s is not supported. You can only choose kokkos", backendstr);
36220c24465SJunchao Zhang   #endif
36320c24465SJunchao Zhang   }
364c2a741eeSJunchao Zhang #endif
365dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf, SetFromOptions, PetscOptionsObject);
366d0609cedSBarry Smith   PetscOptionsEnd();
36795fce210SBarry Smith   PetscFunctionReturn(0);
36895fce210SBarry Smith }
36995fce210SBarry Smith 
37029046d53SLisandro Dalcin /*@
37195fce210SBarry Smith    PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order
37295fce210SBarry Smith 
37395fce210SBarry Smith    Logically Collective
37495fce210SBarry Smith 
3754165533cSJose E. Roman    Input Parameters:
37695fce210SBarry Smith +  sf - star forest
37795fce210SBarry Smith -  flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic)
37895fce210SBarry Smith 
37995fce210SBarry Smith    Level: advanced
38095fce210SBarry Smith 
381db781477SPatrick Sanan .seealso: `PetscSFGatherBegin()`, `PetscSFScatterBegin()`
38295fce210SBarry Smith @*/
383*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetRankOrder(PetscSF sf, PetscBool flg)
384*d71ae5a4SJacob Faibussowitsch {
38595fce210SBarry Smith   PetscFunctionBegin;
38695fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
38795fce210SBarry Smith   PetscValidLogicalCollectiveBool(sf, flg, 2);
38828b400f6SJacob Faibussowitsch   PetscCheck(!sf->multi, PetscObjectComm((PetscObject)sf), PETSC_ERR_ARG_WRONGSTATE, "Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
38995fce210SBarry Smith   sf->rankorder = flg;
39095fce210SBarry Smith   PetscFunctionReturn(0);
39195fce210SBarry Smith }
39295fce210SBarry Smith 
3938dbb0df6SBarry Smith /*@C
39495fce210SBarry Smith    PetscSFSetGraph - Set a parallel star forest
39595fce210SBarry Smith 
39695fce210SBarry Smith    Collective
39795fce210SBarry Smith 
3984165533cSJose E. Roman    Input Parameters:
39995fce210SBarry Smith +  sf - star forest
40095fce210SBarry Smith .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
40195fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
402c4e6a40aSLawrence Mitchell .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage (locations must be >= 0, enforced
403c4e6a40aSLawrence Mitchell during setup in debug mode)
40495fce210SBarry Smith .  localmode - copy mode for ilocal
405c4e6a40aSLawrence Mitchell .  iremote - remote locations of root vertices for each leaf on the current process (locations must be >= 0, enforced
406c4e6a40aSLawrence Mitchell during setup in debug mode)
40795fce210SBarry Smith -  remotemode - copy mode for iremote
40895fce210SBarry Smith 
40995fce210SBarry Smith    Level: intermediate
41095fce210SBarry Smith 
41195452b02SPatrick Sanan    Notes:
412db2b9530SVaclav Hapla    Leaf indices in ilocal must be unique, otherwise an error occurs.
41338ab3f8aSBarry Smith 
414db2b9530SVaclav Hapla    Input arrays ilocal and iremote follow the PetscCopyMode semantics.
415db2b9530SVaclav Hapla    In particular, if localmode/remotemode is PETSC_OWN_POINTER or PETSC_USE_POINTER,
416db2b9530SVaclav Hapla    PETSc might modify the respective array;
417db2b9530SVaclav Hapla    if PETSC_USE_POINTER, the user must delete the array after PetscSFDestroy().
418e97d3de3SVaclav Hapla    Only if PETSC_COPY_VALUES is used, the respective array is guaranteed to stay intact and a const array can be passed (but a cast to non-const is needed).
419db2b9530SVaclav Hapla 
420db2b9530SVaclav Hapla    Fortran Notes:
421db2b9530SVaclav Hapla    In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode.
422c4e6a40aSLawrence Mitchell 
42372bf8598SVaclav Hapla    Developer Notes:
424db2b9530SVaclav Hapla    We sort leaves to check for duplicates and contiguousness and to find minleaf/maxleaf.
425db2b9530SVaclav Hapla    This also allows to compare leaf sets of two SFs easily.
42672bf8598SVaclav Hapla 
427db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFView()`, `PetscSFGetGraph()`
42895fce210SBarry Smith @*/
429*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraph(PetscSF sf, PetscInt nroots, PetscInt nleaves, PetscInt *ilocal, PetscCopyMode localmode, PetscSFNode *iremote, PetscCopyMode remotemode)
430*d71ae5a4SJacob Faibussowitsch {
431db2b9530SVaclav Hapla   PetscBool unique, contiguous;
43295fce210SBarry Smith 
43395fce210SBarry Smith   PetscFunctionBegin;
43495fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
43529046d53SLisandro Dalcin   if (nleaves > 0 && ilocal) PetscValidIntPointer(ilocal, 4);
43629046d53SLisandro Dalcin   if (nleaves > 0) PetscValidPointer(iremote, 6);
43708401ef6SPierre Jolivet   PetscCheck(nroots >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nroots %" PetscInt_FMT ", cannot be negative", nroots);
43808401ef6SPierre Jolivet   PetscCheck(nleaves >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nleaves %" PetscInt_FMT ", cannot be negative", nleaves);
4398da24d32SBarry Smith   /* enums may be handled as unsigned by some compilers, NVHPC for example, the int cast
4408da24d32SBarry Smith    * below is to prevent NVHPC from warning about meaningless comparison of unsigned with zero */
4418da24d32SBarry Smith   PetscCheck((int)localmode >= PETSC_COPY_VALUES && localmode <= PETSC_USE_POINTER, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wrong localmode %d", localmode);
4428da24d32SBarry Smith   PetscCheck((int)remotemode >= PETSC_COPY_VALUES && remotemode <= PETSC_USE_POINTER, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wrong remotemode %d", remotemode);
44329046d53SLisandro Dalcin 
4442a67d2daSStefano Zampini   if (sf->nroots >= 0) { /* Reset only if graph already set */
4459566063dSJacob Faibussowitsch     PetscCall(PetscSFReset(sf));
4462a67d2daSStefano Zampini   }
4472a67d2daSStefano Zampini 
4489566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_SetGraph, sf, 0, 0, 0));
44929046d53SLisandro Dalcin 
45095fce210SBarry Smith   sf->nroots  = nroots;
45195fce210SBarry Smith   sf->nleaves = nleaves;
45229046d53SLisandro Dalcin 
453db2b9530SVaclav Hapla   if (localmode == PETSC_COPY_VALUES && ilocal) {
454db2b9530SVaclav Hapla     PetscInt *tlocal = NULL;
455db2b9530SVaclav Hapla 
4569566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &tlocal));
4579566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tlocal, ilocal, nleaves));
458db2b9530SVaclav Hapla     ilocal = tlocal;
459db2b9530SVaclav Hapla   }
460db2b9530SVaclav Hapla   if (remotemode == PETSC_COPY_VALUES) {
461db2b9530SVaclav Hapla     PetscSFNode *tremote = NULL;
462db2b9530SVaclav Hapla 
4639566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &tremote));
4649566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tremote, iremote, nleaves));
465db2b9530SVaclav Hapla     iremote = tremote;
466db2b9530SVaclav Hapla   }
467db2b9530SVaclav Hapla 
46829046d53SLisandro Dalcin   if (nleaves && ilocal) {
469db2b9530SVaclav Hapla     PetscSFNode work;
470db2b9530SVaclav Hapla 
4719566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithDataArray(nleaves, ilocal, iremote, sizeof(PetscSFNode), &work));
4729566063dSJacob Faibussowitsch     PetscCall(PetscSortedCheckDupsInt(nleaves, ilocal, &unique));
473db2b9530SVaclav Hapla     unique = PetscNot(unique);
474db2b9530SVaclav Hapla     PetscCheck(sf->allow_multi_leaves || unique, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input ilocal has duplicate entries which is not allowed for this PetscSF");
475db2b9530SVaclav Hapla     sf->minleaf = ilocal[0];
476db2b9530SVaclav Hapla     sf->maxleaf = ilocal[nleaves - 1];
477db2b9530SVaclav Hapla     contiguous  = (PetscBool)(unique && ilocal[0] == 0 && ilocal[nleaves - 1] == nleaves - 1);
47829046d53SLisandro Dalcin   } else {
47929046d53SLisandro Dalcin     sf->minleaf = 0;
48029046d53SLisandro Dalcin     sf->maxleaf = nleaves - 1;
481db2b9530SVaclav Hapla     unique      = PETSC_TRUE;
482db2b9530SVaclav Hapla     contiguous  = PETSC_TRUE;
48329046d53SLisandro Dalcin   }
48429046d53SLisandro Dalcin 
485db2b9530SVaclav Hapla   if (contiguous) {
486db2b9530SVaclav Hapla     if (localmode == PETSC_USE_POINTER) {
487db2b9530SVaclav Hapla       ilocal = NULL;
488db2b9530SVaclav Hapla     } else {
4899566063dSJacob Faibussowitsch       PetscCall(PetscFree(ilocal));
490db2b9530SVaclav Hapla     }
491db2b9530SVaclav Hapla   }
492db2b9530SVaclav Hapla   sf->mine = ilocal;
493db2b9530SVaclav Hapla   if (localmode == PETSC_USE_POINTER) {
49429046d53SLisandro Dalcin     sf->mine_alloc = NULL;
495db2b9530SVaclav Hapla   } else {
496db2b9530SVaclav Hapla     sf->mine_alloc = ilocal;
49795fce210SBarry Smith   }
498db2b9530SVaclav Hapla   sf->remote = iremote;
499db2b9530SVaclav Hapla   if (remotemode == PETSC_USE_POINTER) {
50029046d53SLisandro Dalcin     sf->remote_alloc = NULL;
501db2b9530SVaclav Hapla   } else {
502db2b9530SVaclav Hapla     sf->remote_alloc = iremote;
50395fce210SBarry Smith   }
5049566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_SetGraph, sf, 0, 0, 0));
50529046d53SLisandro Dalcin   sf->graphset = PETSC_TRUE;
50695fce210SBarry Smith   PetscFunctionReturn(0);
50795fce210SBarry Smith }
50895fce210SBarry Smith 
50929046d53SLisandro Dalcin /*@
510dd5b3ca6SJunchao Zhang   PetscSFSetGraphWithPattern - Sets the graph of an SF with a specific pattern
511dd5b3ca6SJunchao Zhang 
512dd5b3ca6SJunchao Zhang   Collective
513dd5b3ca6SJunchao Zhang 
514dd5b3ca6SJunchao Zhang   Input Parameters:
515dd5b3ca6SJunchao Zhang + sf      - The PetscSF
516dd5b3ca6SJunchao Zhang . map     - Layout of roots over all processes (insignificant when pattern is PETSCSF_PATTERN_ALLTOALL)
517dd5b3ca6SJunchao Zhang - pattern - One of PETSCSF_PATTERN_ALLGATHER, PETSCSF_PATTERN_GATHER, PETSCSF_PATTERN_ALLTOALL
518dd5b3ca6SJunchao Zhang 
519dd5b3ca6SJunchao Zhang   Notes:
520dd5b3ca6SJunchao Zhang   It is easier to explain PetscSFPattern using vectors. Suppose we have an MPI vector x and its layout is map.
521dd5b3ca6SJunchao Zhang   n and N are local and global sizes of x respectively.
522dd5b3ca6SJunchao Zhang 
523dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLGATHER, the routine creates a graph that if one does Bcast on it, it will copy x to
524dd5b3ca6SJunchao Zhang   sequential vectors y on all ranks.
525dd5b3ca6SJunchao Zhang 
526dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_GATHER, the routine creates a graph that if one does Bcast on it, it will copy x to a
527dd5b3ca6SJunchao Zhang   sequential vector y on rank 0.
528dd5b3ca6SJunchao Zhang 
529dd5b3ca6SJunchao Zhang   In above cases, entries of x are roots and entries of y are leaves.
530dd5b3ca6SJunchao Zhang 
531dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLTOALL, map is insignificant. Suppose NP is size of sf's communicator. The routine
532dd5b3ca6SJunchao Zhang   creates a graph that every rank has NP leaves and NP roots. On rank i, its leaf j is connected to root i
533dd5b3ca6SJunchao 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
534dd5b3ca6SJunchao Zhang   not mean one can not send multiple items. One just needs to create a new MPI datatype for the mulptiple data
535dd5b3ca6SJunchao Zhang   items with MPI_Type_contiguous() and use that as the <unit> argument in SF routines.
536dd5b3ca6SJunchao Zhang 
537dd5b3ca6SJunchao Zhang   In this case, roots and leaves are symmetric.
538dd5b3ca6SJunchao Zhang 
539dd5b3ca6SJunchao Zhang   Level: intermediate
540dd5b3ca6SJunchao Zhang  @*/
541*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphWithPattern(PetscSF sf, PetscLayout map, PetscSFPattern pattern)
542*d71ae5a4SJacob Faibussowitsch {
543dd5b3ca6SJunchao Zhang   MPI_Comm    comm;
544dd5b3ca6SJunchao Zhang   PetscInt    n, N, res[2];
545dd5b3ca6SJunchao Zhang   PetscMPIInt rank, size;
546dd5b3ca6SJunchao Zhang   PetscSFType type;
547dd5b3ca6SJunchao Zhang 
548dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
5492abc8c78SJacob Faibussowitsch   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
5502abc8c78SJacob Faibussowitsch   if (pattern != PETSCSF_PATTERN_ALLTOALL) PetscValidPointer(map, 2);
5519566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
5522c71b3e2SJacob Faibussowitsch   PetscCheck(pattern >= PETSCSF_PATTERN_ALLGATHER && pattern <= PETSCSF_PATTERN_ALLTOALL, comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PetscSFPattern %d", pattern);
5539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
555dd5b3ca6SJunchao Zhang 
556dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLTOALL) {
557dd5b3ca6SJunchao Zhang     type = PETSCSFALLTOALL;
5589566063dSJacob Faibussowitsch     PetscCall(PetscLayoutCreate(comm, &sf->map));
5599566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetLocalSize(sf->map, size));
5609566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetSize(sf->map, ((PetscInt)size) * size));
5619566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetUp(sf->map));
562dd5b3ca6SJunchao Zhang   } else {
5639566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetLocalSize(map, &n));
5649566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetSize(map, &N));
565dd5b3ca6SJunchao Zhang     res[0] = n;
566dd5b3ca6SJunchao Zhang     res[1] = -n;
567dd5b3ca6SJunchao Zhang     /* Check if n are same over all ranks so that we can optimize it */
5681c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, res, 2, MPIU_INT, MPI_MAX, comm));
569dd5b3ca6SJunchao Zhang     if (res[0] == -res[1]) { /* same n */
570dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHER : PETSCSFGATHER;
571dd5b3ca6SJunchao Zhang     } else {
572dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHERV : PETSCSFGATHERV;
573dd5b3ca6SJunchao Zhang     }
5749566063dSJacob Faibussowitsch     PetscCall(PetscLayoutReference(map, &sf->map));
575dd5b3ca6SJunchao Zhang   }
5769566063dSJacob Faibussowitsch   PetscCall(PetscSFSetType(sf, type));
577dd5b3ca6SJunchao Zhang 
578dd5b3ca6SJunchao Zhang   sf->pattern = pattern;
579dd5b3ca6SJunchao Zhang   sf->mine    = NULL; /* Contiguous */
580dd5b3ca6SJunchao Zhang 
581dd5b3ca6SJunchao Zhang   /* Set nleaves, nroots here in case user calls PetscSFGetGraph, which is legal to call even before PetscSFSetUp is called.
582dd5b3ca6SJunchao Zhang      Also set other easy stuff.
583dd5b3ca6SJunchao Zhang    */
584dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLGATHER) {
585dd5b3ca6SJunchao Zhang     sf->nleaves = N;
586dd5b3ca6SJunchao Zhang     sf->nroots  = n;
587dd5b3ca6SJunchao Zhang     sf->nranks  = size;
588dd5b3ca6SJunchao Zhang     sf->minleaf = 0;
589dd5b3ca6SJunchao Zhang     sf->maxleaf = N - 1;
590dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_GATHER) {
591dd5b3ca6SJunchao Zhang     sf->nleaves = rank ? 0 : N;
592dd5b3ca6SJunchao Zhang     sf->nroots  = n;
593dd5b3ca6SJunchao Zhang     sf->nranks  = rank ? 0 : size;
594dd5b3ca6SJunchao Zhang     sf->minleaf = 0;
595dd5b3ca6SJunchao Zhang     sf->maxleaf = rank ? -1 : N - 1;
596dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_ALLTOALL) {
597dd5b3ca6SJunchao Zhang     sf->nleaves = size;
598dd5b3ca6SJunchao Zhang     sf->nroots  = size;
599dd5b3ca6SJunchao Zhang     sf->nranks  = size;
600dd5b3ca6SJunchao Zhang     sf->minleaf = 0;
601dd5b3ca6SJunchao Zhang     sf->maxleaf = size - 1;
602dd5b3ca6SJunchao Zhang   }
603dd5b3ca6SJunchao Zhang   sf->ndranks  = 0; /* We do not need to separate out distinguished ranks for patterned graphs to improve communication performance */
604dd5b3ca6SJunchao Zhang   sf->graphset = PETSC_TRUE;
605dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
606dd5b3ca6SJunchao Zhang }
607dd5b3ca6SJunchao Zhang 
608dd5b3ca6SJunchao Zhang /*@
60995fce210SBarry Smith    PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map
61095fce210SBarry Smith 
61195fce210SBarry Smith    Collective
61295fce210SBarry Smith 
6134165533cSJose E. Roman    Input Parameter:
61495fce210SBarry Smith .  sf - star forest to invert
61595fce210SBarry Smith 
6164165533cSJose E. Roman    Output Parameter:
61795fce210SBarry Smith .  isf - inverse of sf
6184165533cSJose E. Roman 
61995fce210SBarry Smith    Level: advanced
62095fce210SBarry Smith 
62195fce210SBarry Smith    Notes:
62295fce210SBarry Smith    All roots must have degree 1.
62395fce210SBarry Smith 
62495fce210SBarry Smith    The local space may be a permutation, but cannot be sparse.
62595fce210SBarry Smith 
626db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`
62795fce210SBarry Smith @*/
628*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateInverseSF(PetscSF sf, PetscSF *isf)
629*d71ae5a4SJacob Faibussowitsch {
63095fce210SBarry Smith   PetscMPIInt     rank;
63195fce210SBarry Smith   PetscInt        i, nroots, nleaves, maxlocal, count, *newilocal;
63295fce210SBarry Smith   const PetscInt *ilocal;
63395fce210SBarry Smith   PetscSFNode    *roots, *leaves;
63495fce210SBarry Smith 
63595fce210SBarry Smith   PetscFunctionBegin;
63629046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
63729046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf, 1);
63829046d53SLisandro Dalcin   PetscValidPointer(isf, 2);
63929046d53SLisandro Dalcin 
6409566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL));
64129046d53SLisandro Dalcin   maxlocal = sf->maxleaf + 1; /* TODO: We should use PetscSFGetLeafRange() */
64229046d53SLisandro Dalcin 
6439566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
6449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nroots, &roots, maxlocal, &leaves));
645ae9aee6dSMatthew G. Knepley   for (i = 0; i < maxlocal; i++) {
64695fce210SBarry Smith     leaves[i].rank  = rank;
64795fce210SBarry Smith     leaves[i].index = i;
64895fce210SBarry Smith   }
64995fce210SBarry Smith   for (i = 0; i < nroots; i++) {
65095fce210SBarry Smith     roots[i].rank  = -1;
65195fce210SBarry Smith     roots[i].index = -1;
65295fce210SBarry Smith   }
6539566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf, MPIU_2INT, leaves, roots, MPI_REPLACE));
6549566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sf, MPIU_2INT, leaves, roots, MPI_REPLACE));
65595fce210SBarry Smith 
65695fce210SBarry Smith   /* Check whether our leaves are sparse */
6579371c9d4SSatish Balay   for (i = 0, count = 0; i < nroots; i++)
6589371c9d4SSatish Balay     if (roots[i].rank >= 0) count++;
65995fce210SBarry Smith   if (count == nroots) newilocal = NULL;
6609371c9d4SSatish Balay   else { /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */ PetscCall(PetscMalloc1(count, &newilocal));
66195fce210SBarry Smith     for (i = 0, count = 0; i < nroots; i++) {
66295fce210SBarry Smith       if (roots[i].rank >= 0) {
66395fce210SBarry Smith         newilocal[count]   = i;
66495fce210SBarry Smith         roots[count].rank  = roots[i].rank;
66595fce210SBarry Smith         roots[count].index = roots[i].index;
66695fce210SBarry Smith         count++;
66795fce210SBarry Smith       }
66895fce210SBarry Smith     }
66995fce210SBarry Smith   }
67095fce210SBarry Smith 
6719566063dSJacob Faibussowitsch   PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_CONFONLY, isf));
6729566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*isf, maxlocal, count, newilocal, PETSC_OWN_POINTER, roots, PETSC_COPY_VALUES));
6739566063dSJacob Faibussowitsch   PetscCall(PetscFree2(roots, leaves));
67495fce210SBarry Smith   PetscFunctionReturn(0);
67595fce210SBarry Smith }
67695fce210SBarry Smith 
67795fce210SBarry Smith /*@
67895fce210SBarry Smith    PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph
67995fce210SBarry Smith 
68095fce210SBarry Smith    Collective
68195fce210SBarry Smith 
6824165533cSJose E. Roman    Input Parameters:
68395fce210SBarry Smith +  sf - communication object to duplicate
68495fce210SBarry Smith -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)
68595fce210SBarry Smith 
6864165533cSJose E. Roman    Output Parameter:
68795fce210SBarry Smith .  newsf - new communication object
68895fce210SBarry Smith 
68995fce210SBarry Smith    Level: beginner
69095fce210SBarry Smith 
691db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFSetType()`, `PetscSFSetGraph()`
69295fce210SBarry Smith @*/
693*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFDuplicate(PetscSF sf, PetscSFDuplicateOption opt, PetscSF *newsf)
694*d71ae5a4SJacob Faibussowitsch {
69529046d53SLisandro Dalcin   PetscSFType  type;
69697929ea7SJunchao Zhang   MPI_Datatype dtype = MPIU_SCALAR;
69795fce210SBarry Smith 
69895fce210SBarry Smith   PetscFunctionBegin;
69929046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
70029046d53SLisandro Dalcin   PetscValidLogicalCollectiveEnum(sf, opt, 2);
70129046d53SLisandro Dalcin   PetscValidPointer(newsf, 3);
7029566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sf), newsf));
7039566063dSJacob Faibussowitsch   PetscCall(PetscSFGetType(sf, &type));
7049566063dSJacob Faibussowitsch   if (type) PetscCall(PetscSFSetType(*newsf, type));
70554a24607SJunchao Zhang   (*newsf)->allow_multi_leaves = sf->allow_multi_leaves; /* Dup this flag ealier since PetscSFSetGraph() below checks on this flag */
70695fce210SBarry Smith   if (opt == PETSCSF_DUPLICATE_GRAPH) {
707dd5b3ca6SJunchao Zhang     PetscSFCheckGraphSet(sf, 1);
708dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
70995fce210SBarry Smith       PetscInt           nroots, nleaves;
71095fce210SBarry Smith       const PetscInt    *ilocal;
71195fce210SBarry Smith       const PetscSFNode *iremote;
7129566063dSJacob Faibussowitsch       PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, &iremote));
7139566063dSJacob Faibussowitsch       PetscCall(PetscSFSetGraph(*newsf, nroots, nleaves, (PetscInt *)ilocal, PETSC_COPY_VALUES, (PetscSFNode *)iremote, PETSC_COPY_VALUES));
714dd5b3ca6SJunchao Zhang     } else {
7159566063dSJacob Faibussowitsch       PetscCall(PetscSFSetGraphWithPattern(*newsf, sf->map, sf->pattern));
716dd5b3ca6SJunchao Zhang     }
71795fce210SBarry Smith   }
71897929ea7SJunchao Zhang   /* Since oldtype is committed, so is newtype, according to MPI */
7199566063dSJacob Faibussowitsch   if (sf->vscat.bs > 1) PetscCallMPI(MPI_Type_dup(sf->vscat.unit, &dtype));
72097929ea7SJunchao Zhang   (*newsf)->vscat.bs     = sf->vscat.bs;
72197929ea7SJunchao Zhang   (*newsf)->vscat.unit   = dtype;
72297929ea7SJunchao Zhang   (*newsf)->vscat.to_n   = sf->vscat.to_n;
72397929ea7SJunchao Zhang   (*newsf)->vscat.from_n = sf->vscat.from_n;
72497929ea7SJunchao Zhang   /* Do not copy lsf. Build it on demand since it is rarely used */
72597929ea7SJunchao Zhang 
72620c24465SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
72720c24465SJunchao Zhang   (*newsf)->backend              = sf->backend;
72871438e86SJunchao Zhang   (*newsf)->unknown_input_stream = sf->unknown_input_stream;
72920c24465SJunchao Zhang   (*newsf)->use_gpu_aware_mpi    = sf->use_gpu_aware_mpi;
73020c24465SJunchao Zhang   (*newsf)->use_stream_aware_mpi = sf->use_stream_aware_mpi;
73120c24465SJunchao Zhang #endif
732dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf, Duplicate, opt, *newsf);
73320c24465SJunchao Zhang   /* Don't do PetscSFSetUp() since the new sf's graph might have not been set. */
73495fce210SBarry Smith   PetscFunctionReturn(0);
73595fce210SBarry Smith }
73695fce210SBarry Smith 
73795fce210SBarry Smith /*@C
73895fce210SBarry Smith    PetscSFGetGraph - Get the graph specifying a parallel star forest
73995fce210SBarry Smith 
74095fce210SBarry Smith    Not Collective
74195fce210SBarry Smith 
7424165533cSJose E. Roman    Input Parameter:
74395fce210SBarry Smith .  sf - star forest
74495fce210SBarry Smith 
7454165533cSJose E. Roman    Output Parameters:
74695fce210SBarry Smith +  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
74795fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
748bc6585dcSJunchao Zhang .  ilocal - locations of leaves in leafdata buffers (if returned value is NULL, it means leaves are in contiguous storage)
74995fce210SBarry Smith -  iremote - remote locations of root vertices for each leaf on the current process
75095fce210SBarry Smith 
751373e0d91SLisandro Dalcin    Notes:
752373e0d91SLisandro Dalcin      We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set yet
753373e0d91SLisandro Dalcin 
754db2b9530SVaclav Hapla      The returned ilocal/iremote might contain values in different order than the input ones in PetscSFSetGraph(),
755db2b9530SVaclav Hapla      see its manpage for details.
756db2b9530SVaclav Hapla 
7578dbb0df6SBarry Smith    Fortran Notes:
7588dbb0df6SBarry Smith      The returned iremote array is a copy and must be deallocated after use. Consequently, if you
7598dbb0df6SBarry Smith      want to update the graph, you must call PetscSFSetGraph() after modifying the iremote array.
7608dbb0df6SBarry Smith 
7618dbb0df6SBarry Smith      To check for a NULL ilocal use
7628dbb0df6SBarry Smith $      if (loc(ilocal) == loc(PETSC_NULL_INTEGER)) then
763ca797d7aSLawrence Mitchell 
76495fce210SBarry Smith    Level: intermediate
76595fce210SBarry Smith 
766db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`
76795fce210SBarry Smith @*/
768*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGraph(PetscSF sf, PetscInt *nroots, PetscInt *nleaves, const PetscInt **ilocal, const PetscSFNode **iremote)
769*d71ae5a4SJacob Faibussowitsch {
77095fce210SBarry Smith   PetscFunctionBegin;
77195fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
772b8dee149SJunchao Zhang   if (sf->ops->GetGraph) {
7739566063dSJacob Faibussowitsch     PetscCall((sf->ops->GetGraph)(sf, nroots, nleaves, ilocal, iremote));
774b8dee149SJunchao Zhang   } else {
77595fce210SBarry Smith     if (nroots) *nroots = sf->nroots;
77695fce210SBarry Smith     if (nleaves) *nleaves = sf->nleaves;
77795fce210SBarry Smith     if (ilocal) *ilocal = sf->mine;
77895fce210SBarry Smith     if (iremote) *iremote = sf->remote;
779b8dee149SJunchao Zhang   }
78095fce210SBarry Smith   PetscFunctionReturn(0);
78195fce210SBarry Smith }
78295fce210SBarry Smith 
78329046d53SLisandro Dalcin /*@
78495fce210SBarry Smith    PetscSFGetLeafRange - Get the active leaf ranges
78595fce210SBarry Smith 
78695fce210SBarry Smith    Not Collective
78795fce210SBarry Smith 
7884165533cSJose E. Roman    Input Parameter:
78995fce210SBarry Smith .  sf - star forest
79095fce210SBarry Smith 
7914165533cSJose E. Roman    Output Parameters:
792dd5b3ca6SJunchao Zhang +  minleaf - minimum active leaf on this process. Return 0 if there are no leaves.
793dd5b3ca6SJunchao Zhang -  maxleaf - maximum active leaf on this process. Return -1 if there are no leaves.
79495fce210SBarry Smith 
79595fce210SBarry Smith    Level: developer
79695fce210SBarry Smith 
797db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
79895fce210SBarry Smith @*/
799*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetLeafRange(PetscSF sf, PetscInt *minleaf, PetscInt *maxleaf)
800*d71ae5a4SJacob Faibussowitsch {
80195fce210SBarry Smith   PetscFunctionBegin;
80295fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
80329046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf, 1);
80495fce210SBarry Smith   if (minleaf) *minleaf = sf->minleaf;
80595fce210SBarry Smith   if (maxleaf) *maxleaf = sf->maxleaf;
80695fce210SBarry Smith   PetscFunctionReturn(0);
80795fce210SBarry Smith }
80895fce210SBarry Smith 
80995fce210SBarry Smith /*@C
810fe2efc57SMark    PetscSFViewFromOptions - View from Options
811fe2efc57SMark 
812fe2efc57SMark    Collective on PetscSF
813fe2efc57SMark 
814fe2efc57SMark    Input Parameters:
815fe2efc57SMark +  A - the star forest
816736c3998SJose E. Roman .  obj - Optional object
817736c3998SJose E. Roman -  name - command line option
818fe2efc57SMark 
819fe2efc57SMark    Level: intermediate
820db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFView`, `PetscObjectViewFromOptions()`, `PetscSFCreate()`
821fe2efc57SMark @*/
822*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFViewFromOptions(PetscSF A, PetscObject obj, const char name[])
823*d71ae5a4SJacob Faibussowitsch {
824fe2efc57SMark   PetscFunctionBegin;
825fe2efc57SMark   PetscValidHeaderSpecific(A, PETSCSF_CLASSID, 1);
8269566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
827fe2efc57SMark   PetscFunctionReturn(0);
828fe2efc57SMark }
829fe2efc57SMark 
830fe2efc57SMark /*@C
83195fce210SBarry Smith    PetscSFView - view a star forest
83295fce210SBarry Smith 
83395fce210SBarry Smith    Collective
83495fce210SBarry Smith 
8354165533cSJose E. Roman    Input Parameters:
83695fce210SBarry Smith +  sf - star forest
83795fce210SBarry Smith -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD
83895fce210SBarry Smith 
83995fce210SBarry Smith    Level: beginner
84095fce210SBarry Smith 
841db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFSetGraph()`
84295fce210SBarry Smith @*/
843*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFView(PetscSF sf, PetscViewer viewer)
844*d71ae5a4SJacob Faibussowitsch {
84595fce210SBarry Smith   PetscBool         iascii;
84695fce210SBarry Smith   PetscViewerFormat format;
84795fce210SBarry Smith 
84895fce210SBarry Smith   PetscFunctionBegin;
84995fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
8509566063dSJacob Faibussowitsch   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf), &viewer));
85195fce210SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
85295fce210SBarry Smith   PetscCheckSameComm(sf, 1, viewer, 2);
8539566063dSJacob Faibussowitsch   if (sf->graphset) PetscCall(PetscSFSetUp(sf));
8549566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
85553dd6d7dSJunchao Zhang   if (iascii && viewer->format != PETSC_VIEWER_ASCII_MATLAB) {
85695fce210SBarry Smith     PetscMPIInt rank;
85781bfa7aaSJed Brown     PetscInt    ii, i, j;
85895fce210SBarry Smith 
8599566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)sf, viewer));
8609566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
861dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
86280153354SVaclav Hapla       if (!sf->graphset) {
8639566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "PetscSFSetGraph() has not been called yet\n"));
8649566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
86580153354SVaclav Hapla         PetscFunctionReturn(0);
86680153354SVaclav Hapla       }
8679566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
8689566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
8699566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Number of roots=%" PetscInt_FMT ", leaves=%" PetscInt_FMT ", remote ranks=%" PetscInt_FMT "\n", rank, sf->nroots, sf->nleaves, sf->nranks));
87048a46eb9SPierre Jolivet       for (i = 0; i < sf->nleaves; i++) PetscCall(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));
8719566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
8729566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetFormat(viewer, &format));
87395fce210SBarry Smith       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
87481bfa7aaSJed Brown         PetscMPIInt *tmpranks, *perm;
8759566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(sf->nranks, &tmpranks, sf->nranks, &perm));
8769566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(tmpranks, sf->ranks, sf->nranks));
87781bfa7aaSJed Brown         for (i = 0; i < sf->nranks; i++) perm[i] = i;
8789566063dSJacob Faibussowitsch         PetscCall(PetscSortMPIIntWithArray(sf->nranks, tmpranks, perm));
8799566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Roots referenced by my leaves, by rank\n", rank));
88081bfa7aaSJed Brown         for (ii = 0; ii < sf->nranks; ii++) {
88181bfa7aaSJed Brown           i = perm[ii];
8829566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %d: %" PetscInt_FMT " edges\n", rank, sf->ranks[i], sf->roffset[i + 1] - sf->roffset[i]));
88348a46eb9SPierre Jolivet           for (j = sf->roffset[i]; j < sf->roffset[i + 1]; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d]    %" PetscInt_FMT " <- %" PetscInt_FMT "\n", rank, sf->rmine[j], sf->rremote[j]));
88495fce210SBarry Smith         }
8859566063dSJacob Faibussowitsch         PetscCall(PetscFree2(tmpranks, perm));
88695fce210SBarry Smith       }
8879566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
8889566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
889dd5b3ca6SJunchao Zhang     }
8909566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
89195fce210SBarry Smith   }
892dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf, View, viewer);
89395fce210SBarry Smith   PetscFunctionReturn(0);
89495fce210SBarry Smith }
89595fce210SBarry Smith 
89695fce210SBarry Smith /*@C
897dec1416fSJunchao Zhang    PetscSFGetRootRanks - Get root ranks and number of vertices referenced by leaves on this process
89895fce210SBarry Smith 
89995fce210SBarry Smith    Not Collective
90095fce210SBarry Smith 
9014165533cSJose E. Roman    Input Parameter:
90295fce210SBarry Smith .  sf - star forest
90395fce210SBarry Smith 
9044165533cSJose E. Roman    Output Parameters:
90595fce210SBarry Smith +  nranks - number of ranks referenced by local part
90695fce210SBarry Smith .  ranks - array of ranks
90795fce210SBarry Smith .  roffset - offset in rmine/rremote for each rank (length nranks+1)
90895fce210SBarry Smith .  rmine - concatenated array holding local indices referencing each remote rank
90995fce210SBarry Smith -  rremote - concatenated array holding remote indices referenced for each remote rank
91095fce210SBarry Smith 
91195fce210SBarry Smith    Level: developer
91295fce210SBarry Smith 
913db781477SPatrick Sanan .seealso: `PetscSFGetLeafRanks()`
91495fce210SBarry Smith @*/
915*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetRootRanks(PetscSF sf, PetscInt *nranks, const PetscMPIInt **ranks, const PetscInt **roffset, const PetscInt **rmine, const PetscInt **rremote)
916*d71ae5a4SJacob Faibussowitsch {
91795fce210SBarry Smith   PetscFunctionBegin;
91895fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
91928b400f6SJacob Faibussowitsch   PetscCheck(sf->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFSetUp() before obtaining ranks");
920dec1416fSJunchao Zhang   if (sf->ops->GetRootRanks) {
9219566063dSJacob Faibussowitsch     PetscCall((sf->ops->GetRootRanks)(sf, nranks, ranks, roffset, rmine, rremote));
922dec1416fSJunchao Zhang   } else {
923dec1416fSJunchao Zhang     /* The generic implementation */
92495fce210SBarry Smith     if (nranks) *nranks = sf->nranks;
92595fce210SBarry Smith     if (ranks) *ranks = sf->ranks;
92695fce210SBarry Smith     if (roffset) *roffset = sf->roffset;
92795fce210SBarry Smith     if (rmine) *rmine = sf->rmine;
92895fce210SBarry Smith     if (rremote) *rremote = sf->rremote;
929dec1416fSJunchao Zhang   }
93095fce210SBarry Smith   PetscFunctionReturn(0);
93195fce210SBarry Smith }
93295fce210SBarry Smith 
9338750ddebSJunchao Zhang /*@C
9348750ddebSJunchao Zhang    PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process
9358750ddebSJunchao Zhang 
9368750ddebSJunchao Zhang    Not Collective
9378750ddebSJunchao Zhang 
9384165533cSJose E. Roman    Input Parameter:
9398750ddebSJunchao Zhang .  sf - star forest
9408750ddebSJunchao Zhang 
9414165533cSJose E. Roman    Output Parameters:
9428750ddebSJunchao Zhang +  niranks - number of leaf ranks referencing roots on this process
9438750ddebSJunchao Zhang .  iranks - array of ranks
9448750ddebSJunchao Zhang .  ioffset - offset in irootloc for each rank (length niranks+1)
9458750ddebSJunchao Zhang -  irootloc - concatenated array holding local indices of roots referenced by each leaf rank
9468750ddebSJunchao Zhang 
9478750ddebSJunchao Zhang    Level: developer
9488750ddebSJunchao Zhang 
949db781477SPatrick Sanan .seealso: `PetscSFGetRootRanks()`
9508750ddebSJunchao Zhang @*/
951*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetLeafRanks(PetscSF sf, PetscInt *niranks, const PetscMPIInt **iranks, const PetscInt **ioffset, const PetscInt **irootloc)
952*d71ae5a4SJacob Faibussowitsch {
9538750ddebSJunchao Zhang   PetscFunctionBegin;
9548750ddebSJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
95528b400f6SJacob Faibussowitsch   PetscCheck(sf->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFSetUp() before obtaining ranks");
9568750ddebSJunchao Zhang   if (sf->ops->GetLeafRanks) {
9579566063dSJacob Faibussowitsch     PetscCall((sf->ops->GetLeafRanks)(sf, niranks, iranks, ioffset, irootloc));
9588750ddebSJunchao Zhang   } else {
9598750ddebSJunchao Zhang     PetscSFType type;
9609566063dSJacob Faibussowitsch     PetscCall(PetscSFGetType(sf, &type));
96198921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type);
9628750ddebSJunchao Zhang   }
9638750ddebSJunchao Zhang   PetscFunctionReturn(0);
9648750ddebSJunchao Zhang }
9658750ddebSJunchao Zhang 
966*d71ae5a4SJacob Faibussowitsch static PetscBool InList(PetscMPIInt needle, PetscMPIInt n, const PetscMPIInt *list)
967*d71ae5a4SJacob Faibussowitsch {
968b5a8e515SJed Brown   PetscInt i;
969b5a8e515SJed Brown   for (i = 0; i < n; i++) {
970b5a8e515SJed Brown     if (needle == list[i]) return PETSC_TRUE;
971b5a8e515SJed Brown   }
972b5a8e515SJed Brown   return PETSC_FALSE;
973b5a8e515SJed Brown }
974b5a8e515SJed Brown 
97595fce210SBarry Smith /*@C
97621c688dcSJed Brown    PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations.
97721c688dcSJed Brown 
97821c688dcSJed Brown    Collective
97921c688dcSJed Brown 
9804165533cSJose E. Roman    Input Parameters:
981b5a8e515SJed Brown +  sf - PetscSF to set up; PetscSFSetGraph() must have been called
982b5a8e515SJed Brown -  dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange)
98321c688dcSJed Brown 
98421c688dcSJed Brown    Level: developer
98521c688dcSJed Brown 
986db781477SPatrick Sanan .seealso: `PetscSFGetRootRanks()`
98721c688dcSJed Brown @*/
988*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetUpRanks(PetscSF sf, MPI_Group dgroup)
989*d71ae5a4SJacob Faibussowitsch {
99021c688dcSJed Brown   PetscTable         table;
99121c688dcSJed Brown   PetscTablePosition pos;
992b5a8e515SJed Brown   PetscMPIInt        size, groupsize, *groupranks;
993247e8311SStefano Zampini   PetscInt          *rcount, *ranks;
994247e8311SStefano Zampini   PetscInt           i, irank = -1, orank = -1;
99521c688dcSJed Brown 
99621c688dcSJed Brown   PetscFunctionBegin;
99721c688dcSJed Brown   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
99829046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf, 1);
9999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf), &size));
10009566063dSJacob Faibussowitsch   PetscCall(PetscTableCreate(10, size, &table));
100121c688dcSJed Brown   for (i = 0; i < sf->nleaves; i++) {
100221c688dcSJed Brown     /* Log 1-based rank */
10039566063dSJacob Faibussowitsch     PetscCall(PetscTableAdd(table, sf->remote[i].rank + 1, 1, ADD_VALUES));
100421c688dcSJed Brown   }
10059566063dSJacob Faibussowitsch   PetscCall(PetscTableGetCount(table, &sf->nranks));
10069566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(sf->nranks, &sf->ranks, sf->nranks + 1, &sf->roffset, sf->nleaves, &sf->rmine, sf->nleaves, &sf->rremote));
10079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(sf->nranks, &rcount, sf->nranks, &ranks));
10089566063dSJacob Faibussowitsch   PetscCall(PetscTableGetHeadPosition(table, &pos));
100921c688dcSJed Brown   for (i = 0; i < sf->nranks; i++) {
10109566063dSJacob Faibussowitsch     PetscCall(PetscTableGetNext(table, &pos, &ranks[i], &rcount[i]));
101121c688dcSJed Brown     ranks[i]--; /* Convert back to 0-based */
101221c688dcSJed Brown   }
10139566063dSJacob Faibussowitsch   PetscCall(PetscTableDestroy(&table));
1014b5a8e515SJed Brown 
1015b5a8e515SJed Brown   /* We expect that dgroup is reliably "small" while nranks could be large */
1016b5a8e515SJed Brown   {
10177fb8a5e4SKarl Rupp     MPI_Group    group = MPI_GROUP_NULL;
1018b5a8e515SJed Brown     PetscMPIInt *dgroupranks;
10199566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)sf), &group));
10209566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_size(dgroup, &groupsize));
10219566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(groupsize, &dgroupranks));
10229566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(groupsize, &groupranks));
1023b5a8e515SJed Brown     for (i = 0; i < groupsize; i++) dgroupranks[i] = i;
10249566063dSJacob Faibussowitsch     if (groupsize) PetscCallMPI(MPI_Group_translate_ranks(dgroup, groupsize, dgroupranks, group, groupranks));
10259566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_free(&group));
10269566063dSJacob Faibussowitsch     PetscCall(PetscFree(dgroupranks));
1027b5a8e515SJed Brown   }
1028b5a8e515SJed Brown 
1029b5a8e515SJed Brown   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
1030b5a8e515SJed Brown   for (sf->ndranks = 0, i = sf->nranks; sf->ndranks < i;) {
1031b5a8e515SJed Brown     for (i--; sf->ndranks < i; i--) { /* Scan i backward looking for distinguished rank */
1032b5a8e515SJed Brown       if (InList(ranks[i], groupsize, groupranks)) break;
1033b5a8e515SJed Brown     }
1034b5a8e515SJed Brown     for (; sf->ndranks <= i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
1035b5a8e515SJed Brown       if (!InList(ranks[sf->ndranks], groupsize, groupranks)) break;
1036b5a8e515SJed Brown     }
1037b5a8e515SJed Brown     if (sf->ndranks < i) { /* Swap ranks[sf->ndranks] with ranks[i] */
1038b5a8e515SJed Brown       PetscInt tmprank, tmpcount;
1039247e8311SStefano Zampini 
1040b5a8e515SJed Brown       tmprank             = ranks[i];
1041b5a8e515SJed Brown       tmpcount            = rcount[i];
1042b5a8e515SJed Brown       ranks[i]            = ranks[sf->ndranks];
1043b5a8e515SJed Brown       rcount[i]           = rcount[sf->ndranks];
1044b5a8e515SJed Brown       ranks[sf->ndranks]  = tmprank;
1045b5a8e515SJed Brown       rcount[sf->ndranks] = tmpcount;
1046b5a8e515SJed Brown       sf->ndranks++;
1047b5a8e515SJed Brown     }
1048b5a8e515SJed Brown   }
10499566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupranks));
10509566063dSJacob Faibussowitsch   PetscCall(PetscSortIntWithArray(sf->ndranks, ranks, rcount));
10519566063dSJacob Faibussowitsch   PetscCall(PetscSortIntWithArray(sf->nranks - sf->ndranks, ranks + sf->ndranks, rcount + sf->ndranks));
105221c688dcSJed Brown   sf->roffset[0] = 0;
105321c688dcSJed Brown   for (i = 0; i < sf->nranks; i++) {
10549566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(ranks[i], sf->ranks + i));
105521c688dcSJed Brown     sf->roffset[i + 1] = sf->roffset[i] + rcount[i];
105621c688dcSJed Brown     rcount[i]          = 0;
105721c688dcSJed Brown   }
1058247e8311SStefano Zampini   for (i = 0, irank = -1, orank = -1; i < sf->nleaves; i++) {
1059247e8311SStefano Zampini     /* short circuit */
1060247e8311SStefano Zampini     if (orank != sf->remote[i].rank) {
106121c688dcSJed Brown       /* Search for index of iremote[i].rank in sf->ranks */
10629566063dSJacob Faibussowitsch       PetscCall(PetscFindMPIInt(sf->remote[i].rank, sf->ndranks, sf->ranks, &irank));
1063b5a8e515SJed Brown       if (irank < 0) {
10649566063dSJacob Faibussowitsch         PetscCall(PetscFindMPIInt(sf->remote[i].rank, sf->nranks - sf->ndranks, sf->ranks + sf->ndranks, &irank));
1065b5a8e515SJed Brown         if (irank >= 0) irank += sf->ndranks;
106621c688dcSJed Brown       }
1067247e8311SStefano Zampini       orank = sf->remote[i].rank;
1068247e8311SStefano Zampini     }
106908401ef6SPierre Jolivet     PetscCheck(irank >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find rank %" PetscInt_FMT " in array", sf->remote[i].rank);
107021c688dcSJed Brown     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
107121c688dcSJed Brown     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
107221c688dcSJed Brown     rcount[irank]++;
107321c688dcSJed Brown   }
10749566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rcount, ranks));
107521c688dcSJed Brown   PetscFunctionReturn(0);
107621c688dcSJed Brown }
107721c688dcSJed Brown 
107821c688dcSJed Brown /*@C
107995fce210SBarry Smith    PetscSFGetGroups - gets incoming and outgoing process groups
108095fce210SBarry Smith 
108195fce210SBarry Smith    Collective
108295fce210SBarry Smith 
10834165533cSJose E. Roman    Input Parameter:
108495fce210SBarry Smith .  sf - star forest
108595fce210SBarry Smith 
10864165533cSJose E. Roman    Output Parameters:
108795fce210SBarry Smith +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
108895fce210SBarry Smith -  outgoing - group of destination processes for outgoing edges (roots that I reference)
108995fce210SBarry Smith 
109095fce210SBarry Smith    Level: developer
109195fce210SBarry Smith 
1092db781477SPatrick Sanan .seealso: `PetscSFGetWindow()`, `PetscSFRestoreWindow()`
109395fce210SBarry Smith @*/
1094*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGroups(PetscSF sf, MPI_Group *incoming, MPI_Group *outgoing)
1095*d71ae5a4SJacob Faibussowitsch {
10967fb8a5e4SKarl Rupp   MPI_Group group = MPI_GROUP_NULL;
109795fce210SBarry Smith 
109895fce210SBarry Smith   PetscFunctionBegin;
109908401ef6SPierre Jolivet   PetscCheck(sf->nranks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFSetUpRanks() before obtaining groups");
110095fce210SBarry Smith   if (sf->ingroup == MPI_GROUP_NULL) {
110195fce210SBarry Smith     PetscInt        i;
110295fce210SBarry Smith     const PetscInt *indegree;
110395fce210SBarry Smith     PetscMPIInt     rank, *outranks, *inranks;
110495fce210SBarry Smith     PetscSFNode    *remote;
110595fce210SBarry Smith     PetscSF         bgcount;
110695fce210SBarry Smith 
110795fce210SBarry Smith     /* Compute the number of incoming ranks */
11089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(sf->nranks, &remote));
110995fce210SBarry Smith     for (i = 0; i < sf->nranks; i++) {
111095fce210SBarry Smith       remote[i].rank  = sf->ranks[i];
111195fce210SBarry Smith       remote[i].index = 0;
111295fce210SBarry Smith     }
11139566063dSJacob Faibussowitsch     PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_CONFONLY, &bgcount));
11149566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(bgcount, 1, sf->nranks, NULL, PETSC_COPY_VALUES, remote, PETSC_OWN_POINTER));
11159566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(bgcount, &indegree));
11169566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(bgcount, &indegree));
111795fce210SBarry Smith     /* Enumerate the incoming ranks */
11189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(indegree[0], &inranks, sf->nranks, &outranks));
11199566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
112095fce210SBarry Smith     for (i = 0; i < sf->nranks; i++) outranks[i] = rank;
11219566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherBegin(bgcount, MPI_INT, outranks, inranks));
11229566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherEnd(bgcount, MPI_INT, outranks, inranks));
11239566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)sf), &group));
11249566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_incl(group, indegree[0], inranks, &sf->ingroup));
11259566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_free(&group));
11269566063dSJacob Faibussowitsch     PetscCall(PetscFree2(inranks, outranks));
11279566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&bgcount));
112895fce210SBarry Smith   }
112995fce210SBarry Smith   *incoming = sf->ingroup;
113095fce210SBarry Smith 
113195fce210SBarry Smith   if (sf->outgroup == MPI_GROUP_NULL) {
11329566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)sf), &group));
11339566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_incl(group, sf->nranks, sf->ranks, &sf->outgroup));
11349566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_free(&group));
113595fce210SBarry Smith   }
113695fce210SBarry Smith   *outgoing = sf->outgroup;
113795fce210SBarry Smith   PetscFunctionReturn(0);
113895fce210SBarry Smith }
113995fce210SBarry Smith 
114029046d53SLisandro Dalcin /*@
11413b8d980fSPierre Jolivet    PetscSFGetMultiSF - gets the inner SF implementing gathers and scatters
114295fce210SBarry Smith 
114395fce210SBarry Smith    Collective
114495fce210SBarry Smith 
11454165533cSJose E. Roman    Input Parameter:
114695fce210SBarry Smith .  sf - star forest that may contain roots with 0 or with more than 1 vertex
114795fce210SBarry Smith 
11484165533cSJose E. Roman    Output Parameter:
114995fce210SBarry Smith .  multi - star forest with split roots, such that each root has degree exactly 1
115095fce210SBarry Smith 
115195fce210SBarry Smith    Level: developer
115295fce210SBarry Smith 
115395fce210SBarry Smith    Notes:
115495fce210SBarry Smith 
115595fce210SBarry Smith    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
115695fce210SBarry Smith    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
115795fce210SBarry Smith    edge, it is a candidate for future optimization that might involve its removal.
115895fce210SBarry Smith 
1159db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`, `PetscSFGatherBegin()`, `PetscSFScatterBegin()`, `PetscSFComputeMultiRootOriginalNumbering()`
116095fce210SBarry Smith @*/
1161*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetMultiSF(PetscSF sf, PetscSF *multi)
1162*d71ae5a4SJacob Faibussowitsch {
116395fce210SBarry Smith   PetscFunctionBegin;
116495fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
116595fce210SBarry Smith   PetscValidPointer(multi, 2);
116695fce210SBarry Smith   if (sf->nroots < 0) { /* Graph has not been set yet; why do we need this? */
11679566063dSJacob Faibussowitsch     PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, &sf->multi));
116895fce210SBarry Smith     *multi           = sf->multi;
1169013b3241SStefano Zampini     sf->multi->multi = sf->multi;
117095fce210SBarry Smith     PetscFunctionReturn(0);
117195fce210SBarry Smith   }
117295fce210SBarry Smith   if (!sf->multi) {
117395fce210SBarry Smith     const PetscInt *indegree;
11749837ea96SMatthew G. Knepley     PetscInt        i, *inoffset, *outones, *outoffset, maxlocal;
117595fce210SBarry Smith     PetscSFNode    *remote;
117629046d53SLisandro Dalcin     maxlocal = sf->maxleaf + 1; /* TODO: We should use PetscSFGetLeafRange() */
11779566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf, &indegree));
11789566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf, &indegree));
11799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(sf->nroots + 1, &inoffset, maxlocal, &outones, maxlocal, &outoffset));
118095fce210SBarry Smith     inoffset[0] = 0;
118195fce210SBarry Smith     for (i = 0; i < sf->nroots; i++) inoffset[i + 1] = inoffset[i] + indegree[i];
11829837ea96SMatthew G. Knepley     for (i = 0; i < maxlocal; i++) outones[i] = 1;
11839566063dSJacob Faibussowitsch     PetscCall(PetscSFFetchAndOpBegin(sf, MPIU_INT, inoffset, outones, outoffset, MPI_SUM));
11849566063dSJacob Faibussowitsch     PetscCall(PetscSFFetchAndOpEnd(sf, MPIU_INT, inoffset, outones, outoffset, MPI_SUM));
118595fce210SBarry Smith     for (i = 0; i < sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
118676bd3646SJed Brown     if (PetscDefined(USE_DEBUG)) {                               /* Check that the expected number of increments occurred */
1187ad540459SPierre Jolivet       for (i = 0; i < sf->nroots; i++) PetscCheck(inoffset[i] + indegree[i] == inoffset[i + 1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrect result after PetscSFFetchAndOp");
118876bd3646SJed Brown     }
11899566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(sf->nleaves, &remote));
119095fce210SBarry Smith     for (i = 0; i < sf->nleaves; i++) {
119195fce210SBarry Smith       remote[i].rank  = sf->remote[i].rank;
119238e7336fSToby Isaac       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
119395fce210SBarry Smith     }
11949566063dSJacob Faibussowitsch     PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, &sf->multi));
1195013b3241SStefano Zampini     sf->multi->multi = sf->multi;
11969566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(sf->multi, inoffset[sf->nroots], sf->nleaves, sf->mine, PETSC_COPY_VALUES, remote, PETSC_OWN_POINTER));
119795fce210SBarry Smith     if (sf->rankorder) { /* Sort the ranks */
119895fce210SBarry Smith       PetscMPIInt  rank;
119995fce210SBarry Smith       PetscInt    *inranks, *newoffset, *outranks, *newoutoffset, *tmpoffset, maxdegree;
120095fce210SBarry Smith       PetscSFNode *newremote;
12019566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
120295fce210SBarry Smith       for (i = 0, maxdegree = 0; i < sf->nroots; i++) maxdegree = PetscMax(maxdegree, indegree[i]);
12039566063dSJacob Faibussowitsch       PetscCall(PetscMalloc5(sf->multi->nroots, &inranks, sf->multi->nroots, &newoffset, maxlocal, &outranks, maxlocal, &newoutoffset, maxdegree, &tmpoffset));
12049837ea96SMatthew G. Knepley       for (i = 0; i < maxlocal; i++) outranks[i] = rank;
12059566063dSJacob Faibussowitsch       PetscCall(PetscSFReduceBegin(sf->multi, MPIU_INT, outranks, inranks, MPI_REPLACE));
12069566063dSJacob Faibussowitsch       PetscCall(PetscSFReduceEnd(sf->multi, MPIU_INT, outranks, inranks, MPI_REPLACE));
120795fce210SBarry Smith       /* Sort the incoming ranks at each vertex, build the inverse map */
120895fce210SBarry Smith       for (i = 0; i < sf->nroots; i++) {
120995fce210SBarry Smith         PetscInt j;
121095fce210SBarry Smith         for (j = 0; j < indegree[i]; j++) tmpoffset[j] = j;
12119566063dSJacob Faibussowitsch         PetscCall(PetscSortIntWithArray(indegree[i], inranks + inoffset[i], tmpoffset));
121295fce210SBarry Smith         for (j = 0; j < indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
121395fce210SBarry Smith       }
12149566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(sf->multi, MPIU_INT, newoffset, newoutoffset, MPI_REPLACE));
12159566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(sf->multi, MPIU_INT, newoffset, newoutoffset, MPI_REPLACE));
12169566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(sf->nleaves, &newremote));
121795fce210SBarry Smith       for (i = 0; i < sf->nleaves; i++) {
121895fce210SBarry Smith         newremote[i].rank  = sf->remote[i].rank;
121901365b40SToby Isaac         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
122095fce210SBarry Smith       }
12219566063dSJacob Faibussowitsch       PetscCall(PetscSFSetGraph(sf->multi, inoffset[sf->nroots], sf->nleaves, sf->mine, PETSC_COPY_VALUES, newremote, PETSC_OWN_POINTER));
12229566063dSJacob Faibussowitsch       PetscCall(PetscFree5(inranks, newoffset, outranks, newoutoffset, tmpoffset));
122395fce210SBarry Smith     }
12249566063dSJacob Faibussowitsch     PetscCall(PetscFree3(inoffset, outones, outoffset));
122595fce210SBarry Smith   }
122695fce210SBarry Smith   *multi = sf->multi;
122795fce210SBarry Smith   PetscFunctionReturn(0);
122895fce210SBarry Smith }
122995fce210SBarry Smith 
123095fce210SBarry Smith /*@C
123172502a1fSJunchao Zhang    PetscSFCreateEmbeddedRootSF - removes edges from all but the selected roots, does not remap indices
123295fce210SBarry Smith 
123395fce210SBarry Smith    Collective
123495fce210SBarry Smith 
12354165533cSJose E. Roman    Input Parameters:
123695fce210SBarry Smith +  sf - original star forest
1237ba2a7774SJunchao Zhang .  nselected  - number of selected roots on this process
1238ba2a7774SJunchao Zhang -  selected   - indices of the selected roots on this process
123995fce210SBarry Smith 
12404165533cSJose E. Roman    Output Parameter:
1241cd620004SJunchao Zhang .  esf - new star forest
124295fce210SBarry Smith 
124395fce210SBarry Smith    Level: advanced
124495fce210SBarry Smith 
124595fce210SBarry Smith    Note:
124695fce210SBarry Smith    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
124795fce210SBarry Smith    be done by calling PetscSFGetGraph().
124895fce210SBarry Smith 
1249db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`, `PetscSFGetGraph()`
125095fce210SBarry Smith @*/
1251*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateEmbeddedRootSF(PetscSF sf, PetscInt nselected, const PetscInt *selected, PetscSF *esf)
1252*d71ae5a4SJacob Faibussowitsch {
1253cd620004SJunchao Zhang   PetscInt           i, j, n, nroots, nleaves, esf_nleaves, *new_ilocal, minleaf, maxleaf, maxlocal;
1254cd620004SJunchao Zhang   const PetscInt    *ilocal;
1255cd620004SJunchao Zhang   signed char       *rootdata, *leafdata, *leafmem;
1256ba2a7774SJunchao Zhang   const PetscSFNode *iremote;
1257f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1258f659e5c7SJunchao Zhang   MPI_Comm           comm;
125995fce210SBarry Smith 
126095fce210SBarry Smith   PetscFunctionBegin;
126195fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
126229046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf, 1);
1263dadcf809SJacob Faibussowitsch   if (nselected) PetscValidIntPointer(selected, 3);
1264cd620004SJunchao Zhang   PetscValidPointer(esf, 4);
12650511a646SMatthew G. Knepley 
12669566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
12679566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_EmbedSF, sf, 0, 0, 0));
12689566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
12699566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, &iremote));
1270cd620004SJunchao Zhang 
127176bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) { /* Error out if selected[] has dups or  out of range indices */
1272cd620004SJunchao Zhang     PetscBool dups;
12739566063dSJacob Faibussowitsch     PetscCall(PetscCheckDupsInt(nselected, selected, &dups));
127428b400f6SJacob Faibussowitsch     PetscCheck(!dups, comm, PETSC_ERR_ARG_WRONG, "selected[] has dups");
12759371c9d4SSatish Balay     for (i = 0; i < nselected; i++) PetscCheck(selected[i] >= 0 && selected[i] < nroots, comm, PETSC_ERR_ARG_OUTOFRANGE, "selected root indice %" PetscInt_FMT " is out of [0,%" PetscInt_FMT ")", selected[i], nroots);
1276cd620004SJunchao Zhang   }
1277f659e5c7SJunchao Zhang 
1278dbbe0bcdSBarry Smith   if (sf->ops->CreateEmbeddedRootSF) PetscUseTypeMethod(sf, CreateEmbeddedRootSF, nselected, selected, esf);
1279dbbe0bcdSBarry Smith   else {
1280cd620004SJunchao Zhang     /* A generic version of creating embedded sf */
12819566063dSJacob Faibussowitsch     PetscCall(PetscSFGetLeafRange(sf, &minleaf, &maxleaf));
1282cd620004SJunchao Zhang     maxlocal = maxleaf - minleaf + 1;
12839566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(nroots, &rootdata, maxlocal, &leafmem));
1284cd620004SJunchao Zhang     leafdata = leafmem - minleaf;
1285cd620004SJunchao Zhang     /* Tag selected roots and bcast to leaves */
1286cd620004SJunchao Zhang     for (i = 0; i < nselected; i++) rootdata[selected[i]] = 1;
12879566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPI_SIGNED_CHAR, rootdata, leafdata, MPI_REPLACE));
12889566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPI_SIGNED_CHAR, rootdata, leafdata, MPI_REPLACE));
1289ba2a7774SJunchao Zhang 
1290cd620004SJunchao Zhang     /* Build esf with leaves that are still connected */
1291cd620004SJunchao Zhang     esf_nleaves = 0;
1292cd620004SJunchao Zhang     for (i = 0; i < nleaves; i++) {
1293cd620004SJunchao Zhang       j = ilocal ? ilocal[i] : i;
1294cd620004SJunchao Zhang       /* esf_nleaves += leafdata[j] should work in theory, but failed with SFWindow bugs
1295cd620004SJunchao Zhang          with PetscSFBcast. See https://gitlab.com/petsc/petsc/issues/555
1296cd620004SJunchao Zhang       */
1297cd620004SJunchao Zhang       esf_nleaves += (leafdata[j] ? 1 : 0);
1298cd620004SJunchao Zhang     }
12999566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(esf_nleaves, &new_ilocal));
13009566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(esf_nleaves, &new_iremote));
1301cd620004SJunchao Zhang     for (i = n = 0; i < nleaves; i++) {
1302cd620004SJunchao Zhang       j = ilocal ? ilocal[i] : i;
1303cd620004SJunchao Zhang       if (leafdata[j]) {
1304cd620004SJunchao Zhang         new_ilocal[n]        = j;
1305cd620004SJunchao Zhang         new_iremote[n].rank  = iremote[i].rank;
1306cd620004SJunchao Zhang         new_iremote[n].index = iremote[i].index;
1307fc1ede2bSMatthew G. Knepley         ++n;
130895fce210SBarry Smith       }
130995fce210SBarry Smith     }
13109566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(comm, esf));
13119566063dSJacob Faibussowitsch     PetscCall(PetscSFSetFromOptions(*esf));
13129566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*esf, nroots, esf_nleaves, new_ilocal, PETSC_OWN_POINTER, new_iremote, PETSC_OWN_POINTER));
13139566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafmem));
1314f659e5c7SJunchao Zhang   }
13159566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_EmbedSF, sf, 0, 0, 0));
131695fce210SBarry Smith   PetscFunctionReturn(0);
131795fce210SBarry Smith }
131895fce210SBarry Smith 
13192f5fb4c2SMatthew G. Knepley /*@C
13202f5fb4c2SMatthew G. Knepley   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
13212f5fb4c2SMatthew G. Knepley 
13222f5fb4c2SMatthew G. Knepley   Collective
13232f5fb4c2SMatthew G. Knepley 
13244165533cSJose E. Roman   Input Parameters:
13252f5fb4c2SMatthew G. Knepley + sf - original star forest
1326f659e5c7SJunchao Zhang . nselected  - number of selected leaves on this process
1327f659e5c7SJunchao Zhang - selected   - indices of the selected leaves on this process
13282f5fb4c2SMatthew G. Knepley 
13294165533cSJose E. Roman   Output Parameter:
13302f5fb4c2SMatthew G. Knepley .  newsf - new star forest
13312f5fb4c2SMatthew G. Knepley 
13322f5fb4c2SMatthew G. Knepley   Level: advanced
13332f5fb4c2SMatthew G. Knepley 
1334db781477SPatrick Sanan .seealso: `PetscSFCreateEmbeddedRootSF()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
13352f5fb4c2SMatthew G. Knepley @*/
1336*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nselected, const PetscInt *selected, PetscSF *newsf)
1337*d71ae5a4SJacob Faibussowitsch {
1338f659e5c7SJunchao Zhang   const PetscSFNode *iremote;
1339f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1340f659e5c7SJunchao Zhang   const PetscInt    *ilocal;
1341f659e5c7SJunchao Zhang   PetscInt           i, nroots, *leaves, *new_ilocal;
1342f659e5c7SJunchao Zhang   MPI_Comm           comm;
13432f5fb4c2SMatthew G. Knepley 
13442f5fb4c2SMatthew G. Knepley   PetscFunctionBegin;
13452f5fb4c2SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
134629046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf, 1);
1347dadcf809SJacob Faibussowitsch   if (nselected) PetscValidIntPointer(selected, 3);
13482f5fb4c2SMatthew G. Knepley   PetscValidPointer(newsf, 4);
13492f5fb4c2SMatthew G. Knepley 
1350f659e5c7SJunchao Zhang   /* Uniq selected[] and put results in leaves[] */
13519566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
13529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nselected, &leaves));
13539566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(leaves, selected, nselected));
13549566063dSJacob Faibussowitsch   PetscCall(PetscSortedRemoveDupsInt(&nselected, leaves));
135508401ef6SPierre Jolivet   PetscCheck(!nselected || !(leaves[0] < 0 || leaves[nselected - 1] >= sf->nleaves), 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);
1356f659e5c7SJunchao Zhang 
1357f659e5c7SJunchao Zhang   /* Optimize the routine only when sf is setup and hence we can reuse sf's communication pattern */
1358dbbe0bcdSBarry Smith   if (sf->setupcalled && sf->ops->CreateEmbeddedLeafSF) PetscUseTypeMethod(sf, CreateEmbeddedLeafSF, nselected, leaves, newsf);
1359dbbe0bcdSBarry Smith   else {
13609566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sf, &nroots, NULL, &ilocal, &iremote));
13619566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nselected, &new_ilocal));
13629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nselected, &new_iremote));
1363f659e5c7SJunchao Zhang     for (i = 0; i < nselected; ++i) {
1364f659e5c7SJunchao Zhang       const PetscInt l     = leaves[i];
1365f659e5c7SJunchao Zhang       new_ilocal[i]        = ilocal ? ilocal[l] : l;
1366f659e5c7SJunchao Zhang       new_iremote[i].rank  = iremote[l].rank;
1367f659e5c7SJunchao Zhang       new_iremote[i].index = iremote[l].index;
13682f5fb4c2SMatthew G. Knepley     }
13699566063dSJacob Faibussowitsch     PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_CONFONLY, newsf));
13709566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*newsf, nroots, nselected, new_ilocal, PETSC_OWN_POINTER, new_iremote, PETSC_OWN_POINTER));
1371f659e5c7SJunchao Zhang   }
13729566063dSJacob Faibussowitsch   PetscCall(PetscFree(leaves));
13732f5fb4c2SMatthew G. Knepley   PetscFunctionReturn(0);
13742f5fb4c2SMatthew G. Knepley }
13752f5fb4c2SMatthew G. Knepley 
137695fce210SBarry Smith /*@C
1377ad227feaSJunchao Zhang    PetscSFBcastBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastEnd()
13783482bfa8SJunchao Zhang 
13793482bfa8SJunchao Zhang    Collective on PetscSF
13803482bfa8SJunchao Zhang 
13814165533cSJose E. Roman    Input Parameters:
13823482bfa8SJunchao Zhang +  sf - star forest on which to communicate
13833482bfa8SJunchao Zhang .  unit - data type associated with each node
13843482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
13853482bfa8SJunchao Zhang -  op - operation to use for reduction
13863482bfa8SJunchao Zhang 
13874165533cSJose E. Roman    Output Parameter:
13883482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
13893482bfa8SJunchao Zhang 
13903482bfa8SJunchao Zhang    Level: intermediate
13913482bfa8SJunchao Zhang 
1392d0295fc0SJunchao Zhang    Notes:
1393d0295fc0SJunchao Zhang     When petsc is configured with device support, it will use its own mechanism to figure out whether the given data pointers
1394d0295fc0SJunchao Zhang     are host pointers or device pointers, which may incur a noticable cost. If you already knew the info, you should
1395ad227feaSJunchao Zhang     use PetscSFBcastWithMemTypeBegin() instead.
1396db781477SPatrick Sanan .seealso: `PetscSFBcastEnd()`, `PetscSFBcastWithMemTypeBegin()`
13973482bfa8SJunchao Zhang @*/
1398*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFBcastBegin(PetscSF sf, MPI_Datatype unit, const void *rootdata, void *leafdata, MPI_Op op)
1399*d71ae5a4SJacob Faibussowitsch {
1400eb02082bSJunchao Zhang   PetscMemType rootmtype, leafmtype;
14013482bfa8SJunchao Zhang 
14023482bfa8SJunchao Zhang   PetscFunctionBegin;
14033482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
14049566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
14059566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_BcastBegin, sf, 0, 0, 0));
14069566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(rootdata, &rootmtype));
14079566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafdata, &leafmtype));
1408dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, BcastBegin, unit, rootmtype, rootdata, leafmtype, leafdata, op);
14099566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_BcastBegin, sf, 0, 0, 0));
14103482bfa8SJunchao Zhang   PetscFunctionReturn(0);
14113482bfa8SJunchao Zhang }
14123482bfa8SJunchao Zhang 
14133482bfa8SJunchao Zhang /*@C
1414ad227feaSJunchao Zhang    PetscSFBcastWithMemTypeBegin - begin pointwise broadcast with root value being reduced to leaf value with explicit memory types, to be concluded with call to PetscSFBcastEnd()
1415d0295fc0SJunchao Zhang 
1416d0295fc0SJunchao Zhang    Collective on PetscSF
1417d0295fc0SJunchao Zhang 
14184165533cSJose E. Roman    Input Parameters:
1419d0295fc0SJunchao Zhang +  sf - star forest on which to communicate
1420d0295fc0SJunchao Zhang .  unit - data type associated with each node
1421d0295fc0SJunchao Zhang .  rootmtype - memory type of rootdata
1422d0295fc0SJunchao Zhang .  rootdata - buffer to broadcast
1423d0295fc0SJunchao Zhang .  leafmtype - memory type of leafdata
1424d0295fc0SJunchao Zhang -  op - operation to use for reduction
1425d0295fc0SJunchao Zhang 
14264165533cSJose E. Roman    Output Parameter:
1427d0295fc0SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
1428d0295fc0SJunchao Zhang 
1429d0295fc0SJunchao Zhang    Level: intermediate
1430d0295fc0SJunchao Zhang 
1431db781477SPatrick Sanan .seealso: `PetscSFBcastEnd()`, `PetscSFBcastBegin()`
1432d0295fc0SJunchao Zhang @*/
1433*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFBcastWithMemTypeBegin(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, const void *rootdata, PetscMemType leafmtype, void *leafdata, MPI_Op op)
1434*d71ae5a4SJacob Faibussowitsch {
1435d0295fc0SJunchao Zhang   PetscFunctionBegin;
1436d0295fc0SJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
14379566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
14389566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_BcastBegin, sf, 0, 0, 0));
1439dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, BcastBegin, unit, rootmtype, rootdata, leafmtype, leafdata, op);
14409566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_BcastBegin, sf, 0, 0, 0));
1441d0295fc0SJunchao Zhang   PetscFunctionReturn(0);
1442d0295fc0SJunchao Zhang }
1443d0295fc0SJunchao Zhang 
1444d0295fc0SJunchao Zhang /*@C
1445ad227feaSJunchao Zhang    PetscSFBcastEnd - end a broadcast & reduce operation started with PetscSFBcastBegin()
14463482bfa8SJunchao Zhang 
14473482bfa8SJunchao Zhang    Collective
14483482bfa8SJunchao Zhang 
14494165533cSJose E. Roman    Input Parameters:
14503482bfa8SJunchao Zhang +  sf - star forest
14513482bfa8SJunchao Zhang .  unit - data type
14523482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
14533482bfa8SJunchao Zhang -  op - operation to use for reduction
14543482bfa8SJunchao Zhang 
14554165533cSJose E. Roman    Output Parameter:
14563482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
14573482bfa8SJunchao Zhang 
14583482bfa8SJunchao Zhang    Level: intermediate
14593482bfa8SJunchao Zhang 
1460db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`, `PetscSFReduceEnd()`
14613482bfa8SJunchao Zhang @*/
1462*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFBcastEnd(PetscSF sf, MPI_Datatype unit, const void *rootdata, void *leafdata, MPI_Op op)
1463*d71ae5a4SJacob Faibussowitsch {
14643482bfa8SJunchao Zhang   PetscFunctionBegin;
14653482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
14669566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_BcastEnd, sf, 0, 0, 0));
1467dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, BcastEnd, unit, rootdata, leafdata, op);
14689566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_BcastEnd, sf, 0, 0, 0));
14693482bfa8SJunchao Zhang   PetscFunctionReturn(0);
14703482bfa8SJunchao Zhang }
14713482bfa8SJunchao Zhang 
14723482bfa8SJunchao Zhang /*@C
147395fce210SBarry Smith    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
147495fce210SBarry Smith 
147595fce210SBarry Smith    Collective
147695fce210SBarry Smith 
14774165533cSJose E. Roman    Input Parameters:
147895fce210SBarry Smith +  sf - star forest
147995fce210SBarry Smith .  unit - data type
148095fce210SBarry Smith .  leafdata - values to reduce
148195fce210SBarry Smith -  op - reduction operation
148295fce210SBarry Smith 
14834165533cSJose E. Roman    Output Parameter:
148495fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
148595fce210SBarry Smith 
148695fce210SBarry Smith    Level: intermediate
148795fce210SBarry Smith 
1488d0295fc0SJunchao Zhang    Notes:
1489d0295fc0SJunchao Zhang     When petsc is configured with device support, it will use its own mechanism to figure out whether the given data pointers
1490d0295fc0SJunchao Zhang     are host pointers or device pointers, which may incur a noticable cost. If you already knew the info, you should
1491d0295fc0SJunchao Zhang     use PetscSFReduceWithMemTypeBegin() instead.
1492d0295fc0SJunchao Zhang 
1493db781477SPatrick Sanan .seealso: `PetscSFBcastBegin()`, `PetscSFReduceWithMemTypeBegin()`
149495fce210SBarry Smith @*/
1495*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFReduceBegin(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *rootdata, MPI_Op op)
1496*d71ae5a4SJacob Faibussowitsch {
1497eb02082bSJunchao Zhang   PetscMemType rootmtype, leafmtype;
149895fce210SBarry Smith 
149995fce210SBarry Smith   PetscFunctionBegin;
150095fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
15019566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
15029566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_ReduceBegin, sf, 0, 0, 0));
15039566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(rootdata, &rootmtype));
15049566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafdata, &leafmtype));
15059566063dSJacob Faibussowitsch   PetscCall((sf->ops->ReduceBegin)(sf, unit, leafmtype, leafdata, rootmtype, rootdata, op));
15069566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_ReduceBegin, sf, 0, 0, 0));
150795fce210SBarry Smith   PetscFunctionReturn(0);
150895fce210SBarry Smith }
150995fce210SBarry Smith 
151095fce210SBarry Smith /*@C
1511d0295fc0SJunchao Zhang    PetscSFReduceWithMemTypeBegin - begin reduction of leafdata into rootdata with explicit memory types, to be completed with call to PetscSFReduceEnd()
1512d0295fc0SJunchao Zhang 
1513d0295fc0SJunchao Zhang    Collective
1514d0295fc0SJunchao Zhang 
15154165533cSJose E. Roman    Input Parameters:
1516d0295fc0SJunchao Zhang +  sf - star forest
1517d0295fc0SJunchao Zhang .  unit - data type
1518d0295fc0SJunchao Zhang .  leafmtype - memory type of leafdata
1519d0295fc0SJunchao Zhang .  leafdata - values to reduce
1520d0295fc0SJunchao Zhang .  rootmtype - memory type of rootdata
1521d0295fc0SJunchao Zhang -  op - reduction operation
1522d0295fc0SJunchao Zhang 
15234165533cSJose E. Roman    Output Parameter:
1524d0295fc0SJunchao Zhang .  rootdata - result of reduction of values from all leaves of each root
1525d0295fc0SJunchao Zhang 
1526d0295fc0SJunchao Zhang    Level: intermediate
1527d0295fc0SJunchao Zhang 
1528db781477SPatrick Sanan .seealso: `PetscSFBcastBegin()`, `PetscSFReduceBegin()`
1529d0295fc0SJunchao Zhang @*/
1530*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFReduceWithMemTypeBegin(PetscSF sf, MPI_Datatype unit, PetscMemType leafmtype, const void *leafdata, PetscMemType rootmtype, void *rootdata, MPI_Op op)
1531*d71ae5a4SJacob Faibussowitsch {
1532d0295fc0SJunchao Zhang   PetscFunctionBegin;
1533d0295fc0SJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
15349566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
15359566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_ReduceBegin, sf, 0, 0, 0));
15369566063dSJacob Faibussowitsch   PetscCall((sf->ops->ReduceBegin)(sf, unit, leafmtype, leafdata, rootmtype, rootdata, op));
15379566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_ReduceBegin, sf, 0, 0, 0));
1538d0295fc0SJunchao Zhang   PetscFunctionReturn(0);
1539d0295fc0SJunchao Zhang }
1540d0295fc0SJunchao Zhang 
1541d0295fc0SJunchao Zhang /*@C
154295fce210SBarry Smith    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
154395fce210SBarry Smith 
154495fce210SBarry Smith    Collective
154595fce210SBarry Smith 
15464165533cSJose E. Roman    Input Parameters:
154795fce210SBarry Smith +  sf - star forest
154895fce210SBarry Smith .  unit - data type
154995fce210SBarry Smith .  leafdata - values to reduce
155095fce210SBarry Smith -  op - reduction operation
155195fce210SBarry Smith 
15524165533cSJose E. Roman    Output Parameter:
155395fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
155495fce210SBarry Smith 
155595fce210SBarry Smith    Level: intermediate
155695fce210SBarry Smith 
1557db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`, `PetscSFBcastEnd()`
155895fce210SBarry Smith @*/
1559*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFReduceEnd(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *rootdata, MPI_Op op)
1560*d71ae5a4SJacob Faibussowitsch {
156195fce210SBarry Smith   PetscFunctionBegin;
156295fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
15639566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_ReduceEnd, sf, 0, 0, 0));
1564dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, ReduceEnd, unit, leafdata, rootdata, op);
15659566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_ReduceEnd, sf, 0, 0, 0));
156695fce210SBarry Smith   PetscFunctionReturn(0);
156795fce210SBarry Smith }
156895fce210SBarry Smith 
156995fce210SBarry Smith /*@C
1570a1729e3fSJunchao Zhang    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
1571a1729e3fSJunchao Zhang 
1572a1729e3fSJunchao Zhang    Collective
1573a1729e3fSJunchao Zhang 
15744165533cSJose E. Roman    Input Parameters:
1575a1729e3fSJunchao Zhang +  sf - star forest
1576a1729e3fSJunchao Zhang .  unit - data type
1577a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1578a1729e3fSJunchao Zhang -  op - operation to use for reduction
1579a1729e3fSJunchao Zhang 
15804165533cSJose E. Roman    Output Parameters:
1581a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1582a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1583a1729e3fSJunchao Zhang 
1584a1729e3fSJunchao Zhang    Level: advanced
1585a1729e3fSJunchao Zhang 
1586a1729e3fSJunchao Zhang    Note:
1587a1729e3fSJunchao Zhang    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
1588a1729e3fSJunchao Zhang    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
1589a1729e3fSJunchao Zhang    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
1590a1729e3fSJunchao Zhang    integers.
1591a1729e3fSJunchao Zhang 
1592db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeBegin()`, `PetscSFReduceBegin()`, `PetscSFSetGraph()`
1593a1729e3fSJunchao Zhang @*/
1594*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf, MPI_Datatype unit, void *rootdata, const void *leafdata, void *leafupdate, MPI_Op op)
1595*d71ae5a4SJacob Faibussowitsch {
1596eb02082bSJunchao Zhang   PetscMemType rootmtype, leafmtype, leafupdatemtype;
1597a1729e3fSJunchao Zhang 
1598a1729e3fSJunchao Zhang   PetscFunctionBegin;
1599a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
16009566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
16019566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_FetchAndOpBegin, sf, 0, 0, 0));
16029566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(rootdata, &rootmtype));
16039566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafdata, &leafmtype));
16049566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafupdate, &leafupdatemtype));
160508401ef6SPierre Jolivet   PetscCheck(leafmtype == leafupdatemtype, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for leafdata and leafupdate in different memory types");
1606dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, FetchAndOpBegin, unit, rootmtype, rootdata, leafmtype, leafdata, leafupdate, op);
16079566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_FetchAndOpBegin, sf, 0, 0, 0));
1608a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1609a1729e3fSJunchao Zhang }
1610a1729e3fSJunchao Zhang 
1611a1729e3fSJunchao Zhang /*@C
1612d3b3e55cSJunchao Zhang    PetscSFFetchAndOpWithMemTypeBegin - begin operation with explicit memory types that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
1613d3b3e55cSJunchao Zhang 
1614d3b3e55cSJunchao Zhang    Collective
1615d3b3e55cSJunchao Zhang 
1616d3b3e55cSJunchao Zhang    Input Parameters:
1617d3b3e55cSJunchao Zhang +  sf - star forest
1618d3b3e55cSJunchao Zhang .  unit - data type
1619d3b3e55cSJunchao Zhang .  rootmtype - memory type of rootdata
1620d3b3e55cSJunchao Zhang .  leafmtype - memory type of leafdata
1621d3b3e55cSJunchao Zhang .  leafdata - leaf values to use in reduction
1622d3b3e55cSJunchao Zhang .  leafupdatemtype - memory type of leafupdate
1623d3b3e55cSJunchao Zhang -  op - operation to use for reduction
1624d3b3e55cSJunchao Zhang 
1625d3b3e55cSJunchao Zhang    Output Parameters:
1626d3b3e55cSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1627d3b3e55cSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1628d3b3e55cSJunchao Zhang 
1629d3b3e55cSJunchao Zhang    Level: advanced
1630d3b3e55cSJunchao Zhang 
1631d3b3e55cSJunchao Zhang    Note: See PetscSFFetchAndOpBegin() for more details.
1632d3b3e55cSJunchao Zhang 
1633c2e3fba1SPatrick Sanan .seealso: `PetscSFFetchAndOpBegin()`, `PetscSFComputeDegreeBegin()`, `PetscSFReduceBegin()`, `PetscSFSetGraph()`
1634d3b3e55cSJunchao Zhang @*/
1635*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFFetchAndOpWithMemTypeBegin(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, void *rootdata, PetscMemType leafmtype, const void *leafdata, PetscMemType leafupdatemtype, void *leafupdate, MPI_Op op)
1636*d71ae5a4SJacob Faibussowitsch {
1637d3b3e55cSJunchao Zhang   PetscFunctionBegin;
1638d3b3e55cSJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
16399566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
16409566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_FetchAndOpBegin, sf, 0, 0, 0));
164108401ef6SPierre Jolivet   PetscCheck(leafmtype == leafupdatemtype, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for leafdata and leafupdate in different memory types");
1642dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, FetchAndOpBegin, unit, rootmtype, rootdata, leafmtype, leafdata, leafupdate, op);
16439566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_FetchAndOpBegin, sf, 0, 0, 0));
1644d3b3e55cSJunchao Zhang   PetscFunctionReturn(0);
1645d3b3e55cSJunchao Zhang }
1646d3b3e55cSJunchao Zhang 
1647d3b3e55cSJunchao Zhang /*@C
1648a1729e3fSJunchao 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
1649a1729e3fSJunchao Zhang 
1650a1729e3fSJunchao Zhang    Collective
1651a1729e3fSJunchao Zhang 
16524165533cSJose E. Roman    Input Parameters:
1653a1729e3fSJunchao Zhang +  sf - star forest
1654a1729e3fSJunchao Zhang .  unit - data type
1655a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1656a1729e3fSJunchao Zhang -  op - operation to use for reduction
1657a1729e3fSJunchao Zhang 
16584165533cSJose E. Roman    Output Parameters:
1659a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1660a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1661a1729e3fSJunchao Zhang 
1662a1729e3fSJunchao Zhang    Level: advanced
1663a1729e3fSJunchao Zhang 
1664db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeEnd()`, `PetscSFReduceEnd()`, `PetscSFSetGraph()`
1665a1729e3fSJunchao Zhang @*/
1666*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf, MPI_Datatype unit, void *rootdata, const void *leafdata, void *leafupdate, MPI_Op op)
1667*d71ae5a4SJacob Faibussowitsch {
1668a1729e3fSJunchao Zhang   PetscFunctionBegin;
1669a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
16709566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_FetchAndOpEnd, sf, 0, 0, 0));
1671dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, FetchAndOpEnd, unit, rootdata, leafdata, leafupdate, op);
16729566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_FetchAndOpEnd, sf, 0, 0, 0));
1673a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1674a1729e3fSJunchao Zhang }
1675a1729e3fSJunchao Zhang 
1676a1729e3fSJunchao Zhang /*@C
167795fce210SBarry Smith    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
167895fce210SBarry Smith 
167995fce210SBarry Smith    Collective
168095fce210SBarry Smith 
16814165533cSJose E. Roman    Input Parameter:
168295fce210SBarry Smith .  sf - star forest
168395fce210SBarry Smith 
16844165533cSJose E. Roman    Output Parameter:
168595fce210SBarry Smith .  degree - degree of each root vertex
168695fce210SBarry Smith 
168795fce210SBarry Smith    Level: advanced
168895fce210SBarry Smith 
1689ffe67aa5SVáclav Hapla    Notes:
1690ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1691ffe67aa5SVáclav Hapla 
1692db781477SPatrick Sanan .seealso: `PetscSFGatherBegin()`
169395fce210SBarry Smith @*/
1694*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf, const PetscInt **degree)
1695*d71ae5a4SJacob Faibussowitsch {
169695fce210SBarry Smith   PetscFunctionBegin;
169795fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
169895fce210SBarry Smith   PetscSFCheckGraphSet(sf, 1);
169995fce210SBarry Smith   PetscValidPointer(degree, 2);
1700803bd9e8SMatthew G. Knepley   if (!sf->degreeknown) {
17015b0d146aSStefano Zampini     PetscInt i, nroots = sf->nroots, maxlocal;
170228b400f6SJacob Faibussowitsch     PetscCheck(!sf->degree, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
17035b0d146aSStefano Zampini     maxlocal = sf->maxleaf - sf->minleaf + 1;
17049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nroots, &sf->degree));
17059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PetscMax(maxlocal, 1), &sf->degreetmp)); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
170629046d53SLisandro Dalcin     for (i = 0; i < nroots; i++) sf->degree[i] = 0;
17079837ea96SMatthew G. Knepley     for (i = 0; i < maxlocal; i++) sf->degreetmp[i] = 1;
17089566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf, MPIU_INT, sf->degreetmp - sf->minleaf, sf->degree, MPI_SUM));
170995fce210SBarry Smith   }
171095fce210SBarry Smith   *degree = NULL;
171195fce210SBarry Smith   PetscFunctionReturn(0);
171295fce210SBarry Smith }
171395fce210SBarry Smith 
171495fce210SBarry Smith /*@C
171595fce210SBarry Smith    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
171695fce210SBarry Smith 
171795fce210SBarry Smith    Collective
171895fce210SBarry Smith 
17194165533cSJose E. Roman    Input Parameter:
172095fce210SBarry Smith .  sf - star forest
172195fce210SBarry Smith 
17224165533cSJose E. Roman    Output Parameter:
172395fce210SBarry Smith .  degree - degree of each root vertex
172495fce210SBarry Smith 
172595fce210SBarry Smith    Level: developer
172695fce210SBarry Smith 
1727ffe67aa5SVáclav Hapla    Notes:
1728ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1729ffe67aa5SVáclav Hapla 
173095fce210SBarry Smith .seealso:
173195fce210SBarry Smith @*/
1732*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf, const PetscInt **degree)
1733*d71ae5a4SJacob Faibussowitsch {
173495fce210SBarry Smith   PetscFunctionBegin;
173595fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
173695fce210SBarry Smith   PetscSFCheckGraphSet(sf, 1);
173729046d53SLisandro Dalcin   PetscValidPointer(degree, 2);
173895fce210SBarry Smith   if (!sf->degreeknown) {
173928b400f6SJacob Faibussowitsch     PetscCheck(sf->degreetmp, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
17409566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf, MPIU_INT, sf->degreetmp - sf->minleaf, sf->degree, MPI_SUM));
17419566063dSJacob Faibussowitsch     PetscCall(PetscFree(sf->degreetmp));
174295fce210SBarry Smith     sf->degreeknown = PETSC_TRUE;
174395fce210SBarry Smith   }
174495fce210SBarry Smith   *degree = sf->degree;
174595fce210SBarry Smith   PetscFunctionReturn(0);
174695fce210SBarry Smith }
174795fce210SBarry Smith 
1748673100f5SVaclav Hapla /*@C
174966dfcd1aSVaclav Hapla    PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()).
175066dfcd1aSVaclav Hapla    Each multi-root is assigned index of the corresponding original root.
1751673100f5SVaclav Hapla 
1752673100f5SVaclav Hapla    Collective
1753673100f5SVaclav Hapla 
17544165533cSJose E. Roman    Input Parameters:
1755673100f5SVaclav Hapla +  sf - star forest
1756673100f5SVaclav Hapla -  degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd()
1757673100f5SVaclav Hapla 
17584165533cSJose E. Roman    Output Parameters:
175966dfcd1aSVaclav Hapla +  nMultiRoots - (optional) number of multi-roots (roots of multi-SF)
176066dfcd1aSVaclav Hapla -  multiRootsOrigNumbering - original indices of multi-roots; length of this array is nMultiRoots
1761673100f5SVaclav Hapla 
1762673100f5SVaclav Hapla    Level: developer
1763673100f5SVaclav Hapla 
1764ffe67aa5SVáclav Hapla    Notes:
1765ffe67aa5SVáclav Hapla    The returned array multiRootsOrigNumbering is newly allocated and should be destroyed with PetscFree() when no longer needed.
1766ffe67aa5SVáclav Hapla 
1767db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeBegin()`, `PetscSFComputeDegreeEnd()`, `PetscSFGetMultiSF()`
1768673100f5SVaclav Hapla @*/
1769*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[])
1770*d71ae5a4SJacob Faibussowitsch {
1771673100f5SVaclav Hapla   PetscSF  msf;
1772673100f5SVaclav Hapla   PetscInt i, j, k, nroots, nmroots;
1773673100f5SVaclav Hapla 
1774673100f5SVaclav Hapla   PetscFunctionBegin;
1775673100f5SVaclav Hapla   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
17769566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
177729328920SVaclav Hapla   if (nroots) PetscValidIntPointer(degree, 2);
177866dfcd1aSVaclav Hapla   if (nMultiRoots) PetscValidIntPointer(nMultiRoots, 3);
177966dfcd1aSVaclav Hapla   PetscValidPointer(multiRootsOrigNumbering, 4);
17809566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf, &msf));
17819566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL));
17829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nmroots, multiRootsOrigNumbering));
1783673100f5SVaclav Hapla   for (i = 0, j = 0, k = 0; i < nroots; i++) {
1784673100f5SVaclav Hapla     if (!degree[i]) continue;
1785ad540459SPierre Jolivet     for (j = 0; j < degree[i]; j++, k++) (*multiRootsOrigNumbering)[k] = i;
1786673100f5SVaclav Hapla   }
178708401ef6SPierre Jolivet   PetscCheck(k == nmroots, PETSC_COMM_SELF, PETSC_ERR_PLIB, "sanity check fail");
178866dfcd1aSVaclav Hapla   if (nMultiRoots) *nMultiRoots = nmroots;
1789673100f5SVaclav Hapla   PetscFunctionReturn(0);
1790673100f5SVaclav Hapla }
1791673100f5SVaclav Hapla 
179295fce210SBarry Smith /*@C
179395fce210SBarry Smith    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
179495fce210SBarry Smith 
179595fce210SBarry Smith    Collective
179695fce210SBarry Smith 
17974165533cSJose E. Roman    Input Parameters:
179895fce210SBarry Smith +  sf - star forest
179995fce210SBarry Smith .  unit - data type
180095fce210SBarry Smith -  leafdata - leaf data to gather to roots
180195fce210SBarry Smith 
18024165533cSJose E. Roman    Output Parameter:
180395fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
180495fce210SBarry Smith 
180595fce210SBarry Smith    Level: intermediate
180695fce210SBarry Smith 
1807db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeBegin()`, `PetscSFScatterBegin()`
180895fce210SBarry Smith @*/
1809*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGatherBegin(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *multirootdata)
1810*d71ae5a4SJacob Faibussowitsch {
1811a5526d50SJunchao Zhang   PetscSF multi = NULL;
181295fce210SBarry Smith 
181395fce210SBarry Smith   PetscFunctionBegin;
181495fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
18159566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
18169566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf, &multi));
18179566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(multi, unit, leafdata, multirootdata, MPI_REPLACE));
181895fce210SBarry Smith   PetscFunctionReturn(0);
181995fce210SBarry Smith }
182095fce210SBarry Smith 
182195fce210SBarry Smith /*@C
182295fce210SBarry Smith    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
182395fce210SBarry Smith 
182495fce210SBarry Smith    Collective
182595fce210SBarry Smith 
18264165533cSJose E. Roman    Input Parameters:
182795fce210SBarry Smith +  sf - star forest
182895fce210SBarry Smith .  unit - data type
182995fce210SBarry Smith -  leafdata - leaf data to gather to roots
183095fce210SBarry Smith 
18314165533cSJose E. Roman    Output Parameter:
183295fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
183395fce210SBarry Smith 
183495fce210SBarry Smith    Level: intermediate
183595fce210SBarry Smith 
1836db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeEnd()`, `PetscSFScatterEnd()`
183795fce210SBarry Smith @*/
1838*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGatherEnd(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *multirootdata)
1839*d71ae5a4SJacob Faibussowitsch {
1840a5526d50SJunchao Zhang   PetscSF multi = NULL;
184195fce210SBarry Smith 
184295fce210SBarry Smith   PetscFunctionBegin;
184395fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
18449566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf, &multi));
18459566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(multi, unit, leafdata, multirootdata, MPI_REPLACE));
184695fce210SBarry Smith   PetscFunctionReturn(0);
184795fce210SBarry Smith }
184895fce210SBarry Smith 
184995fce210SBarry Smith /*@C
185095fce210SBarry Smith    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
185195fce210SBarry Smith 
185295fce210SBarry Smith    Collective
185395fce210SBarry Smith 
18544165533cSJose E. Roman    Input Parameters:
185595fce210SBarry Smith +  sf - star forest
185695fce210SBarry Smith .  unit - data type
185795fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
185895fce210SBarry Smith 
18594165533cSJose E. Roman    Output Parameter:
186095fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
186195fce210SBarry Smith 
186295fce210SBarry Smith    Level: intermediate
186395fce210SBarry Smith 
1864db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeBegin()`, `PetscSFScatterBegin()`
186595fce210SBarry Smith @*/
1866*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFScatterBegin(PetscSF sf, MPI_Datatype unit, const void *multirootdata, void *leafdata)
1867*d71ae5a4SJacob Faibussowitsch {
1868a5526d50SJunchao Zhang   PetscSF multi = NULL;
186995fce210SBarry Smith 
187095fce210SBarry Smith   PetscFunctionBegin;
187195fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
18729566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
18739566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf, &multi));
18749566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(multi, unit, multirootdata, leafdata, MPI_REPLACE));
187595fce210SBarry Smith   PetscFunctionReturn(0);
187695fce210SBarry Smith }
187795fce210SBarry Smith 
187895fce210SBarry Smith /*@C
187995fce210SBarry Smith    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
188095fce210SBarry Smith 
188195fce210SBarry Smith    Collective
188295fce210SBarry Smith 
18834165533cSJose E. Roman    Input Parameters:
188495fce210SBarry Smith +  sf - star forest
188595fce210SBarry Smith .  unit - data type
188695fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
188795fce210SBarry Smith 
18884165533cSJose E. Roman    Output Parameter:
188995fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
189095fce210SBarry Smith 
189195fce210SBarry Smith    Level: intermediate
189295fce210SBarry Smith 
1893db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeEnd()`, `PetscSFScatterEnd()`
189495fce210SBarry Smith @*/
1895*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFScatterEnd(PetscSF sf, MPI_Datatype unit, const void *multirootdata, void *leafdata)
1896*d71ae5a4SJacob Faibussowitsch {
1897a5526d50SJunchao Zhang   PetscSF multi = NULL;
189895fce210SBarry Smith 
189995fce210SBarry Smith   PetscFunctionBegin;
190095fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
19019566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf, &multi));
19029566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(multi, unit, multirootdata, leafdata, MPI_REPLACE));
190395fce210SBarry Smith   PetscFunctionReturn(0);
190495fce210SBarry Smith }
1905a7b3aa13SAta Mesgarnejad 
1906*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFCheckLeavesUnique_Private(PetscSF sf)
1907*d71ae5a4SJacob Faibussowitsch {
1908a072220fSLawrence Mitchell   PetscInt        i, n, nleaves;
1909a072220fSLawrence Mitchell   const PetscInt *ilocal = NULL;
1910a072220fSLawrence Mitchell   PetscHSetI      seen;
1911a072220fSLawrence Mitchell 
1912a072220fSLawrence Mitchell   PetscFunctionBegin;
1913b458e8f1SJose E. Roman   if (PetscDefined(USE_DEBUG)) {
19149566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, &ilocal, NULL));
19159566063dSJacob Faibussowitsch     PetscCall(PetscHSetICreate(&seen));
1916a072220fSLawrence Mitchell     for (i = 0; i < nleaves; i++) {
1917a072220fSLawrence Mitchell       const PetscInt leaf = ilocal ? ilocal[i] : i;
19189566063dSJacob Faibussowitsch       PetscCall(PetscHSetIAdd(seen, leaf));
1919a072220fSLawrence Mitchell     }
19209566063dSJacob Faibussowitsch     PetscCall(PetscHSetIGetSize(seen, &n));
192108401ef6SPierre Jolivet     PetscCheck(n == nleaves, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Provided leaves have repeated values: all leaves must be unique");
19229566063dSJacob Faibussowitsch     PetscCall(PetscHSetIDestroy(&seen));
1923b458e8f1SJose E. Roman   }
1924a072220fSLawrence Mitchell   PetscFunctionReturn(0);
1925a072220fSLawrence Mitchell }
192654729392SStefano Zampini 
1927a7b3aa13SAta Mesgarnejad /*@
192804c0ada0SJunchao Zhang   PetscSFCompose - Compose a new PetscSF by putting the second SF under the first one in a top (roots) down (leaves) view
1929a7b3aa13SAta Mesgarnejad 
1930a7b3aa13SAta Mesgarnejad   Input Parameters:
193154729392SStefano Zampini + sfA - The first PetscSF
193254729392SStefano Zampini - sfB - The second PetscSF
1933a7b3aa13SAta Mesgarnejad 
1934a7b3aa13SAta Mesgarnejad   Output Parameters:
193554729392SStefano Zampini . sfBA - The composite SF
1936a7b3aa13SAta Mesgarnejad 
1937a7b3aa13SAta Mesgarnejad   Level: developer
1938a7b3aa13SAta Mesgarnejad 
1939a072220fSLawrence Mitchell   Notes:
194054729392SStefano Zampini   Currently, the two SFs must be defined on congruent communicators and they must be true star
194154729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots.
194254729392SStefano Zampini 
194354729392SStefano Zampini   sfA's leaf space and sfB's root space might be partially overlapped. The composition builds
194454729392SStefano Zampini   a graph with sfA's roots and sfB's leaves only when there is a path between them. Unconnected
194554729392SStefano Zampini   nodes (roots or leaves) are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a
194654729392SStefano Zampini   Bcast on sfA, then a Bcast on sfB, on connected nodes.
1947a072220fSLawrence Mitchell 
1948db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFComposeInverse()`, `PetscSFGetGraph()`, `PetscSFSetGraph()`
1949a7b3aa13SAta Mesgarnejad @*/
1950*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA)
1951*d71ae5a4SJacob Faibussowitsch {
1952a7b3aa13SAta Mesgarnejad   const PetscSFNode *remotePointsA, *remotePointsB;
1953d41018fbSJunchao Zhang   PetscSFNode       *remotePointsBA = NULL, *reorderedRemotePointsA = NULL, *leafdataB;
195454729392SStefano Zampini   const PetscInt    *localPointsA, *localPointsB;
195554729392SStefano Zampini   PetscInt          *localPointsBA;
195654729392SStefano Zampini   PetscInt           i, numRootsA, numLeavesA, numRootsB, numLeavesB, minleaf, maxleaf, numLeavesBA;
195754729392SStefano Zampini   PetscBool          denseB;
1958a7b3aa13SAta Mesgarnejad 
1959a7b3aa13SAta Mesgarnejad   PetscFunctionBegin;
1960a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1);
196129046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfA, 1);
196229046d53SLisandro Dalcin   PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 2);
196329046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfB, 2);
196454729392SStefano Zampini   PetscCheckSameComm(sfA, 1, sfB, 2);
196529046d53SLisandro Dalcin   PetscValidPointer(sfBA, 3);
19669566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckLeavesUnique_Private(sfA));
19679566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckLeavesUnique_Private(sfB));
196854729392SStefano Zampini 
19699566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA));
19709566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB));
19710ea77edaSksagiyam   /* Make sure that PetscSFBcast{Begin, End}(sfB, ...) works with root data of size       */
19720ea77edaSksagiyam   /* numRootsB; otherwise, garbage will be broadcasted.                                   */
19730ea77edaSksagiyam   /* Example (comm size = 1):                                                             */
19740ea77edaSksagiyam   /* sfA: 0 <- (0, 0)                                                                     */
19750ea77edaSksagiyam   /* sfB: 100 <- (0, 0)                                                                   */
19760ea77edaSksagiyam   /*      101 <- (0, 1)                                                                   */
19770ea77edaSksagiyam   /* Here, we have remotePointsA = [(0, 0)], but for remotePointsA to be a valid tartget  */
19780ea77edaSksagiyam   /* of sfB, it has to be recasted as [(0, 0), (-1, -1)] so that points 100 and 101 would */
19790ea77edaSksagiyam   /* receive (0, 0) and (-1, -1), respectively, when PetscSFBcast(sfB, ...) is called on  */
19800ea77edaSksagiyam   /* remotePointsA; if not recasted, point 101 would receive a garbage value.             */
19819566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numRootsB, &reorderedRemotePointsA));
198254729392SStefano Zampini   for (i = 0; i < numRootsB; i++) {
198354729392SStefano Zampini     reorderedRemotePointsA[i].rank  = -1;
198454729392SStefano Zampini     reorderedRemotePointsA[i].index = -1;
198554729392SStefano Zampini   }
198654729392SStefano Zampini   for (i = 0; i < numLeavesA; i++) {
19870ea77edaSksagiyam     PetscInt localp = localPointsA ? localPointsA[i] : i;
19880ea77edaSksagiyam 
19890ea77edaSksagiyam     if (localp >= numRootsB) continue;
19900ea77edaSksagiyam     reorderedRemotePointsA[localp] = remotePointsA[i];
199154729392SStefano Zampini   }
1992d41018fbSJunchao Zhang   remotePointsA = reorderedRemotePointsA;
19939566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(sfB, &minleaf, &maxleaf));
19949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxleaf - minleaf + 1, &leafdataB));
19950ea77edaSksagiyam   for (i = 0; i < maxleaf - minleaf + 1; i++) {
19960ea77edaSksagiyam     leafdataB[i].rank  = -1;
19970ea77edaSksagiyam     leafdataB[i].index = -1;
19980ea77edaSksagiyam   }
19999566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, leafdataB - minleaf, MPI_REPLACE));
20009566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, leafdataB - minleaf, MPI_REPLACE));
20019566063dSJacob Faibussowitsch   PetscCall(PetscFree(reorderedRemotePointsA));
2002d41018fbSJunchao Zhang 
200354729392SStefano Zampini   denseB = (PetscBool)!localPointsB;
200454729392SStefano Zampini   for (i = 0, numLeavesBA = 0; i < numLeavesB; i++) {
200554729392SStefano Zampini     if (leafdataB[localPointsB ? localPointsB[i] - minleaf : i].rank == -1) denseB = PETSC_FALSE;
200654729392SStefano Zampini     else numLeavesBA++;
200754729392SStefano Zampini   }
200854729392SStefano Zampini   if (denseB) {
2009d41018fbSJunchao Zhang     localPointsBA  = NULL;
2010d41018fbSJunchao Zhang     remotePointsBA = leafdataB;
2011d41018fbSJunchao Zhang   } else {
20129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numLeavesBA, &localPointsBA));
20139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numLeavesBA, &remotePointsBA));
201454729392SStefano Zampini     for (i = 0, numLeavesBA = 0; i < numLeavesB; i++) {
201554729392SStefano Zampini       const PetscInt l = localPointsB ? localPointsB[i] : i;
201654729392SStefano Zampini 
201754729392SStefano Zampini       if (leafdataB[l - minleaf].rank == -1) continue;
201854729392SStefano Zampini       remotePointsBA[numLeavesBA] = leafdataB[l - minleaf];
201954729392SStefano Zampini       localPointsBA[numLeavesBA]  = l;
202054729392SStefano Zampini       numLeavesBA++;
202154729392SStefano Zampini     }
20229566063dSJacob Faibussowitsch     PetscCall(PetscFree(leafdataB));
2023d41018fbSJunchao Zhang   }
20249566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfA), sfBA));
20259566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(*sfBA));
20269566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sfBA, numRootsA, numLeavesBA, localPointsBA, PETSC_OWN_POINTER, remotePointsBA, PETSC_OWN_POINTER));
2027a7b3aa13SAta Mesgarnejad   PetscFunctionReturn(0);
2028a7b3aa13SAta Mesgarnejad }
20291c6ba672SJunchao Zhang 
203004c0ada0SJunchao Zhang /*@
203154729392SStefano Zampini   PetscSFComposeInverse - Compose a new PetscSF by putting the inverse of the second SF under the first one
203204c0ada0SJunchao Zhang 
203304c0ada0SJunchao Zhang   Input Parameters:
203454729392SStefano Zampini + sfA - The first PetscSF
203554729392SStefano Zampini - sfB - The second PetscSF
203604c0ada0SJunchao Zhang 
203704c0ada0SJunchao Zhang   Output Parameters:
203854729392SStefano Zampini . sfBA - The composite SF.
203904c0ada0SJunchao Zhang 
204004c0ada0SJunchao Zhang   Level: developer
204104c0ada0SJunchao Zhang 
204254729392SStefano Zampini   Notes:
204354729392SStefano Zampini   Currently, the two SFs must be defined on congruent communicators and they must be true star
204454729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots. Even more, all roots of the
204554729392SStefano Zampini   second SF must have a degree of 1, i.e., no roots have more than one leaf connected.
204654729392SStefano Zampini 
204754729392SStefano Zampini   sfA's leaf space and sfB's leaf space might be partially overlapped. The composition builds
204854729392SStefano Zampini   a graph with sfA's roots and sfB's roots only when there is a path between them. Unconnected
204954729392SStefano Zampini   roots are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a Bcast on sfA, then
205054729392SStefano Zampini   a Reduce on sfB, on connected roots.
205154729392SStefano Zampini 
2052db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFCompose()`, `PetscSFGetGraph()`, `PetscSFSetGraph()`, `PetscSFCreateInverseSF()`
205304c0ada0SJunchao Zhang @*/
2054*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFComposeInverse(PetscSF sfA, PetscSF sfB, PetscSF *sfBA)
2055*d71ae5a4SJacob Faibussowitsch {
205604c0ada0SJunchao Zhang   const PetscSFNode *remotePointsA, *remotePointsB;
205704c0ada0SJunchao Zhang   PetscSFNode       *remotePointsBA;
205804c0ada0SJunchao Zhang   const PetscInt    *localPointsA, *localPointsB;
205954729392SStefano Zampini   PetscSFNode       *reorderedRemotePointsA = NULL;
206054729392SStefano Zampini   PetscInt           i, numRootsA, numLeavesA, numLeavesBA, numRootsB, numLeavesB, minleaf, maxleaf, *localPointsBA;
20615b0d146aSStefano Zampini   MPI_Op             op;
20625b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
20635b0d146aSStefano Zampini   PetscBool iswin;
20645b0d146aSStefano Zampini #endif
206504c0ada0SJunchao Zhang 
206604c0ada0SJunchao Zhang   PetscFunctionBegin;
206704c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1);
206804c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfA, 1);
206904c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 2);
207004c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfB, 2);
207154729392SStefano Zampini   PetscCheckSameComm(sfA, 1, sfB, 2);
207204c0ada0SJunchao Zhang   PetscValidPointer(sfBA, 3);
20739566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckLeavesUnique_Private(sfA));
20749566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckLeavesUnique_Private(sfB));
207554729392SStefano Zampini 
20769566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA));
20779566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB));
20785b0d146aSStefano Zampini 
20795b0d146aSStefano Zampini   /* TODO: Check roots of sfB have degree of 1 */
20805b0d146aSStefano Zampini   /* Once we implement it, we can replace the MPI_MAXLOC
208183df288dSJunchao Zhang      with MPI_REPLACE. In that case, MPI_MAXLOC and MPI_REPLACE have the same effect.
20825b0d146aSStefano Zampini      We use MPI_MAXLOC only to have a deterministic output from this routine if
20835b0d146aSStefano Zampini      the root condition is not meet.
20845b0d146aSStefano Zampini    */
20855b0d146aSStefano Zampini   op = MPI_MAXLOC;
20865b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
20875b0d146aSStefano Zampini   /* we accept a non-deterministic output (if any) with PETSCSFWINDOW, since MPI_MAXLOC cannot operate on MPIU_2INT with MPI_Accumulate */
20889566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)sfB, PETSCSFWINDOW, &iswin));
208983df288dSJunchao Zhang   if (iswin) op = MPI_REPLACE;
20905b0d146aSStefano Zampini #endif
20915b0d146aSStefano Zampini 
20929566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(sfB, &minleaf, &maxleaf));
20939566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxleaf - minleaf + 1, &reorderedRemotePointsA));
209454729392SStefano Zampini   for (i = 0; i < maxleaf - minleaf + 1; i++) {
209554729392SStefano Zampini     reorderedRemotePointsA[i].rank  = -1;
209654729392SStefano Zampini     reorderedRemotePointsA[i].index = -1;
209754729392SStefano Zampini   }
209854729392SStefano Zampini   if (localPointsA) {
209954729392SStefano Zampini     for (i = 0; i < numLeavesA; i++) {
210054729392SStefano Zampini       if (localPointsA[i] > maxleaf || localPointsA[i] < minleaf) continue;
210154729392SStefano Zampini       reorderedRemotePointsA[localPointsA[i] - minleaf] = remotePointsA[i];
210254729392SStefano Zampini     }
210354729392SStefano Zampini   } else {
210454729392SStefano Zampini     for (i = 0; i < numLeavesA; i++) {
210554729392SStefano Zampini       if (i > maxleaf || i < minleaf) continue;
210654729392SStefano Zampini       reorderedRemotePointsA[i - minleaf] = remotePointsA[i];
210754729392SStefano Zampini     }
210854729392SStefano Zampini   }
210954729392SStefano Zampini 
21109566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numRootsB, &localPointsBA));
21119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numRootsB, &remotePointsBA));
211254729392SStefano Zampini   for (i = 0; i < numRootsB; i++) {
211354729392SStefano Zampini     remotePointsBA[i].rank  = -1;
211454729392SStefano Zampini     remotePointsBA[i].index = -1;
211554729392SStefano Zampini   }
211654729392SStefano Zampini 
21179566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sfB, MPIU_2INT, reorderedRemotePointsA - minleaf, remotePointsBA, op));
21189566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sfB, MPIU_2INT, reorderedRemotePointsA - minleaf, remotePointsBA, op));
21199566063dSJacob Faibussowitsch   PetscCall(PetscFree(reorderedRemotePointsA));
212054729392SStefano Zampini   for (i = 0, numLeavesBA = 0; i < numRootsB; i++) {
212154729392SStefano Zampini     if (remotePointsBA[i].rank == -1) continue;
212254729392SStefano Zampini     remotePointsBA[numLeavesBA].rank  = remotePointsBA[i].rank;
212354729392SStefano Zampini     remotePointsBA[numLeavesBA].index = remotePointsBA[i].index;
212454729392SStefano Zampini     localPointsBA[numLeavesBA]        = i;
212554729392SStefano Zampini     numLeavesBA++;
212654729392SStefano Zampini   }
21279566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfA), sfBA));
21289566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(*sfBA));
21299566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sfBA, numRootsA, numLeavesBA, localPointsBA, PETSC_OWN_POINTER, remotePointsBA, PETSC_OWN_POINTER));
213004c0ada0SJunchao Zhang   PetscFunctionReturn(0);
213104c0ada0SJunchao Zhang }
213204c0ada0SJunchao Zhang 
21331c6ba672SJunchao Zhang /*
21341c6ba672SJunchao Zhang   PetscSFCreateLocalSF_Private - Creates a local PetscSF that only has intra-process edges of the global PetscSF
21351c6ba672SJunchao Zhang 
21361c6ba672SJunchao Zhang   Input Parameters:
21371c6ba672SJunchao Zhang . sf - The global PetscSF
21381c6ba672SJunchao Zhang 
21391c6ba672SJunchao Zhang   Output Parameters:
21401c6ba672SJunchao Zhang . out - The local PetscSF
21411c6ba672SJunchao Zhang  */
2142*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateLocalSF_Private(PetscSF sf, PetscSF *out)
2143*d71ae5a4SJacob Faibussowitsch {
21441c6ba672SJunchao Zhang   MPI_Comm           comm;
21451c6ba672SJunchao Zhang   PetscMPIInt        myrank;
21461c6ba672SJunchao Zhang   const PetscInt    *ilocal;
21471c6ba672SJunchao Zhang   const PetscSFNode *iremote;
21481c6ba672SJunchao Zhang   PetscInt           i, j, nroots, nleaves, lnleaves, *lilocal;
21491c6ba672SJunchao Zhang   PetscSFNode       *liremote;
21501c6ba672SJunchao Zhang   PetscSF            lsf;
21511c6ba672SJunchao Zhang 
21521c6ba672SJunchao Zhang   PetscFunctionBegin;
21531c6ba672SJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
2154dbbe0bcdSBarry Smith   if (sf->ops->CreateLocalSF) PetscUseTypeMethod(sf, CreateLocalSF, out);
2155dbbe0bcdSBarry Smith   else {
21561c6ba672SJunchao Zhang     /* Could use PetscSFCreateEmbeddedLeafSF, but since we know the comm is PETSC_COMM_SELF, we can make it fast */
21579566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
21589566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &myrank));
21591c6ba672SJunchao Zhang 
21601c6ba672SJunchao Zhang     /* Find out local edges and build a local SF */
21619566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, &iremote));
21629371c9d4SSatish Balay     for (i = lnleaves = 0; i < nleaves; i++) {
21639371c9d4SSatish Balay       if (iremote[i].rank == (PetscInt)myrank) lnleaves++;
21649371c9d4SSatish Balay     }
21659566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lnleaves, &lilocal));
21669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lnleaves, &liremote));
21671c6ba672SJunchao Zhang 
21681c6ba672SJunchao Zhang     for (i = j = 0; i < nleaves; i++) {
21691c6ba672SJunchao Zhang       if (iremote[i].rank == (PetscInt)myrank) {
21701c6ba672SJunchao Zhang         lilocal[j]        = ilocal ? ilocal[i] : i; /* ilocal=NULL for contiguous storage */
21711c6ba672SJunchao Zhang         liremote[j].rank  = 0;                      /* rank in PETSC_COMM_SELF */
21721c6ba672SJunchao Zhang         liremote[j].index = iremote[i].index;
21731c6ba672SJunchao Zhang         j++;
21741c6ba672SJunchao Zhang       }
21751c6ba672SJunchao Zhang     }
21769566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PETSC_COMM_SELF, &lsf));
21779566063dSJacob Faibussowitsch     PetscCall(PetscSFSetFromOptions(lsf));
21789566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(lsf, nroots, lnleaves, lilocal, PETSC_OWN_POINTER, liremote, PETSC_OWN_POINTER));
21799566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(lsf));
21801c6ba672SJunchao Zhang     *out = lsf;
21811c6ba672SJunchao Zhang   }
21821c6ba672SJunchao Zhang   PetscFunctionReturn(0);
21831c6ba672SJunchao Zhang }
2184dd5b3ca6SJunchao Zhang 
2185dd5b3ca6SJunchao Zhang /* Similar to PetscSFBcast, but only Bcast to leaves on rank 0 */
2186*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFBcastToZero_Private(PetscSF sf, MPI_Datatype unit, const void *rootdata, void *leafdata)
2187*d71ae5a4SJacob Faibussowitsch {
2188eb02082bSJunchao Zhang   PetscMemType rootmtype, leafmtype;
2189dd5b3ca6SJunchao Zhang 
2190dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
2191dd5b3ca6SJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
21929566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
21939566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_BcastBegin, sf, 0, 0, 0));
21949566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(rootdata, &rootmtype));
21959566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafdata, &leafmtype));
2196dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, BcastToZero, unit, rootmtype, rootdata, leafmtype, leafdata);
21979566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_BcastBegin, sf, 0, 0, 0));
2198dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
2199dd5b3ca6SJunchao Zhang }
2200dd5b3ca6SJunchao Zhang 
2201157edd7aSVaclav Hapla /*@
2202157edd7aSVaclav Hapla   PetscSFConcatenate - concatenate multiple SFs into one
2203157edd7aSVaclav Hapla 
2204157edd7aSVaclav Hapla   Input Parameters:
2205157edd7aSVaclav Hapla + comm - the communicator
2206157edd7aSVaclav Hapla . nsfs - the number of input PetscSF
2207157edd7aSVaclav Hapla . sfs  - the array of input PetscSF
2208157edd7aSVaclav Hapla . shareRoots - the flag whether roots of input PetscSFs are taken as shared (PETSC_TRUE), or separate and concatenated (PETSC_FALSE)
2209157edd7aSVaclav Hapla - leafOffsets - the array of local leaf offsets, one for each input PetscSF, or NULL for contiguous storage
2210157edd7aSVaclav Hapla 
2211157edd7aSVaclav Hapla   Output Parameters:
2212157edd7aSVaclav Hapla . newsf - The resulting PetscSF
2213157edd7aSVaclav Hapla 
2214157edd7aSVaclav Hapla   Level: developer
2215157edd7aSVaclav Hapla 
2216157edd7aSVaclav Hapla   Notes:
2217157edd7aSVaclav Hapla   The communicator of all SFs in sfs must be comm.
2218157edd7aSVaclav Hapla 
2219157edd7aSVaclav Hapla   The offsets in leafOffsets are added to the original leaf indices.
2220157edd7aSVaclav Hapla 
2221157edd7aSVaclav Hapla   If all input SFs use contiguous leaf storage (ilocal = NULL), leafOffsets can be passed as NULL as well.
2222157edd7aSVaclav Hapla   In this case, NULL is also passed as ilocal to the resulting SF.
2223157edd7aSVaclav Hapla 
2224157edd7aSVaclav Hapla   If any input SF has non-null ilocal, leafOffsets is needed to distinguish leaves from different input SFs.
2225157edd7aSVaclav Hapla   In this case, user is responsible to provide correct offsets so that the resulting leaves are unique (otherwise an error occurs).
2226157edd7aSVaclav Hapla 
2227db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFCompose()`, `PetscSFGetGraph()`, `PetscSFSetGraph()`
2228157edd7aSVaclav Hapla @*/
2229*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFConcatenate(MPI_Comm comm, PetscInt nsfs, PetscSF sfs[], PetscBool shareRoots, PetscInt leafOffsets[], PetscSF *newsf)
2230*d71ae5a4SJacob Faibussowitsch {
2231157edd7aSVaclav Hapla   PetscInt     i, s, nLeaves, nRoots;
2232157edd7aSVaclav Hapla   PetscInt    *leafArrayOffsets;
2233157edd7aSVaclav Hapla   PetscInt    *ilocal_new;
2234157edd7aSVaclav Hapla   PetscSFNode *iremote_new;
2235157edd7aSVaclav Hapla   PetscInt    *rootOffsets;
2236157edd7aSVaclav Hapla   PetscBool    all_ilocal_null = PETSC_FALSE;
2237157edd7aSVaclav Hapla   PetscMPIInt  rank;
2238157edd7aSVaclav Hapla 
2239157edd7aSVaclav Hapla   PetscFunctionBegin;
2240157edd7aSVaclav Hapla   {
2241157edd7aSVaclav Hapla     PetscSF dummy; /* just to have a PetscObject on comm for input validation */
2242157edd7aSVaclav Hapla 
22439566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(comm, &dummy));
2244157edd7aSVaclav Hapla     PetscValidLogicalCollectiveInt(dummy, nsfs, 2);
2245157edd7aSVaclav Hapla     PetscValidPointer(sfs, 3);
2246157edd7aSVaclav Hapla     for (i = 0; i < nsfs; i++) {
2247157edd7aSVaclav Hapla       PetscValidHeaderSpecific(sfs[i], PETSCSF_CLASSID, 3);
2248157edd7aSVaclav Hapla       PetscCheckSameComm(dummy, 1, sfs[i], 3);
2249157edd7aSVaclav Hapla     }
2250157edd7aSVaclav Hapla     PetscValidLogicalCollectiveBool(dummy, shareRoots, 4);
2251157edd7aSVaclav Hapla     if (leafOffsets) PetscValidIntPointer(leafOffsets, 5);
2252157edd7aSVaclav Hapla     PetscValidPointer(newsf, 6);
22539566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&dummy));
2254157edd7aSVaclav Hapla   }
2255157edd7aSVaclav Hapla   if (!nsfs) {
22569566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(comm, newsf));
22579566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*newsf, 0, 0, NULL, PETSC_OWN_POINTER, NULL, PETSC_OWN_POINTER));
2258157edd7aSVaclav Hapla     PetscFunctionReturn(0);
2259157edd7aSVaclav Hapla   }
22609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2261157edd7aSVaclav Hapla 
22629566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nsfs + 1, &rootOffsets));
2263157edd7aSVaclav Hapla   if (shareRoots) {
22649566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfs[0], &nRoots, NULL, NULL, NULL));
2265157edd7aSVaclav Hapla     if (PetscDefined(USE_DEBUG)) {
2266157edd7aSVaclav Hapla       for (s = 1; s < nsfs; s++) {
2267157edd7aSVaclav Hapla         PetscInt nr;
2268157edd7aSVaclav Hapla 
22699566063dSJacob Faibussowitsch         PetscCall(PetscSFGetGraph(sfs[s], &nr, NULL, NULL, NULL));
2270157edd7aSVaclav Hapla         PetscCheck(nr == nRoots, comm, PETSC_ERR_ARG_SIZ, "shareRoots = PETSC_TRUE but sfs[%" PetscInt_FMT "] has a different number of roots (%" PetscInt_FMT ") than sfs[0] (%" PetscInt_FMT ")", s, nr, nRoots);
2271157edd7aSVaclav Hapla       }
2272157edd7aSVaclav Hapla     }
2273157edd7aSVaclav Hapla   } else {
2274157edd7aSVaclav Hapla     for (s = 0; s < nsfs; s++) {
2275157edd7aSVaclav Hapla       PetscInt nr;
2276157edd7aSVaclav Hapla 
22779566063dSJacob Faibussowitsch       PetscCall(PetscSFGetGraph(sfs[s], &nr, NULL, NULL, NULL));
2278157edd7aSVaclav Hapla       rootOffsets[s + 1] = rootOffsets[s] + nr;
2279157edd7aSVaclav Hapla     }
2280157edd7aSVaclav Hapla     nRoots = rootOffsets[nsfs];
2281157edd7aSVaclav Hapla   }
2282157edd7aSVaclav Hapla 
2283157edd7aSVaclav Hapla   /* Calculate leaf array offsets and automatic root offsets */
22849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsfs + 1, &leafArrayOffsets));
2285157edd7aSVaclav Hapla   leafArrayOffsets[0] = 0;
2286157edd7aSVaclav Hapla   for (s = 0; s < nsfs; s++) {
2287157edd7aSVaclav Hapla     PetscInt nl;
2288157edd7aSVaclav Hapla 
22899566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfs[s], NULL, &nl, NULL, NULL));
2290157edd7aSVaclav Hapla     leafArrayOffsets[s + 1] = leafArrayOffsets[s] + nl;
2291157edd7aSVaclav Hapla   }
2292157edd7aSVaclav Hapla   nLeaves = leafArrayOffsets[nsfs];
2293157edd7aSVaclav Hapla 
2294157edd7aSVaclav Hapla   if (!leafOffsets) {
2295157edd7aSVaclav Hapla     all_ilocal_null = PETSC_TRUE;
2296157edd7aSVaclav Hapla     for (s = 0; s < nsfs; s++) {
2297157edd7aSVaclav Hapla       const PetscInt *ilocal;
2298157edd7aSVaclav Hapla 
22999566063dSJacob Faibussowitsch       PetscCall(PetscSFGetGraph(sfs[s], NULL, NULL, &ilocal, NULL));
2300157edd7aSVaclav Hapla       if (ilocal) {
2301157edd7aSVaclav Hapla         all_ilocal_null = PETSC_FALSE;
2302157edd7aSVaclav Hapla         break;
2303157edd7aSVaclav Hapla       }
2304157edd7aSVaclav Hapla     }
2305157edd7aSVaclav Hapla     PetscCheck(all_ilocal_null, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "leafOffsets can be passed as NULL only if all SFs have ilocal = NULL");
2306157edd7aSVaclav Hapla   }
2307157edd7aSVaclav Hapla 
2308157edd7aSVaclav Hapla   /* Renumber and concatenate local leaves */
2309157edd7aSVaclav Hapla   ilocal_new = NULL;
2310157edd7aSVaclav Hapla   if (!all_ilocal_null) {
23119566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nLeaves, &ilocal_new));
2312157edd7aSVaclav Hapla     for (i = 0; i < nLeaves; i++) ilocal_new[i] = -1;
2313157edd7aSVaclav Hapla     for (s = 0; s < nsfs; s++) {
2314157edd7aSVaclav Hapla       const PetscInt *ilocal;
2315157edd7aSVaclav Hapla       PetscInt       *ilocal_l = &ilocal_new[leafArrayOffsets[s]];
2316157edd7aSVaclav Hapla       PetscInt        i, nleaves_l;
2317157edd7aSVaclav Hapla 
23189566063dSJacob Faibussowitsch       PetscCall(PetscSFGetGraph(sfs[s], NULL, &nleaves_l, &ilocal, NULL));
2319157edd7aSVaclav Hapla       for (i = 0; i < nleaves_l; i++) ilocal_l[i] = (ilocal ? ilocal[i] : i) + leafOffsets[s];
2320157edd7aSVaclav Hapla     }
2321157edd7aSVaclav Hapla   }
2322157edd7aSVaclav Hapla 
2323157edd7aSVaclav Hapla   /* Renumber and concatenate remote roots */
23249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nLeaves, &iremote_new));
2325157edd7aSVaclav Hapla   for (i = 0; i < nLeaves; i++) {
2326157edd7aSVaclav Hapla     iremote_new[i].rank  = -1;
2327157edd7aSVaclav Hapla     iremote_new[i].index = -1;
2328157edd7aSVaclav Hapla   }
2329157edd7aSVaclav Hapla   for (s = 0; s < nsfs; s++) {
2330157edd7aSVaclav Hapla     PetscInt           i, nl, nr;
2331157edd7aSVaclav Hapla     PetscSF            tmp_sf;
2332157edd7aSVaclav Hapla     const PetscSFNode *iremote;
2333157edd7aSVaclav Hapla     PetscSFNode       *tmp_rootdata;
2334157edd7aSVaclav Hapla     PetscSFNode       *tmp_leafdata = &iremote_new[leafArrayOffsets[s]];
2335157edd7aSVaclav Hapla 
23369566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfs[s], &nr, &nl, NULL, &iremote));
23379566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(comm, &tmp_sf));
2338157edd7aSVaclav Hapla     /* create helper SF with contiguous leaves */
23399566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(tmp_sf, nr, nl, NULL, PETSC_USE_POINTER, (PetscSFNode *)iremote, PETSC_COPY_VALUES));
23409566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(tmp_sf));
23419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nr, &tmp_rootdata));
2342157edd7aSVaclav Hapla     for (i = 0; i < nr; i++) {
2343157edd7aSVaclav Hapla       tmp_rootdata[i].index = i + rootOffsets[s];
2344157edd7aSVaclav Hapla       tmp_rootdata[i].rank  = (PetscInt)rank;
2345157edd7aSVaclav Hapla     }
23469566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(tmp_sf, MPIU_2INT, tmp_rootdata, tmp_leafdata, MPI_REPLACE));
23479566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(tmp_sf, MPIU_2INT, tmp_rootdata, tmp_leafdata, MPI_REPLACE));
23489566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&tmp_sf));
23499566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp_rootdata));
2350157edd7aSVaclav Hapla   }
2351157edd7aSVaclav Hapla 
2352157edd7aSVaclav Hapla   /* Build the new SF */
23539566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, newsf));
23549566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*newsf, nRoots, nLeaves, ilocal_new, PETSC_OWN_POINTER, iremote_new, PETSC_OWN_POINTER));
23559566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(*newsf));
23569566063dSJacob Faibussowitsch   PetscCall(PetscFree(rootOffsets));
23579566063dSJacob Faibussowitsch   PetscCall(PetscFree(leafArrayOffsets));
2358157edd7aSVaclav Hapla   PetscFunctionReturn(0);
2359157edd7aSVaclav Hapla }
2360