xref: /petsc/src/vec/is/sf/interface/sf.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
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>
4eec179cfSJacob Faibussowitsch #include <petsc/private/hashmapi.h>
595fce210SBarry Smith 
67fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
77fd2d3dbSJunchao Zhang   #include <cuda_runtime.h>
8715b587bSJunchao Zhang   #include <petscdevice_cuda.h>
97fd2d3dbSJunchao Zhang #endif
107fd2d3dbSJunchao Zhang 
117fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_HIP)
127fd2d3dbSJunchao Zhang   #include <hip/hip_runtime.h>
137fd2d3dbSJunchao Zhang #endif
147fd2d3dbSJunchao Zhang 
152abc8c78SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
164bf303faSJacob Faibussowitsch extern void PetscSFCheckGraphSet(PetscSF, int);
172abc8c78SJacob Faibussowitsch #else
1895fce210SBarry Smith   #if defined(PETSC_USE_DEBUG)
19a8f51744SPierre Jolivet     #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)
2095fce210SBarry Smith   #else
219371c9d4SSatish Balay     #define PetscSFCheckGraphSet(sf, arg) \
229371c9d4SSatish Balay       do { \
239371c9d4SSatish Balay       } while (0)
2495fce210SBarry Smith   #endif
252abc8c78SJacob Faibussowitsch #endif
2695fce210SBarry Smith 
274c8fdceaSLisandro Dalcin const char *const PetscSFDuplicateOptions[]     = {"CONFONLY", "RANKS", "GRAPH", "PetscSFDuplicateOption", "PETSCSF_DUPLICATE_", NULL};
281f40158dSVaclav Hapla const char *const PetscSFConcatenateRootModes[] = {"local", "shared", "global", "PetscSFConcatenateRootMode", "PETSCSF_CONCATENATE_ROOTMODE_", NULL};
2995fce210SBarry Smith 
308af6ec1cSBarry Smith /*@
3195fce210SBarry Smith   PetscSFCreate - create a star forest communication context
3295fce210SBarry Smith 
33d083f849SBarry Smith   Collective
3495fce210SBarry Smith 
354165533cSJose E. Roman   Input Parameter:
3695fce210SBarry Smith . comm - communicator on which the star forest will operate
3795fce210SBarry Smith 
384165533cSJose E. Roman   Output Parameter:
3995fce210SBarry Smith . sf - new star forest context
4095fce210SBarry Smith 
4120662ed9SBarry Smith   Options Database Key:
426677b1c1SJunchao Zhang + -sf_type basic                 - Use MPI persistent Isend/Irecv for communication (Default)
436677b1c1SJunchao Zhang . -sf_type window                - Use MPI-3 one-sided window for communication
446677b1c1SJunchao Zhang . -sf_type neighbor              - Use MPI-3 neighborhood collectives for communication
456677b1c1SJunchao Zhang - -sf_neighbor_persistent <bool> - If true, use MPI-4 persistent neighborhood collectives for communication (used along with -sf_type neighbor)
46dd5b3ca6SJunchao Zhang 
4795fce210SBarry Smith   Level: intermediate
4895fce210SBarry Smith 
49cab54364SBarry Smith   Note:
50cab54364SBarry Smith   When one knows the communication graph is one of the predefined graph, such as `MPI_Alltoall()`, `MPI_Allgatherv()`,
51cab54364SBarry Smith   `MPI_Gatherv()`, one can create a `PetscSF` and then set its graph with `PetscSFSetGraphWithPattern()`. These special
5220662ed9SBarry Smith   `SF`s are optimized and they have better performance than the general `SF`s.
53dd5b3ca6SJunchao Zhang 
5438b5cf2dSJacob Faibussowitsch .seealso: `PetscSF`, `PetscSFSetType`, `PetscSFSetGraph()`, `PetscSFSetGraphWithPattern()`, `PetscSFDestroy()`
5595fce210SBarry Smith @*/
56d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreate(MPI_Comm comm, PetscSF *sf)
57d71ae5a4SJacob Faibussowitsch {
5895fce210SBarry Smith   PetscSF b;
5995fce210SBarry Smith 
6095fce210SBarry Smith   PetscFunctionBegin;
614f572ea9SToby Isaac   PetscAssertPointer(sf, 2);
629566063dSJacob Faibussowitsch   PetscCall(PetscSFInitializePackage());
6395fce210SBarry Smith 
649566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(b, PETSCSF_CLASSID, "PetscSF", "Star Forest", "PetscSF", comm, PetscSFDestroy, PetscSFView));
6595fce210SBarry Smith 
6695fce210SBarry Smith   b->nroots    = -1;
6795fce210SBarry Smith   b->nleaves   = -1;
6829046d53SLisandro Dalcin   b->minleaf   = PETSC_MAX_INT;
6929046d53SLisandro Dalcin   b->maxleaf   = PETSC_MIN_INT;
7095fce210SBarry Smith   b->nranks    = -1;
7195fce210SBarry Smith   b->rankorder = PETSC_TRUE;
7295fce210SBarry Smith   b->ingroup   = MPI_GROUP_NULL;
7395fce210SBarry Smith   b->outgroup  = MPI_GROUP_NULL;
7495fce210SBarry Smith   b->graphset  = PETSC_FALSE;
7520c24465SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
7620c24465SJunchao Zhang   b->use_gpu_aware_mpi    = use_gpu_aware_mpi;
7720c24465SJunchao Zhang   b->use_stream_aware_mpi = PETSC_FALSE;
7871438e86SJunchao Zhang   b->unknown_input_stream = PETSC_FALSE;
7927f636e8SJunchao Zhang   #if defined(PETSC_HAVE_KOKKOS) /* Prefer kokkos over cuda*/
8020c24465SJunchao Zhang   b->backend = PETSCSF_BACKEND_KOKKOS;
8127f636e8SJunchao Zhang   #elif defined(PETSC_HAVE_CUDA)
8227f636e8SJunchao Zhang   b->backend = PETSCSF_BACKEND_CUDA;
8359af0bd3SScott Kruger   #elif defined(PETSC_HAVE_HIP)
8459af0bd3SScott Kruger   b->backend = PETSCSF_BACKEND_HIP;
8520c24465SJunchao Zhang   #endif
8671438e86SJunchao Zhang 
8771438e86SJunchao Zhang   #if defined(PETSC_HAVE_NVSHMEM)
8871438e86SJunchao Zhang   b->use_nvshmem     = PETSC_FALSE; /* Default is not to try NVSHMEM */
8971438e86SJunchao Zhang   b->use_nvshmem_get = PETSC_FALSE; /* Default is to use nvshmem_put based protocol */
909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_nvshmem", &b->use_nvshmem, NULL));
919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_nvshmem_get", &b->use_nvshmem_get, NULL));
9271438e86SJunchao Zhang   #endif
9320c24465SJunchao Zhang #endif
9460c22052SBarry Smith   b->vscat.from_n = -1;
9560c22052SBarry Smith   b->vscat.to_n   = -1;
9660c22052SBarry Smith   b->vscat.unit   = MPIU_SCALAR;
9795fce210SBarry Smith   *sf             = b;
983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9995fce210SBarry Smith }
10095fce210SBarry Smith 
10129046d53SLisandro Dalcin /*@
10295fce210SBarry Smith   PetscSFReset - Reset a star forest so that different sizes or neighbors can be used
10395fce210SBarry Smith 
10495fce210SBarry Smith   Collective
10595fce210SBarry Smith 
1064165533cSJose E. Roman   Input Parameter:
10795fce210SBarry Smith . sf - star forest
10895fce210SBarry Smith 
10995fce210SBarry Smith   Level: advanced
11095fce210SBarry Smith 
111cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`, `PetscSFSetGraph()`, `PetscSFDestroy()`
11295fce210SBarry Smith @*/
113d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFReset(PetscSF sf)
114d71ae5a4SJacob Faibussowitsch {
11595fce210SBarry Smith   PetscFunctionBegin;
11695fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
117dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf, Reset);
11829046d53SLisandro Dalcin   sf->nroots   = -1;
11929046d53SLisandro Dalcin   sf->nleaves  = -1;
12029046d53SLisandro Dalcin   sf->minleaf  = PETSC_MAX_INT;
12129046d53SLisandro Dalcin   sf->maxleaf  = PETSC_MIN_INT;
12295fce210SBarry Smith   sf->mine     = NULL;
12395fce210SBarry Smith   sf->remote   = NULL;
12429046d53SLisandro Dalcin   sf->graphset = PETSC_FALSE;
1259566063dSJacob Faibussowitsch   PetscCall(PetscFree(sf->mine_alloc));
1269566063dSJacob Faibussowitsch   PetscCall(PetscFree(sf->remote_alloc));
12721c688dcSJed Brown   sf->nranks = -1;
1289566063dSJacob Faibussowitsch   PetscCall(PetscFree4(sf->ranks, sf->roffset, sf->rmine, sf->rremote));
12929046d53SLisandro Dalcin   sf->degreeknown = PETSC_FALSE;
1309566063dSJacob Faibussowitsch   PetscCall(PetscFree(sf->degree));
1319566063dSJacob Faibussowitsch   if (sf->ingroup != MPI_GROUP_NULL) PetscCallMPI(MPI_Group_free(&sf->ingroup));
1329566063dSJacob Faibussowitsch   if (sf->outgroup != MPI_GROUP_NULL) PetscCallMPI(MPI_Group_free(&sf->outgroup));
133013b3241SStefano Zampini   if (sf->multi) sf->multi->multi = NULL;
1349566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf->multi));
1359566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&sf->map));
13671438e86SJunchao Zhang 
13771438e86SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
1389566063dSJacob Faibussowitsch   for (PetscInt i = 0; i < 2; i++) PetscCall(PetscSFFree(sf, PETSC_MEMTYPE_DEVICE, sf->rmine_d[i]));
13971438e86SJunchao Zhang #endif
14071438e86SJunchao Zhang 
14195fce210SBarry Smith   sf->setupcalled = PETSC_FALSE;
1423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14395fce210SBarry Smith }
14495fce210SBarry Smith 
14595fce210SBarry Smith /*@C
146cab54364SBarry Smith   PetscSFSetType - Set the `PetscSF` communication implementation
14795fce210SBarry Smith 
148c3339decSBarry Smith   Collective
14995fce210SBarry Smith 
15095fce210SBarry Smith   Input Parameters:
151cab54364SBarry Smith + sf   - the `PetscSF` context
15295fce210SBarry Smith - type - a known method
153cab54364SBarry Smith .vb
154cab54364SBarry Smith     PETSCSFWINDOW - MPI-2/3 one-sided
155cab54364SBarry Smith     PETSCSFBASIC - basic implementation using MPI-1 two-sided
156cab54364SBarry Smith .ve
15795fce210SBarry Smith 
15895fce210SBarry Smith   Options Database Key:
15920662ed9SBarry Smith . -sf_type <type> - Sets the method; for example `basic` or `window` use -help for a list of available methods
160cab54364SBarry Smith 
161cab54364SBarry Smith   Level: intermediate
16295fce210SBarry Smith 
16395fce210SBarry Smith   Notes:
16420662ed9SBarry Smith   See `PetscSFType` for possible values
16595fce210SBarry Smith 
16620662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFCreate()`
16795fce210SBarry Smith @*/
168d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetType(PetscSF sf, PetscSFType type)
169d71ae5a4SJacob Faibussowitsch {
17095fce210SBarry Smith   PetscBool match;
1715f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(PetscSF);
17295fce210SBarry Smith 
17395fce210SBarry Smith   PetscFunctionBegin;
17495fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
1754f572ea9SToby Isaac   PetscAssertPointer(type, 2);
17695fce210SBarry Smith 
1779566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)sf, type, &match));
1783ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
17995fce210SBarry Smith 
1809566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(PetscSFList, type, &r));
1816adde796SStefano Zampini   PetscCheck(r, PetscObjectComm((PetscObject)sf), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested PetscSF type %s", type);
18229046d53SLisandro Dalcin   /* Destroy the previous PetscSF implementation context */
183dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf, Destroy);
1849566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(sf->ops, sizeof(*sf->ops)));
1859566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)sf, type));
1869566063dSJacob Faibussowitsch   PetscCall((*r)(sf));
1873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18895fce210SBarry Smith }
18995fce210SBarry Smith 
19029046d53SLisandro Dalcin /*@C
191cab54364SBarry Smith   PetscSFGetType - Get the `PetscSF` communication implementation
19229046d53SLisandro Dalcin 
19329046d53SLisandro Dalcin   Not Collective
19429046d53SLisandro Dalcin 
19529046d53SLisandro Dalcin   Input Parameter:
196cab54364SBarry Smith . sf - the `PetscSF` context
19729046d53SLisandro Dalcin 
19829046d53SLisandro Dalcin   Output Parameter:
199cab54364SBarry Smith . type - the `PetscSF` type name
20029046d53SLisandro Dalcin 
20129046d53SLisandro Dalcin   Level: intermediate
20229046d53SLisandro Dalcin 
20320662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFSetType()`, `PetscSFCreate()`
20429046d53SLisandro Dalcin @*/
205d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type)
206d71ae5a4SJacob Faibussowitsch {
20729046d53SLisandro Dalcin   PetscFunctionBegin;
20829046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
2094f572ea9SToby Isaac   PetscAssertPointer(type, 2);
21029046d53SLisandro Dalcin   *type = ((PetscObject)sf)->type_name;
2113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21229046d53SLisandro Dalcin }
21329046d53SLisandro Dalcin 
2141fb7b255SJunchao Zhang /*@C
21520662ed9SBarry Smith   PetscSFDestroy - destroy a star forest
21695fce210SBarry Smith 
21795fce210SBarry Smith   Collective
21895fce210SBarry Smith 
2194165533cSJose E. Roman   Input Parameter:
22095fce210SBarry Smith . sf - address of star forest
22195fce210SBarry Smith 
22295fce210SBarry Smith   Level: intermediate
22395fce210SBarry Smith 
22420662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFCreate()`, `PetscSFReset()`
22595fce210SBarry Smith @*/
226d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFDestroy(PetscSF *sf)
227d71ae5a4SJacob Faibussowitsch {
22895fce210SBarry Smith   PetscFunctionBegin;
2293ba16761SJacob Faibussowitsch   if (!*sf) PetscFunctionReturn(PETSC_SUCCESS);
230*f4f49eeaSPierre Jolivet   PetscValidHeaderSpecific(*sf, PETSCSF_CLASSID, 1);
231*f4f49eeaSPierre Jolivet   if (--((PetscObject)*sf)->refct > 0) {
2329371c9d4SSatish Balay     *sf = NULL;
2333ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2349371c9d4SSatish Balay   }
2359566063dSJacob Faibussowitsch   PetscCall(PetscSFReset(*sf));
236*f4f49eeaSPierre Jolivet   PetscTryTypeMethod(*sf, Destroy);
2379566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&(*sf)->vscat.lsf));
2389566063dSJacob Faibussowitsch   if ((*sf)->vscat.bs > 1) PetscCallMPI(MPI_Type_free(&(*sf)->vscat.unit));
239c02794c0SJunchao Zhang #if defined(PETSC_HAVE_CUDA) && defined(PETSC_HAVE_MPIX_STREAM)
240715b587bSJunchao Zhang   if ((*sf)->use_stream_aware_mpi) {
241715b587bSJunchao Zhang     PetscCallMPI(MPIX_Stream_free(&(*sf)->mpi_stream));
242715b587bSJunchao Zhang     PetscCallMPI(MPI_Comm_free(&(*sf)->stream_comm));
243715b587bSJunchao Zhang   }
244715b587bSJunchao Zhang #endif
2459566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(sf));
2463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24795fce210SBarry Smith }
24895fce210SBarry Smith 
249d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFCheckGraphValid_Private(PetscSF sf)
250d71ae5a4SJacob Faibussowitsch {
251c4e6a40aSLawrence Mitchell   PetscInt           i, nleaves;
252c4e6a40aSLawrence Mitchell   PetscMPIInt        size;
253c4e6a40aSLawrence Mitchell   const PetscInt    *ilocal;
254c4e6a40aSLawrence Mitchell   const PetscSFNode *iremote;
255c4e6a40aSLawrence Mitchell 
256c4e6a40aSLawrence Mitchell   PetscFunctionBegin;
2573ba16761SJacob Faibussowitsch   if (!sf->graphset || !PetscDefined(USE_DEBUG)) PetscFunctionReturn(PETSC_SUCCESS);
2589566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, &ilocal, &iremote));
2599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf), &size));
260c4e6a40aSLawrence Mitchell   for (i = 0; i < nleaves; i++) {
261c4e6a40aSLawrence Mitchell     const PetscInt rank   = iremote[i].rank;
262c4e6a40aSLawrence Mitchell     const PetscInt remote = iremote[i].index;
263c4e6a40aSLawrence Mitchell     const PetscInt leaf   = ilocal ? ilocal[i] : i;
264c9cc58a2SBarry 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);
26508401ef6SPierre 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);
26608401ef6SPierre 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);
267c4e6a40aSLawrence Mitchell   }
2683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
269c4e6a40aSLawrence Mitchell }
270c4e6a40aSLawrence Mitchell 
27195fce210SBarry Smith /*@
27220662ed9SBarry Smith   PetscSFSetUp - set up communication structures for a `PetscSF`, after this is done it may be used to perform communication
27395fce210SBarry Smith 
27495fce210SBarry Smith   Collective
27595fce210SBarry Smith 
2764165533cSJose E. Roman   Input Parameter:
27795fce210SBarry Smith . sf - star forest communication object
27895fce210SBarry Smith 
27995fce210SBarry Smith   Level: beginner
28095fce210SBarry Smith 
28120662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFSetFromOptions()`, `PetscSFSetType()`
28295fce210SBarry Smith @*/
283d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetUp(PetscSF sf)
284d71ae5a4SJacob Faibussowitsch {
28595fce210SBarry Smith   PetscFunctionBegin;
28629046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
28729046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf, 1);
2883ba16761SJacob Faibussowitsch   if (sf->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
2899566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_SetUp, sf, 0, 0, 0));
2909566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckGraphValid_Private(sf));
2919566063dSJacob Faibussowitsch   if (!((PetscObject)sf)->type_name) PetscCall(PetscSFSetType(sf, PETSCSFBASIC)); /* Zero all sf->ops */
292dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf, SetUp);
29320c24465SJunchao Zhang #if defined(PETSC_HAVE_CUDA)
29420c24465SJunchao Zhang   if (sf->backend == PETSCSF_BACKEND_CUDA) {
29571438e86SJunchao Zhang     sf->ops->Malloc = PetscSFMalloc_CUDA;
29671438e86SJunchao Zhang     sf->ops->Free   = PetscSFFree_CUDA;
29720c24465SJunchao Zhang   }
29820c24465SJunchao Zhang #endif
29959af0bd3SScott Kruger #if defined(PETSC_HAVE_HIP)
30059af0bd3SScott Kruger   if (sf->backend == PETSCSF_BACKEND_HIP) {
30159af0bd3SScott Kruger     sf->ops->Malloc = PetscSFMalloc_HIP;
30259af0bd3SScott Kruger     sf->ops->Free   = PetscSFFree_HIP;
30359af0bd3SScott Kruger   }
30459af0bd3SScott Kruger #endif
30520c24465SJunchao Zhang 
30620c24465SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
30720c24465SJunchao Zhang   if (sf->backend == PETSCSF_BACKEND_KOKKOS) {
30820c24465SJunchao Zhang     sf->ops->Malloc = PetscSFMalloc_Kokkos;
30920c24465SJunchao Zhang     sf->ops->Free   = PetscSFFree_Kokkos;
31020c24465SJunchao Zhang   }
31120c24465SJunchao Zhang #endif
3129566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_SetUp, sf, 0, 0, 0));
31395fce210SBarry Smith   sf->setupcalled = PETSC_TRUE;
3143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31595fce210SBarry Smith }
31695fce210SBarry Smith 
3178af6ec1cSBarry Smith /*@
318cab54364SBarry Smith   PetscSFSetFromOptions - set `PetscSF` options using the options database
31995fce210SBarry Smith 
32095fce210SBarry Smith   Logically Collective
32195fce210SBarry Smith 
3224165533cSJose E. Roman   Input Parameter:
32395fce210SBarry Smith . sf - star forest
32495fce210SBarry Smith 
32595fce210SBarry Smith   Options Database Keys:
32620662ed9SBarry Smith + -sf_type                                                                                                         - implementation type, see `PetscSFSetType()`
32751ccb202SJunchao Zhang . -sf_rank_order                                                                                                   - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
32820662ed9SBarry Smith . -sf_use_default_stream                                                                                           - Assume callers of `PetscSF` computed the input root/leafdata with the default CUDA stream. `PetscSF` will also
32920662ed9SBarry Smith                             use the default stream to process data. Therefore, no stream synchronization is needed between `PetscSF` and its caller (default: true).
33020662ed9SBarry Smith                             If true, this option only works with `-use_gpu_aware_mpi 1`.
33120662ed9SBarry Smith . -sf_use_stream_aware_mpi                                                                                         - Assume the underlying MPI is CUDA-stream aware and `PetscSF` won't sync streams for send/recv buffers passed to MPI (default: false).
33220662ed9SBarry Smith                                If true, this option only works with `-use_gpu_aware_mpi 1`.
33395fce210SBarry Smith 
33438b5cf2dSJacob Faibussowitsch - -sf_backend cuda | hip | kokkos -Select the device backend SF uses. Currently `PetscSF` has these backends: cuda - hip and Kokkos.
33559af0bd3SScott Kruger                               On CUDA (HIP) devices, one can choose cuda (hip) or kokkos with the default being kokkos. On other devices,
33620c24465SJunchao Zhang                               the only available is kokkos.
33720c24465SJunchao Zhang 
33895fce210SBarry Smith   Level: intermediate
339cab54364SBarry Smith 
340cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`, `PetscSFSetType()`
34195fce210SBarry Smith @*/
342d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
343d71ae5a4SJacob Faibussowitsch {
34495fce210SBarry Smith   PetscSFType deft;
34595fce210SBarry Smith   char        type[256];
34695fce210SBarry Smith   PetscBool   flg;
34795fce210SBarry Smith 
34895fce210SBarry Smith   PetscFunctionBegin;
34995fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
350d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)sf);
35195fce210SBarry Smith   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
3529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-sf_type", "PetscSF implementation type", "PetscSFSetType", PetscSFList, deft, type, sizeof(type), &flg));
3539566063dSJacob Faibussowitsch   PetscCall(PetscSFSetType(sf, flg ? type : deft));
3549566063dSJacob 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));
3557fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
35620c24465SJunchao Zhang   {
35720c24465SJunchao Zhang     char      backendstr[32] = {0};
35859af0bd3SScott Kruger     PetscBool isCuda = PETSC_FALSE, isHip = PETSC_FALSE, isKokkos = PETSC_FALSE, set;
35920c24465SJunchao Zhang     /* Change the defaults set in PetscSFCreate() with command line options */
360d5b43468SJose E. Roman     PetscCall(PetscOptionsBool("-sf_unknown_input_stream", "SF root/leafdata is computed on arbitrary streams unknown to SF", "PetscSFSetFromOptions", sf->unknown_input_stream, &sf->unknown_input_stream, NULL));
3619566063dSJacob 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));
3629566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-sf_backend", "Select the device backend SF uses", "PetscSFSetFromOptions", NULL, backendstr, sizeof(backendstr), &set));
3639566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp("cuda", backendstr, &isCuda));
3649566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp("kokkos", backendstr, &isKokkos));
3659566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp("hip", backendstr, &isHip));
36659af0bd3SScott Kruger   #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
36720c24465SJunchao Zhang     if (isCuda) sf->backend = PETSCSF_BACKEND_CUDA;
36820c24465SJunchao Zhang     else if (isKokkos) sf->backend = PETSCSF_BACKEND_KOKKOS;
36959af0bd3SScott Kruger     else if (isHip) sf->backend = PETSCSF_BACKEND_HIP;
37028b400f6SJacob 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);
37120c24465SJunchao Zhang   #elif defined(PETSC_HAVE_KOKKOS)
37208401ef6SPierre Jolivet     PetscCheck(!set || isKokkos, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "-sf_backend %s is not supported. You can only choose kokkos", backendstr);
37320c24465SJunchao Zhang   #endif
374715b587bSJunchao Zhang 
375715b587bSJunchao Zhang   #if defined(PETSC_HAVE_CUDA) && defined(PETSC_HAVE_MPIX_STREAM)
376715b587bSJunchao Zhang     if (sf->use_stream_aware_mpi) {
377715b587bSJunchao Zhang       MPI_Info info;
378715b587bSJunchao Zhang 
379715b587bSJunchao Zhang       PetscCallMPI(MPI_Info_create(&info));
380715b587bSJunchao Zhang       PetscCallMPI(MPI_Info_set(info, "type", "cudaStream_t"));
381715b587bSJunchao Zhang       PetscCallMPI(MPIX_Info_set_hex(info, "value", &PetscDefaultCudaStream, sizeof(PetscDefaultCudaStream)));
382715b587bSJunchao Zhang       PetscCallMPI(MPIX_Stream_create(info, &sf->mpi_stream));
383715b587bSJunchao Zhang       PetscCallMPI(MPI_Info_free(&info));
384715b587bSJunchao Zhang       PetscCallMPI(MPIX_Stream_comm_create(PetscObjectComm((PetscObject)sf), sf->mpi_stream, &sf->stream_comm));
385715b587bSJunchao Zhang     }
386715b587bSJunchao Zhang   #endif
38720c24465SJunchao Zhang   }
388c2a741eeSJunchao Zhang #endif
389dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf, SetFromOptions, PetscOptionsObject);
390d0609cedSBarry Smith   PetscOptionsEnd();
3913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39295fce210SBarry Smith }
39395fce210SBarry Smith 
39429046d53SLisandro Dalcin /*@
39595fce210SBarry Smith   PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order
39695fce210SBarry Smith 
39795fce210SBarry Smith   Logically Collective
39895fce210SBarry Smith 
3994165533cSJose E. Roman   Input Parameters:
40095fce210SBarry Smith + sf  - star forest
401cab54364SBarry Smith - flg - `PETSC_TRUE` to sort, `PETSC_FALSE` to skip sorting (lower setup cost, but non-deterministic)
40295fce210SBarry Smith 
40395fce210SBarry Smith   Level: advanced
40495fce210SBarry Smith 
40520662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFGatherBegin()`, `PetscSFScatterBegin()`
40695fce210SBarry Smith @*/
407d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetRankOrder(PetscSF sf, PetscBool flg)
408d71ae5a4SJacob Faibussowitsch {
40995fce210SBarry Smith   PetscFunctionBegin;
41095fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
41195fce210SBarry Smith   PetscValidLogicalCollectiveBool(sf, flg, 2);
41228b400f6SJacob Faibussowitsch   PetscCheck(!sf->multi, PetscObjectComm((PetscObject)sf), PETSC_ERR_ARG_WRONGSTATE, "Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
41395fce210SBarry Smith   sf->rankorder = flg;
4143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41595fce210SBarry Smith }
41695fce210SBarry Smith 
4178dbb0df6SBarry Smith /*@C
41895fce210SBarry Smith   PetscSFSetGraph - Set a parallel star forest
41995fce210SBarry Smith 
42095fce210SBarry Smith   Collective
42195fce210SBarry Smith 
4224165533cSJose E. Roman   Input Parameters:
42395fce210SBarry Smith + sf         - star forest
42495fce210SBarry Smith . nroots     - number of root vertices on the current process (these are possible targets for other process to attach leaves)
42595fce210SBarry Smith . nleaves    - number of leaf vertices on the current process, each of these references a root on any process
42620662ed9SBarry Smith . ilocal     - locations of leaves in leafdata buffers, pass `NULL` for contiguous storage (locations must be >= 0, enforced
427c4e6a40aSLawrence Mitchell during setup in debug mode)
42820662ed9SBarry Smith . localmode  - copy mode for `ilocal`
429c4e6a40aSLawrence Mitchell . iremote    - remote locations of root vertices for each leaf on the current process (locations must be >= 0, enforced
430c4e6a40aSLawrence Mitchell during setup in debug mode)
43120662ed9SBarry Smith - remotemode - copy mode for `iremote`
43295fce210SBarry Smith 
43395fce210SBarry Smith   Level: intermediate
43495fce210SBarry Smith 
43595452b02SPatrick Sanan   Notes:
43620662ed9SBarry Smith   Leaf indices in `ilocal` must be unique, otherwise an error occurs.
43738ab3f8aSBarry Smith 
43820662ed9SBarry Smith   Input arrays `ilocal` and `iremote` follow the `PetscCopyMode` semantics.
43920662ed9SBarry Smith   In particular, if `localmode` or `remotemode` is `PETSC_OWN_POINTER` or `PETSC_USE_POINTER`,
440db2b9530SVaclav Hapla   PETSc might modify the respective array;
44120662ed9SBarry Smith   if `PETSC_USE_POINTER`, the user must delete the array after `PetscSFDestroy()`.
442cab54364SBarry Smith   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).
443db2b9530SVaclav Hapla 
44438b5cf2dSJacob Faibussowitsch   Fortran Notes:
44520662ed9SBarry Smith   In Fortran you must use `PETSC_COPY_VALUES` for `localmode` and `remotemode`.
446c4e6a40aSLawrence Mitchell 
44738b5cf2dSJacob Faibussowitsch   Developer Notes:
448db2b9530SVaclav Hapla   We sort leaves to check for duplicates and contiguousness and to find minleaf/maxleaf.
44920662ed9SBarry Smith   This also allows to compare leaf sets of two `PetscSF`s easily.
45072bf8598SVaclav Hapla 
45120662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFGetGraph()`
45295fce210SBarry Smith @*/
453d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraph(PetscSF sf, PetscInt nroots, PetscInt nleaves, PetscInt *ilocal, PetscCopyMode localmode, PetscSFNode *iremote, PetscCopyMode remotemode)
454d71ae5a4SJacob Faibussowitsch {
455db2b9530SVaclav Hapla   PetscBool unique, contiguous;
45695fce210SBarry Smith 
45795fce210SBarry Smith   PetscFunctionBegin;
45895fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
4594f572ea9SToby Isaac   if (nleaves > 0 && ilocal) PetscAssertPointer(ilocal, 4);
4604f572ea9SToby Isaac   if (nleaves > 0) PetscAssertPointer(iremote, 6);
46108401ef6SPierre Jolivet   PetscCheck(nroots >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nroots %" PetscInt_FMT ", cannot be negative", nroots);
46208401ef6SPierre Jolivet   PetscCheck(nleaves >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nleaves %" PetscInt_FMT ", cannot be negative", nleaves);
4638da24d32SBarry Smith   /* enums may be handled as unsigned by some compilers, NVHPC for example, the int cast
4648da24d32SBarry Smith    * below is to prevent NVHPC from warning about meaningless comparison of unsigned with zero */
4658da24d32SBarry Smith   PetscCheck((int)localmode >= PETSC_COPY_VALUES && localmode <= PETSC_USE_POINTER, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wrong localmode %d", localmode);
4668da24d32SBarry Smith   PetscCheck((int)remotemode >= PETSC_COPY_VALUES && remotemode <= PETSC_USE_POINTER, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wrong remotemode %d", remotemode);
46729046d53SLisandro Dalcin 
4682a67d2daSStefano Zampini   if (sf->nroots >= 0) { /* Reset only if graph already set */
4699566063dSJacob Faibussowitsch     PetscCall(PetscSFReset(sf));
4702a67d2daSStefano Zampini   }
4712a67d2daSStefano Zampini 
4729566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_SetGraph, sf, 0, 0, 0));
47329046d53SLisandro Dalcin 
47495fce210SBarry Smith   sf->nroots  = nroots;
47595fce210SBarry Smith   sf->nleaves = nleaves;
47629046d53SLisandro Dalcin 
477db2b9530SVaclav Hapla   if (localmode == PETSC_COPY_VALUES && ilocal) {
478db2b9530SVaclav Hapla     PetscInt *tlocal = NULL;
479db2b9530SVaclav Hapla 
4809566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &tlocal));
4819566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tlocal, ilocal, nleaves));
482db2b9530SVaclav Hapla     ilocal = tlocal;
483db2b9530SVaclav Hapla   }
484db2b9530SVaclav Hapla   if (remotemode == PETSC_COPY_VALUES) {
485db2b9530SVaclav Hapla     PetscSFNode *tremote = NULL;
486db2b9530SVaclav Hapla 
4879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &tremote));
4889566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tremote, iremote, nleaves));
489db2b9530SVaclav Hapla     iremote = tremote;
490db2b9530SVaclav Hapla   }
491db2b9530SVaclav Hapla 
49229046d53SLisandro Dalcin   if (nleaves && ilocal) {
493db2b9530SVaclav Hapla     PetscSFNode work;
494db2b9530SVaclav Hapla 
4959566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithDataArray(nleaves, ilocal, iremote, sizeof(PetscSFNode), &work));
4969566063dSJacob Faibussowitsch     PetscCall(PetscSortedCheckDupsInt(nleaves, ilocal, &unique));
497db2b9530SVaclav Hapla     unique = PetscNot(unique);
498db2b9530SVaclav 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");
499db2b9530SVaclav Hapla     sf->minleaf = ilocal[0];
500db2b9530SVaclav Hapla     sf->maxleaf = ilocal[nleaves - 1];
501db2b9530SVaclav Hapla     contiguous  = (PetscBool)(unique && ilocal[0] == 0 && ilocal[nleaves - 1] == nleaves - 1);
50229046d53SLisandro Dalcin   } else {
50329046d53SLisandro Dalcin     sf->minleaf = 0;
50429046d53SLisandro Dalcin     sf->maxleaf = nleaves - 1;
505db2b9530SVaclav Hapla     unique      = PETSC_TRUE;
506db2b9530SVaclav Hapla     contiguous  = PETSC_TRUE;
50729046d53SLisandro Dalcin   }
50829046d53SLisandro Dalcin 
509db2b9530SVaclav Hapla   if (contiguous) {
510db2b9530SVaclav Hapla     if (localmode == PETSC_USE_POINTER) {
511db2b9530SVaclav Hapla       ilocal = NULL;
512db2b9530SVaclav Hapla     } else {
5139566063dSJacob Faibussowitsch       PetscCall(PetscFree(ilocal));
514db2b9530SVaclav Hapla     }
515db2b9530SVaclav Hapla   }
516db2b9530SVaclav Hapla   sf->mine = ilocal;
517db2b9530SVaclav Hapla   if (localmode == PETSC_USE_POINTER) {
51829046d53SLisandro Dalcin     sf->mine_alloc = NULL;
519db2b9530SVaclav Hapla   } else {
520db2b9530SVaclav Hapla     sf->mine_alloc = ilocal;
52195fce210SBarry Smith   }
522db2b9530SVaclav Hapla   sf->remote = iremote;
523db2b9530SVaclav Hapla   if (remotemode == PETSC_USE_POINTER) {
52429046d53SLisandro Dalcin     sf->remote_alloc = NULL;
525db2b9530SVaclav Hapla   } else {
526db2b9530SVaclav Hapla     sf->remote_alloc = iremote;
52795fce210SBarry Smith   }
5289566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_SetGraph, sf, 0, 0, 0));
52929046d53SLisandro Dalcin   sf->graphset = PETSC_TRUE;
5303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53195fce210SBarry Smith }
53295fce210SBarry Smith 
53329046d53SLisandro Dalcin /*@
534cab54364SBarry Smith   PetscSFSetGraphWithPattern - Sets the graph of a `PetscSF` with a specific pattern
535dd5b3ca6SJunchao Zhang 
536dd5b3ca6SJunchao Zhang   Collective
537dd5b3ca6SJunchao Zhang 
538dd5b3ca6SJunchao Zhang   Input Parameters:
539cab54364SBarry Smith + sf      - The `PetscSF`
540cab54364SBarry Smith . map     - Layout of roots over all processes (insignificant when pattern is `PETSCSF_PATTERN_ALLTOALL`)
541cab54364SBarry Smith - pattern - One of `PETSCSF_PATTERN_ALLGATHER`, `PETSCSF_PATTERN_GATHER`, `PETSCSF_PATTERN_ALLTOALL`
542cab54364SBarry Smith 
543cab54364SBarry Smith   Level: intermediate
544dd5b3ca6SJunchao Zhang 
545dd5b3ca6SJunchao Zhang   Notes:
54620662ed9SBarry Smith   It is easier to explain `PetscSFPattern` using vectors. Suppose we have an MPI vector `x` and its `PetscLayout` is `map`.
54720662ed9SBarry Smith   `n` and `N` are the local and global sizes of `x` respectively.
548dd5b3ca6SJunchao Zhang 
54920662ed9SBarry Smith   With `PETSCSF_PATTERN_ALLGATHER`, the routine creates a graph that if one does `PetscSFBcastBegin()`/`PetscSFBcastEnd()` on it, it will copy `x` to
55020662ed9SBarry Smith   sequential vectors `y` on all MPI processes.
551dd5b3ca6SJunchao Zhang 
55220662ed9SBarry Smith   With `PETSCSF_PATTERN_GATHER`, the routine creates a graph that if one does `PetscSFBcastBegin()`/`PetscSFBcastEnd()` on it, it will copy `x` to a
55320662ed9SBarry Smith   sequential vector `y` on rank 0.
554dd5b3ca6SJunchao Zhang 
55520662ed9SBarry Smith   In above cases, entries of `x` are roots and entries of `y` are leaves.
556dd5b3ca6SJunchao Zhang 
55720662ed9SBarry Smith   With `PETSCSF_PATTERN_ALLTOALL`, map is insignificant. Suppose NP is size of `sf`'s communicator. The routine
558dd5b3ca6SJunchao Zhang   creates a graph that every rank has NP leaves and NP roots. On rank i, its leaf j is connected to root i
559cab54364SBarry Smith   of rank j. Here 0 <=i,j<NP. It is a kind of `MPI_Alltoall()` with sendcount/recvcount being 1. Note that it does
560dd5b3ca6SJunchao Zhang   not mean one can not send multiple items. One just needs to create a new MPI datatype for the mulptiple data
561cab54364SBarry Smith   items with `MPI_Type_contiguous` and use that as the <unit> argument in SF routines.
562dd5b3ca6SJunchao Zhang 
563dd5b3ca6SJunchao Zhang   In this case, roots and leaves are symmetric.
564dd5b3ca6SJunchao Zhang 
565cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFGetGraph()`
566dd5b3ca6SJunchao Zhang  @*/
567d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphWithPattern(PetscSF sf, PetscLayout map, PetscSFPattern pattern)
568d71ae5a4SJacob Faibussowitsch {
569dd5b3ca6SJunchao Zhang   MPI_Comm    comm;
570dd5b3ca6SJunchao Zhang   PetscInt    n, N, res[2];
571dd5b3ca6SJunchao Zhang   PetscMPIInt rank, size;
572dd5b3ca6SJunchao Zhang   PetscSFType type;
573dd5b3ca6SJunchao Zhang 
574dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
5752abc8c78SJacob Faibussowitsch   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
5764f572ea9SToby Isaac   if (pattern != PETSCSF_PATTERN_ALLTOALL) PetscAssertPointer(map, 2);
5779566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
5782c71b3e2SJacob Faibussowitsch   PetscCheck(pattern >= PETSCSF_PATTERN_ALLGATHER && pattern <= PETSCSF_PATTERN_ALLTOALL, comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PetscSFPattern %d", pattern);
5799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5809566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
581dd5b3ca6SJunchao Zhang 
582dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLTOALL) {
583dd5b3ca6SJunchao Zhang     type = PETSCSFALLTOALL;
5849566063dSJacob Faibussowitsch     PetscCall(PetscLayoutCreate(comm, &sf->map));
5859566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetLocalSize(sf->map, size));
5869566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetSize(sf->map, ((PetscInt)size) * size));
5879566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetUp(sf->map));
588dd5b3ca6SJunchao Zhang   } else {
5899566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetLocalSize(map, &n));
5909566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetSize(map, &N));
591dd5b3ca6SJunchao Zhang     res[0] = n;
592dd5b3ca6SJunchao Zhang     res[1] = -n;
593dd5b3ca6SJunchao Zhang     /* Check if n are same over all ranks so that we can optimize it */
5941c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, res, 2, MPIU_INT, MPI_MAX, comm));
595dd5b3ca6SJunchao Zhang     if (res[0] == -res[1]) { /* same n */
596dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHER : PETSCSFGATHER;
597dd5b3ca6SJunchao Zhang     } else {
598dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHERV : PETSCSFGATHERV;
599dd5b3ca6SJunchao Zhang     }
6009566063dSJacob Faibussowitsch     PetscCall(PetscLayoutReference(map, &sf->map));
601dd5b3ca6SJunchao Zhang   }
6029566063dSJacob Faibussowitsch   PetscCall(PetscSFSetType(sf, type));
603dd5b3ca6SJunchao Zhang 
604dd5b3ca6SJunchao Zhang   sf->pattern = pattern;
605dd5b3ca6SJunchao Zhang   sf->mine    = NULL; /* Contiguous */
606dd5b3ca6SJunchao Zhang 
607dd5b3ca6SJunchao Zhang   /* Set nleaves, nroots here in case user calls PetscSFGetGraph, which is legal to call even before PetscSFSetUp is called.
608dd5b3ca6SJunchao Zhang      Also set other easy stuff.
609dd5b3ca6SJunchao Zhang    */
610dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLGATHER) {
611dd5b3ca6SJunchao Zhang     sf->nleaves = N;
612dd5b3ca6SJunchao Zhang     sf->nroots  = n;
613dd5b3ca6SJunchao Zhang     sf->nranks  = size;
614dd5b3ca6SJunchao Zhang     sf->minleaf = 0;
615dd5b3ca6SJunchao Zhang     sf->maxleaf = N - 1;
616dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_GATHER) {
617dd5b3ca6SJunchao Zhang     sf->nleaves = rank ? 0 : N;
618dd5b3ca6SJunchao Zhang     sf->nroots  = n;
619dd5b3ca6SJunchao Zhang     sf->nranks  = rank ? 0 : size;
620dd5b3ca6SJunchao Zhang     sf->minleaf = 0;
621dd5b3ca6SJunchao Zhang     sf->maxleaf = rank ? -1 : N - 1;
622dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_ALLTOALL) {
623dd5b3ca6SJunchao Zhang     sf->nleaves = size;
624dd5b3ca6SJunchao Zhang     sf->nroots  = size;
625dd5b3ca6SJunchao Zhang     sf->nranks  = size;
626dd5b3ca6SJunchao Zhang     sf->minleaf = 0;
627dd5b3ca6SJunchao Zhang     sf->maxleaf = size - 1;
628dd5b3ca6SJunchao Zhang   }
629dd5b3ca6SJunchao Zhang   sf->ndranks  = 0; /* We do not need to separate out distinguished ranks for patterned graphs to improve communication performance */
630dd5b3ca6SJunchao Zhang   sf->graphset = PETSC_TRUE;
6313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
632dd5b3ca6SJunchao Zhang }
633dd5b3ca6SJunchao Zhang 
634dd5b3ca6SJunchao Zhang /*@
635cab54364SBarry Smith   PetscSFCreateInverseSF - given a `PetscSF` in which all vertices have degree 1, creates the inverse map
63695fce210SBarry Smith 
63795fce210SBarry Smith   Collective
63895fce210SBarry Smith 
6394165533cSJose E. Roman   Input Parameter:
64095fce210SBarry Smith . sf - star forest to invert
64195fce210SBarry Smith 
6424165533cSJose E. Roman   Output Parameter:
64320662ed9SBarry Smith . isf - inverse of `sf`
6444165533cSJose E. Roman 
64595fce210SBarry Smith   Level: advanced
64695fce210SBarry Smith 
64795fce210SBarry Smith   Notes:
64895fce210SBarry Smith   All roots must have degree 1.
64995fce210SBarry Smith 
65095fce210SBarry Smith   The local space may be a permutation, but cannot be sparse.
65195fce210SBarry Smith 
65220662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFSetGraph()`
65395fce210SBarry Smith @*/
654d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateInverseSF(PetscSF sf, PetscSF *isf)
655d71ae5a4SJacob Faibussowitsch {
65695fce210SBarry Smith   PetscMPIInt     rank;
65795fce210SBarry Smith   PetscInt        i, nroots, nleaves, maxlocal, count, *newilocal;
65895fce210SBarry Smith   const PetscInt *ilocal;
65995fce210SBarry Smith   PetscSFNode    *roots, *leaves;
66095fce210SBarry Smith 
66195fce210SBarry Smith   PetscFunctionBegin;
66229046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
66329046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf, 1);
6644f572ea9SToby Isaac   PetscAssertPointer(isf, 2);
66529046d53SLisandro Dalcin 
6669566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL));
66729046d53SLisandro Dalcin   maxlocal = sf->maxleaf + 1; /* TODO: We should use PetscSFGetLeafRange() */
66829046d53SLisandro Dalcin 
6699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
6709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nroots, &roots, maxlocal, &leaves));
671ae9aee6dSMatthew G. Knepley   for (i = 0; i < maxlocal; i++) {
67295fce210SBarry Smith     leaves[i].rank  = rank;
67395fce210SBarry Smith     leaves[i].index = i;
67495fce210SBarry Smith   }
67595fce210SBarry Smith   for (i = 0; i < nroots; i++) {
67695fce210SBarry Smith     roots[i].rank  = -1;
67795fce210SBarry Smith     roots[i].index = -1;
67895fce210SBarry Smith   }
6799566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf, MPIU_2INT, leaves, roots, MPI_REPLACE));
6809566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sf, MPIU_2INT, leaves, roots, MPI_REPLACE));
68195fce210SBarry Smith 
68295fce210SBarry Smith   /* Check whether our leaves are sparse */
6839371c9d4SSatish Balay   for (i = 0, count = 0; i < nroots; i++)
6849371c9d4SSatish Balay     if (roots[i].rank >= 0) count++;
68595fce210SBarry Smith   if (count == nroots) newilocal = NULL;
6869371c9d4SSatish Balay   else { /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */ PetscCall(PetscMalloc1(count, &newilocal));
68795fce210SBarry Smith     for (i = 0, count = 0; i < nroots; i++) {
68895fce210SBarry Smith       if (roots[i].rank >= 0) {
68995fce210SBarry Smith         newilocal[count]   = i;
69095fce210SBarry Smith         roots[count].rank  = roots[i].rank;
69195fce210SBarry Smith         roots[count].index = roots[i].index;
69295fce210SBarry Smith         count++;
69395fce210SBarry Smith       }
69495fce210SBarry Smith     }
69595fce210SBarry Smith   }
69695fce210SBarry Smith 
6979566063dSJacob Faibussowitsch   PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_CONFONLY, isf));
6989566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*isf, maxlocal, count, newilocal, PETSC_OWN_POINTER, roots, PETSC_COPY_VALUES));
6999566063dSJacob Faibussowitsch   PetscCall(PetscFree2(roots, leaves));
7003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70195fce210SBarry Smith }
70295fce210SBarry Smith 
70395fce210SBarry Smith /*@
704cab54364SBarry Smith   PetscSFDuplicate - duplicate a `PetscSF`, optionally preserving rank connectivity and graph
70595fce210SBarry Smith 
70695fce210SBarry Smith   Collective
70795fce210SBarry Smith 
7084165533cSJose E. Roman   Input Parameters:
70995fce210SBarry Smith + sf  - communication object to duplicate
710cab54364SBarry Smith - opt - `PETSCSF_DUPLICATE_CONFONLY`, `PETSCSF_DUPLICATE_RANKS`, or `PETSCSF_DUPLICATE_GRAPH` (see `PetscSFDuplicateOption`)
71195fce210SBarry Smith 
7124165533cSJose E. Roman   Output Parameter:
71395fce210SBarry Smith . newsf - new communication object
71495fce210SBarry Smith 
71595fce210SBarry Smith   Level: beginner
71695fce210SBarry Smith 
71720662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFCreate()`, `PetscSFSetType()`, `PetscSFSetGraph()`
71895fce210SBarry Smith @*/
719d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFDuplicate(PetscSF sf, PetscSFDuplicateOption opt, PetscSF *newsf)
720d71ae5a4SJacob Faibussowitsch {
72129046d53SLisandro Dalcin   PetscSFType  type;
72297929ea7SJunchao Zhang   MPI_Datatype dtype = MPIU_SCALAR;
72395fce210SBarry Smith 
72495fce210SBarry Smith   PetscFunctionBegin;
72529046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
72629046d53SLisandro Dalcin   PetscValidLogicalCollectiveEnum(sf, opt, 2);
7274f572ea9SToby Isaac   PetscAssertPointer(newsf, 3);
7289566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sf), newsf));
7299566063dSJacob Faibussowitsch   PetscCall(PetscSFGetType(sf, &type));
7309566063dSJacob Faibussowitsch   if (type) PetscCall(PetscSFSetType(*newsf, type));
73135cb6cd3SPierre Jolivet   (*newsf)->allow_multi_leaves = sf->allow_multi_leaves; /* Dup this flag earlier since PetscSFSetGraph() below checks on this flag */
73295fce210SBarry Smith   if (opt == PETSCSF_DUPLICATE_GRAPH) {
733dd5b3ca6SJunchao Zhang     PetscSFCheckGraphSet(sf, 1);
734dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
73595fce210SBarry Smith       PetscInt           nroots, nleaves;
73695fce210SBarry Smith       const PetscInt    *ilocal;
73795fce210SBarry Smith       const PetscSFNode *iremote;
7389566063dSJacob Faibussowitsch       PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, &iremote));
7399566063dSJacob Faibussowitsch       PetscCall(PetscSFSetGraph(*newsf, nroots, nleaves, (PetscInt *)ilocal, PETSC_COPY_VALUES, (PetscSFNode *)iremote, PETSC_COPY_VALUES));
740dd5b3ca6SJunchao Zhang     } else {
7419566063dSJacob Faibussowitsch       PetscCall(PetscSFSetGraphWithPattern(*newsf, sf->map, sf->pattern));
742dd5b3ca6SJunchao Zhang     }
74395fce210SBarry Smith   }
74497929ea7SJunchao Zhang   /* Since oldtype is committed, so is newtype, according to MPI */
7459566063dSJacob Faibussowitsch   if (sf->vscat.bs > 1) PetscCallMPI(MPI_Type_dup(sf->vscat.unit, &dtype));
74697929ea7SJunchao Zhang   (*newsf)->vscat.bs     = sf->vscat.bs;
74797929ea7SJunchao Zhang   (*newsf)->vscat.unit   = dtype;
74897929ea7SJunchao Zhang   (*newsf)->vscat.to_n   = sf->vscat.to_n;
74997929ea7SJunchao Zhang   (*newsf)->vscat.from_n = sf->vscat.from_n;
75097929ea7SJunchao Zhang   /* Do not copy lsf. Build it on demand since it is rarely used */
75197929ea7SJunchao Zhang 
75220c24465SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
75320c24465SJunchao Zhang   (*newsf)->backend              = sf->backend;
75471438e86SJunchao Zhang   (*newsf)->unknown_input_stream = sf->unknown_input_stream;
75520c24465SJunchao Zhang   (*newsf)->use_gpu_aware_mpi    = sf->use_gpu_aware_mpi;
75620c24465SJunchao Zhang   (*newsf)->use_stream_aware_mpi = sf->use_stream_aware_mpi;
75720c24465SJunchao Zhang #endif
758dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf, Duplicate, opt, *newsf);
75920c24465SJunchao Zhang   /* Don't do PetscSFSetUp() since the new sf's graph might have not been set. */
7603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76195fce210SBarry Smith }
76295fce210SBarry Smith 
76395fce210SBarry Smith /*@C
76495fce210SBarry Smith   PetscSFGetGraph - Get the graph specifying a parallel star forest
76595fce210SBarry Smith 
76695fce210SBarry Smith   Not Collective
76795fce210SBarry Smith 
7684165533cSJose E. Roman   Input Parameter:
76995fce210SBarry Smith . sf - star forest
77095fce210SBarry Smith 
7714165533cSJose E. Roman   Output Parameters:
77295fce210SBarry Smith + nroots  - number of root vertices on the current process (these are possible targets for other process to attach leaves)
77395fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process
77420662ed9SBarry Smith . ilocal  - locations of leaves in leafdata buffers (if returned value is `NULL`, it means leaves are in contiguous storage)
77595fce210SBarry Smith - iremote - remote locations of root vertices for each leaf on the current process
77695fce210SBarry Smith 
777cab54364SBarry Smith   Level: intermediate
778cab54364SBarry Smith 
779373e0d91SLisandro Dalcin   Notes:
78020662ed9SBarry Smith   We are not currently requiring that the graph is set, thus returning `nroots` = -1 if it has not been set yet
781373e0d91SLisandro Dalcin 
78220662ed9SBarry Smith   The returned `ilocal` and `iremote` might contain values in different order than the input ones in `PetscSFSetGraph()`
783db2b9530SVaclav Hapla 
7848dbb0df6SBarry Smith   Fortran Notes:
78520662ed9SBarry Smith   The returned `iremote` array is a copy and must be deallocated after use. Consequently, if you
78620662ed9SBarry Smith   want to update the graph, you must call `PetscSFSetGraph()` after modifying the `iremote` array.
7878dbb0df6SBarry Smith 
78820662ed9SBarry Smith   To check for a `NULL` `ilocal` use
7898dbb0df6SBarry Smith $      if (loc(ilocal) == loc(PETSC_NULL_INTEGER)) then
790ca797d7aSLawrence Mitchell 
79120662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`
79295fce210SBarry Smith @*/
793d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGraph(PetscSF sf, PetscInt *nroots, PetscInt *nleaves, const PetscInt **ilocal, const PetscSFNode **iremote)
794d71ae5a4SJacob Faibussowitsch {
79595fce210SBarry Smith   PetscFunctionBegin;
79695fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
797b8dee149SJunchao Zhang   if (sf->ops->GetGraph) {
798*f4f49eeaSPierre Jolivet     PetscCall(sf->ops->GetGraph(sf, nroots, nleaves, ilocal, iremote));
799b8dee149SJunchao Zhang   } else {
80095fce210SBarry Smith     if (nroots) *nroots = sf->nroots;
80195fce210SBarry Smith     if (nleaves) *nleaves = sf->nleaves;
80295fce210SBarry Smith     if (ilocal) *ilocal = sf->mine;
80395fce210SBarry Smith     if (iremote) *iremote = sf->remote;
804b8dee149SJunchao Zhang   }
8053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80695fce210SBarry Smith }
80795fce210SBarry Smith 
80829046d53SLisandro Dalcin /*@
80995fce210SBarry Smith   PetscSFGetLeafRange - Get the active leaf ranges
81095fce210SBarry Smith 
81195fce210SBarry Smith   Not Collective
81295fce210SBarry Smith 
8134165533cSJose E. Roman   Input Parameter:
81495fce210SBarry Smith . sf - star forest
81595fce210SBarry Smith 
8164165533cSJose E. Roman   Output Parameters:
81720662ed9SBarry Smith + minleaf - minimum active leaf on this process. Returns 0 if there are no leaves.
81820662ed9SBarry Smith - maxleaf - maximum active leaf on this process. Returns -1 if there are no leaves.
81995fce210SBarry Smith 
82095fce210SBarry Smith   Level: developer
82195fce210SBarry Smith 
82220662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
82395fce210SBarry Smith @*/
824d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetLeafRange(PetscSF sf, PetscInt *minleaf, PetscInt *maxleaf)
825d71ae5a4SJacob Faibussowitsch {
82695fce210SBarry Smith   PetscFunctionBegin;
82795fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
82829046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf, 1);
82995fce210SBarry Smith   if (minleaf) *minleaf = sf->minleaf;
83095fce210SBarry Smith   if (maxleaf) *maxleaf = sf->maxleaf;
8313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
83295fce210SBarry Smith }
83395fce210SBarry Smith 
83495fce210SBarry Smith /*@C
835cab54364SBarry Smith   PetscSFViewFromOptions - View a `PetscSF` based on arguments in the options database
836fe2efc57SMark 
83720f4b53cSBarry Smith   Collective
838fe2efc57SMark 
839fe2efc57SMark   Input Parameters:
840fe2efc57SMark + A    - the star forest
841cab54364SBarry Smith . obj  - Optional object that provides the prefix for the option names
842736c3998SJose E. Roman - name - command line option
843fe2efc57SMark 
844fe2efc57SMark   Level: intermediate
845cab54364SBarry Smith 
84620662ed9SBarry Smith   Note:
84720662ed9SBarry Smith   See `PetscObjectViewFromOptions()` for possible `PetscViewer` and `PetscViewerFormat`
84820662ed9SBarry Smith 
849db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFView`, `PetscObjectViewFromOptions()`, `PetscSFCreate()`
850fe2efc57SMark @*/
851d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFViewFromOptions(PetscSF A, PetscObject obj, const char name[])
852d71ae5a4SJacob Faibussowitsch {
853fe2efc57SMark   PetscFunctionBegin;
854fe2efc57SMark   PetscValidHeaderSpecific(A, PETSCSF_CLASSID, 1);
8559566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
8563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
857fe2efc57SMark }
858fe2efc57SMark 
859fe2efc57SMark /*@C
86095fce210SBarry Smith   PetscSFView - view a star forest
86195fce210SBarry Smith 
86295fce210SBarry Smith   Collective
86395fce210SBarry Smith 
8644165533cSJose E. Roman   Input Parameters:
86595fce210SBarry Smith + sf     - star forest
866cab54364SBarry Smith - viewer - viewer to display graph, for example `PETSC_VIEWER_STDOUT_WORLD`
86795fce210SBarry Smith 
86895fce210SBarry Smith   Level: beginner
86995fce210SBarry Smith 
870cab54364SBarry Smith .seealso: `PetscSF`, `PetscViewer`, `PetscSFCreate()`, `PetscSFSetGraph()`
87195fce210SBarry Smith @*/
872d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFView(PetscSF sf, PetscViewer viewer)
873d71ae5a4SJacob Faibussowitsch {
87495fce210SBarry Smith   PetscBool         iascii;
87595fce210SBarry Smith   PetscViewerFormat format;
87695fce210SBarry Smith 
87795fce210SBarry Smith   PetscFunctionBegin;
87895fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
8799566063dSJacob Faibussowitsch   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf), &viewer));
88095fce210SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
88195fce210SBarry Smith   PetscCheckSameComm(sf, 1, viewer, 2);
8829566063dSJacob Faibussowitsch   if (sf->graphset) PetscCall(PetscSFSetUp(sf));
8839566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
88453dd6d7dSJunchao Zhang   if (iascii && viewer->format != PETSC_VIEWER_ASCII_MATLAB) {
88595fce210SBarry Smith     PetscMPIInt rank;
88681bfa7aaSJed Brown     PetscInt    ii, i, j;
88795fce210SBarry Smith 
8889566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)sf, viewer));
8899566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
890dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
89180153354SVaclav Hapla       if (!sf->graphset) {
8929566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "PetscSFSetGraph() has not been called yet\n"));
8939566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
8943ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
89580153354SVaclav Hapla       }
8969566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
8979566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
8989566063dSJacob 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));
89948a46eb9SPierre 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));
9009566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
9019566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetFormat(viewer, &format));
90295fce210SBarry Smith       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
90381bfa7aaSJed Brown         PetscMPIInt *tmpranks, *perm;
9049566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(sf->nranks, &tmpranks, sf->nranks, &perm));
9059566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(tmpranks, sf->ranks, sf->nranks));
90681bfa7aaSJed Brown         for (i = 0; i < sf->nranks; i++) perm[i] = i;
9079566063dSJacob Faibussowitsch         PetscCall(PetscSortMPIIntWithArray(sf->nranks, tmpranks, perm));
9089566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Roots referenced by my leaves, by rank\n", rank));
90981bfa7aaSJed Brown         for (ii = 0; ii < sf->nranks; ii++) {
91081bfa7aaSJed Brown           i = perm[ii];
9119566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %d: %" PetscInt_FMT " edges\n", rank, sf->ranks[i], sf->roffset[i + 1] - sf->roffset[i]));
91248a46eb9SPierre 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]));
91395fce210SBarry Smith         }
9149566063dSJacob Faibussowitsch         PetscCall(PetscFree2(tmpranks, perm));
91595fce210SBarry Smith       }
9169566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
9179566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
918dd5b3ca6SJunchao Zhang     }
9199566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
92095fce210SBarry Smith   }
921dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf, View, viewer);
9223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92395fce210SBarry Smith }
92495fce210SBarry Smith 
92595fce210SBarry Smith /*@C
926dec1416fSJunchao Zhang   PetscSFGetRootRanks - Get root ranks and number of vertices referenced by leaves on this process
92795fce210SBarry Smith 
92895fce210SBarry Smith   Not Collective
92995fce210SBarry Smith 
9304165533cSJose E. Roman   Input Parameter:
93195fce210SBarry Smith . sf - star forest
93295fce210SBarry Smith 
9334165533cSJose E. Roman   Output Parameters:
93495fce210SBarry Smith + nranks  - number of ranks referenced by local part
93520662ed9SBarry Smith . ranks   - [`nranks`] array of ranks
93620662ed9SBarry Smith . roffset - [`nranks`+1] offset in `rmine`/`rremote` for each rank
93720662ed9SBarry Smith . rmine   - [`roffset`[`nranks`]] concatenated array holding local indices referencing each remote rank
93820662ed9SBarry Smith - rremote - [`roffset`[`nranks`]] concatenated array holding remote indices referenced for each remote rank
93995fce210SBarry Smith 
94095fce210SBarry Smith   Level: developer
94195fce210SBarry Smith 
942cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGetLeafRanks()`
94395fce210SBarry Smith @*/
944d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetRootRanks(PetscSF sf, PetscInt *nranks, const PetscMPIInt **ranks, const PetscInt **roffset, const PetscInt **rmine, const PetscInt **rremote)
945d71ae5a4SJacob Faibussowitsch {
94695fce210SBarry Smith   PetscFunctionBegin;
94795fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
94828b400f6SJacob Faibussowitsch   PetscCheck(sf->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFSetUp() before obtaining ranks");
949dec1416fSJunchao Zhang   if (sf->ops->GetRootRanks) {
9509927e4dfSBarry Smith     PetscUseTypeMethod(sf, GetRootRanks, nranks, ranks, roffset, rmine, rremote);
951dec1416fSJunchao Zhang   } else {
952dec1416fSJunchao Zhang     /* The generic implementation */
95395fce210SBarry Smith     if (nranks) *nranks = sf->nranks;
95495fce210SBarry Smith     if (ranks) *ranks = sf->ranks;
95595fce210SBarry Smith     if (roffset) *roffset = sf->roffset;
95695fce210SBarry Smith     if (rmine) *rmine = sf->rmine;
95795fce210SBarry Smith     if (rremote) *rremote = sf->rremote;
958dec1416fSJunchao Zhang   }
9593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
96095fce210SBarry Smith }
96195fce210SBarry Smith 
9628750ddebSJunchao Zhang /*@C
9638750ddebSJunchao Zhang   PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process
9648750ddebSJunchao Zhang 
9658750ddebSJunchao Zhang   Not Collective
9668750ddebSJunchao Zhang 
9674165533cSJose E. Roman   Input Parameter:
9688750ddebSJunchao Zhang . sf - star forest
9698750ddebSJunchao Zhang 
9704165533cSJose E. Roman   Output Parameters:
9718750ddebSJunchao Zhang + niranks  - number of leaf ranks referencing roots on this process
97220662ed9SBarry Smith . iranks   - [`niranks`] array of ranks
97320662ed9SBarry Smith . ioffset  - [`niranks`+1] offset in `irootloc` for each rank
97420662ed9SBarry Smith - irootloc - [`ioffset`[`niranks`]] concatenated array holding local indices of roots referenced by each leaf rank
9758750ddebSJunchao Zhang 
9768750ddebSJunchao Zhang   Level: developer
9778750ddebSJunchao Zhang 
978cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGetRootRanks()`
9798750ddebSJunchao Zhang @*/
980d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetLeafRanks(PetscSF sf, PetscInt *niranks, const PetscMPIInt **iranks, const PetscInt **ioffset, const PetscInt **irootloc)
981d71ae5a4SJacob Faibussowitsch {
9828750ddebSJunchao Zhang   PetscFunctionBegin;
9838750ddebSJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
98428b400f6SJacob Faibussowitsch   PetscCheck(sf->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFSetUp() before obtaining ranks");
9858750ddebSJunchao Zhang   if (sf->ops->GetLeafRanks) {
9869927e4dfSBarry Smith     PetscUseTypeMethod(sf, GetLeafRanks, niranks, iranks, ioffset, irootloc);
9878750ddebSJunchao Zhang   } else {
9888750ddebSJunchao Zhang     PetscSFType type;
9899566063dSJacob Faibussowitsch     PetscCall(PetscSFGetType(sf, &type));
99098921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type);
9918750ddebSJunchao Zhang   }
9923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9938750ddebSJunchao Zhang }
9948750ddebSJunchao Zhang 
995d71ae5a4SJacob Faibussowitsch static PetscBool InList(PetscMPIInt needle, PetscMPIInt n, const PetscMPIInt *list)
996d71ae5a4SJacob Faibussowitsch {
997b5a8e515SJed Brown   PetscInt i;
998b5a8e515SJed Brown   for (i = 0; i < n; i++) {
999b5a8e515SJed Brown     if (needle == list[i]) return PETSC_TRUE;
1000b5a8e515SJed Brown   }
1001b5a8e515SJed Brown   return PETSC_FALSE;
1002b5a8e515SJed Brown }
1003b5a8e515SJed Brown 
100495fce210SBarry Smith /*@C
1005cab54364SBarry Smith   PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by `PetscSF` implementations.
100621c688dcSJed Brown 
100721c688dcSJed Brown   Collective
100821c688dcSJed Brown 
10094165533cSJose E. Roman   Input Parameters:
1010cab54364SBarry Smith + sf     - `PetscSF` to set up; `PetscSFSetGraph()` must have been called
1011cab54364SBarry Smith - dgroup - `MPI_Group` of ranks to be distinguished (e.g., for self or shared memory exchange)
101221c688dcSJed Brown 
101321c688dcSJed Brown   Level: developer
101421c688dcSJed Brown 
1015cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGetRootRanks()`
101621c688dcSJed Brown @*/
1017d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetUpRanks(PetscSF sf, MPI_Group dgroup)
1018d71ae5a4SJacob Faibussowitsch {
1019eec179cfSJacob Faibussowitsch   PetscHMapI    table;
1020eec179cfSJacob Faibussowitsch   PetscHashIter pos;
1021b5a8e515SJed Brown   PetscMPIInt   size, groupsize, *groupranks;
1022247e8311SStefano Zampini   PetscInt     *rcount, *ranks;
1023247e8311SStefano Zampini   PetscInt      i, irank = -1, orank = -1;
102421c688dcSJed Brown 
102521c688dcSJed Brown   PetscFunctionBegin;
102621c688dcSJed Brown   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
102729046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf, 1);
10289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf), &size));
1029eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapICreateWithSize(10, &table));
103021c688dcSJed Brown   for (i = 0; i < sf->nleaves; i++) {
103121c688dcSJed Brown     /* Log 1-based rank */
1032eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapISetWithMode(table, sf->remote[i].rank + 1, 1, ADD_VALUES));
103321c688dcSJed Brown   }
1034eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapIGetSize(table, &sf->nranks));
10359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(sf->nranks, &sf->ranks, sf->nranks + 1, &sf->roffset, sf->nleaves, &sf->rmine, sf->nleaves, &sf->rremote));
10369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(sf->nranks, &rcount, sf->nranks, &ranks));
1037eec179cfSJacob Faibussowitsch   PetscHashIterBegin(table, pos);
103821c688dcSJed Brown   for (i = 0; i < sf->nranks; i++) {
1039eec179cfSJacob Faibussowitsch     PetscHashIterGetKey(table, pos, ranks[i]);
1040eec179cfSJacob Faibussowitsch     PetscHashIterGetVal(table, pos, rcount[i]);
1041eec179cfSJacob Faibussowitsch     PetscHashIterNext(table, pos);
104221c688dcSJed Brown     ranks[i]--; /* Convert back to 0-based */
104321c688dcSJed Brown   }
1044eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&table));
1045b5a8e515SJed Brown 
1046b5a8e515SJed Brown   /* We expect that dgroup is reliably "small" while nranks could be large */
1047b5a8e515SJed Brown   {
10487fb8a5e4SKarl Rupp     MPI_Group    group = MPI_GROUP_NULL;
1049b5a8e515SJed Brown     PetscMPIInt *dgroupranks;
10509566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)sf), &group));
10519566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_size(dgroup, &groupsize));
10529566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(groupsize, &dgroupranks));
10539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(groupsize, &groupranks));
1054b5a8e515SJed Brown     for (i = 0; i < groupsize; i++) dgroupranks[i] = i;
10559566063dSJacob Faibussowitsch     if (groupsize) PetscCallMPI(MPI_Group_translate_ranks(dgroup, groupsize, dgroupranks, group, groupranks));
10569566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_free(&group));
10579566063dSJacob Faibussowitsch     PetscCall(PetscFree(dgroupranks));
1058b5a8e515SJed Brown   }
1059b5a8e515SJed Brown 
1060b5a8e515SJed Brown   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
1061b5a8e515SJed Brown   for (sf->ndranks = 0, i = sf->nranks; sf->ndranks < i;) {
1062b5a8e515SJed Brown     for (i--; sf->ndranks < i; i--) { /* Scan i backward looking for distinguished rank */
1063b5a8e515SJed Brown       if (InList(ranks[i], groupsize, groupranks)) break;
1064b5a8e515SJed Brown     }
1065b5a8e515SJed Brown     for (; sf->ndranks <= i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
1066b5a8e515SJed Brown       if (!InList(ranks[sf->ndranks], groupsize, groupranks)) break;
1067b5a8e515SJed Brown     }
1068b5a8e515SJed Brown     if (sf->ndranks < i) { /* Swap ranks[sf->ndranks] with ranks[i] */
1069b5a8e515SJed Brown       PetscInt tmprank, tmpcount;
1070247e8311SStefano Zampini 
1071b5a8e515SJed Brown       tmprank             = ranks[i];
1072b5a8e515SJed Brown       tmpcount            = rcount[i];
1073b5a8e515SJed Brown       ranks[i]            = ranks[sf->ndranks];
1074b5a8e515SJed Brown       rcount[i]           = rcount[sf->ndranks];
1075b5a8e515SJed Brown       ranks[sf->ndranks]  = tmprank;
1076b5a8e515SJed Brown       rcount[sf->ndranks] = tmpcount;
1077b5a8e515SJed Brown       sf->ndranks++;
1078b5a8e515SJed Brown     }
1079b5a8e515SJed Brown   }
10809566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupranks));
10819566063dSJacob Faibussowitsch   PetscCall(PetscSortIntWithArray(sf->ndranks, ranks, rcount));
10825c0db29aSPierre Jolivet   if (rcount) PetscCall(PetscSortIntWithArray(sf->nranks - sf->ndranks, ranks + sf->ndranks, rcount + sf->ndranks));
108321c688dcSJed Brown   sf->roffset[0] = 0;
108421c688dcSJed Brown   for (i = 0; i < sf->nranks; i++) {
10859566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(ranks[i], sf->ranks + i));
108621c688dcSJed Brown     sf->roffset[i + 1] = sf->roffset[i] + rcount[i];
108721c688dcSJed Brown     rcount[i]          = 0;
108821c688dcSJed Brown   }
1089247e8311SStefano Zampini   for (i = 0, irank = -1, orank = -1; i < sf->nleaves; i++) {
1090247e8311SStefano Zampini     /* short circuit */
1091247e8311SStefano Zampini     if (orank != sf->remote[i].rank) {
109221c688dcSJed Brown       /* Search for index of iremote[i].rank in sf->ranks */
10939566063dSJacob Faibussowitsch       PetscCall(PetscFindMPIInt(sf->remote[i].rank, sf->ndranks, sf->ranks, &irank));
1094b5a8e515SJed Brown       if (irank < 0) {
10959566063dSJacob Faibussowitsch         PetscCall(PetscFindMPIInt(sf->remote[i].rank, sf->nranks - sf->ndranks, sf->ranks + sf->ndranks, &irank));
1096b5a8e515SJed Brown         if (irank >= 0) irank += sf->ndranks;
109721c688dcSJed Brown       }
1098247e8311SStefano Zampini       orank = sf->remote[i].rank;
1099247e8311SStefano Zampini     }
110008401ef6SPierre Jolivet     PetscCheck(irank >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find rank %" PetscInt_FMT " in array", sf->remote[i].rank);
110121c688dcSJed Brown     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
110221c688dcSJed Brown     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
110321c688dcSJed Brown     rcount[irank]++;
110421c688dcSJed Brown   }
11059566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rcount, ranks));
11063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
110721c688dcSJed Brown }
110821c688dcSJed Brown 
110921c688dcSJed Brown /*@C
111095fce210SBarry Smith   PetscSFGetGroups - gets incoming and outgoing process groups
111195fce210SBarry Smith 
111295fce210SBarry Smith   Collective
111395fce210SBarry Smith 
11144165533cSJose E. Roman   Input Parameter:
111595fce210SBarry Smith . sf - star forest
111695fce210SBarry Smith 
11174165533cSJose E. Roman   Output Parameters:
111895fce210SBarry Smith + incoming - group of origin processes for incoming edges (leaves that reference my roots)
111995fce210SBarry Smith - outgoing - group of destination processes for outgoing edges (roots that I reference)
112095fce210SBarry Smith 
112195fce210SBarry Smith   Level: developer
112295fce210SBarry Smith 
1123cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGetWindow()`, `PetscSFRestoreWindow()`
112495fce210SBarry Smith @*/
1125d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGroups(PetscSF sf, MPI_Group *incoming, MPI_Group *outgoing)
1126d71ae5a4SJacob Faibussowitsch {
11277fb8a5e4SKarl Rupp   MPI_Group group = MPI_GROUP_NULL;
112895fce210SBarry Smith 
112995fce210SBarry Smith   PetscFunctionBegin;
113008401ef6SPierre Jolivet   PetscCheck(sf->nranks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFSetUpRanks() before obtaining groups");
113195fce210SBarry Smith   if (sf->ingroup == MPI_GROUP_NULL) {
113295fce210SBarry Smith     PetscInt        i;
113395fce210SBarry Smith     const PetscInt *indegree;
113495fce210SBarry Smith     PetscMPIInt     rank, *outranks, *inranks;
113595fce210SBarry Smith     PetscSFNode    *remote;
113695fce210SBarry Smith     PetscSF         bgcount;
113795fce210SBarry Smith 
113895fce210SBarry Smith     /* Compute the number of incoming ranks */
11399566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(sf->nranks, &remote));
114095fce210SBarry Smith     for (i = 0; i < sf->nranks; i++) {
114195fce210SBarry Smith       remote[i].rank  = sf->ranks[i];
114295fce210SBarry Smith       remote[i].index = 0;
114395fce210SBarry Smith     }
11449566063dSJacob Faibussowitsch     PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_CONFONLY, &bgcount));
11459566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(bgcount, 1, sf->nranks, NULL, PETSC_COPY_VALUES, remote, PETSC_OWN_POINTER));
11469566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(bgcount, &indegree));
11479566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(bgcount, &indegree));
114895fce210SBarry Smith     /* Enumerate the incoming ranks */
11499566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(indegree[0], &inranks, sf->nranks, &outranks));
11509566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
115195fce210SBarry Smith     for (i = 0; i < sf->nranks; i++) outranks[i] = rank;
11529566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherBegin(bgcount, MPI_INT, outranks, inranks));
11539566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherEnd(bgcount, MPI_INT, outranks, inranks));
11549566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)sf), &group));
11559566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_incl(group, indegree[0], inranks, &sf->ingroup));
11569566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_free(&group));
11579566063dSJacob Faibussowitsch     PetscCall(PetscFree2(inranks, outranks));
11589566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&bgcount));
115995fce210SBarry Smith   }
116095fce210SBarry Smith   *incoming = sf->ingroup;
116195fce210SBarry Smith 
116295fce210SBarry Smith   if (sf->outgroup == MPI_GROUP_NULL) {
11639566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)sf), &group));
11649566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_incl(group, sf->nranks, sf->ranks, &sf->outgroup));
11659566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_free(&group));
116695fce210SBarry Smith   }
116795fce210SBarry Smith   *outgoing = sf->outgroup;
11683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
116995fce210SBarry Smith }
117095fce210SBarry Smith 
117129046d53SLisandro Dalcin /*@
1172cab54364SBarry Smith   PetscSFGetMultiSF - gets the inner `PetscSF` implementing gathers and scatters
117395fce210SBarry Smith 
117495fce210SBarry Smith   Collective
117595fce210SBarry Smith 
11764165533cSJose E. Roman   Input Parameter:
117795fce210SBarry Smith . sf - star forest that may contain roots with 0 or with more than 1 vertex
117895fce210SBarry Smith 
11794165533cSJose E. Roman   Output Parameter:
118095fce210SBarry Smith . multi - star forest with split roots, such that each root has degree exactly 1
118195fce210SBarry Smith 
118295fce210SBarry Smith   Level: developer
118395fce210SBarry Smith 
1184cab54364SBarry Smith   Note:
1185cab54364SBarry Smith   In most cases, users should use `PetscSFGatherBegin()` and `PetscSFScatterBegin()` instead of manipulating multi
118695fce210SBarry Smith   directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
118795fce210SBarry Smith   edge, it is a candidate for future optimization that might involve its removal.
118895fce210SBarry Smith 
1189cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFGatherBegin()`, `PetscSFScatterBegin()`, `PetscSFComputeMultiRootOriginalNumbering()`
119095fce210SBarry Smith @*/
1191d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetMultiSF(PetscSF sf, PetscSF *multi)
1192d71ae5a4SJacob Faibussowitsch {
119395fce210SBarry Smith   PetscFunctionBegin;
119495fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
11954f572ea9SToby Isaac   PetscAssertPointer(multi, 2);
119695fce210SBarry Smith   if (sf->nroots < 0) { /* Graph has not been set yet; why do we need this? */
11979566063dSJacob Faibussowitsch     PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, &sf->multi));
119895fce210SBarry Smith     *multi           = sf->multi;
1199013b3241SStefano Zampini     sf->multi->multi = sf->multi;
12003ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
120195fce210SBarry Smith   }
120295fce210SBarry Smith   if (!sf->multi) {
120395fce210SBarry Smith     const PetscInt *indegree;
12049837ea96SMatthew G. Knepley     PetscInt        i, *inoffset, *outones, *outoffset, maxlocal;
120595fce210SBarry Smith     PetscSFNode    *remote;
120629046d53SLisandro Dalcin     maxlocal = sf->maxleaf + 1; /* TODO: We should use PetscSFGetLeafRange() */
12079566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf, &indegree));
12089566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf, &indegree));
12099566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(sf->nroots + 1, &inoffset, maxlocal, &outones, maxlocal, &outoffset));
121095fce210SBarry Smith     inoffset[0] = 0;
121195fce210SBarry Smith     for (i = 0; i < sf->nroots; i++) inoffset[i + 1] = inoffset[i] + indegree[i];
12129837ea96SMatthew G. Knepley     for (i = 0; i < maxlocal; i++) outones[i] = 1;
12139566063dSJacob Faibussowitsch     PetscCall(PetscSFFetchAndOpBegin(sf, MPIU_INT, inoffset, outones, outoffset, MPI_SUM));
12149566063dSJacob Faibussowitsch     PetscCall(PetscSFFetchAndOpEnd(sf, MPIU_INT, inoffset, outones, outoffset, MPI_SUM));
121595fce210SBarry Smith     for (i = 0; i < sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
121676bd3646SJed Brown     if (PetscDefined(USE_DEBUG)) {                               /* Check that the expected number of increments occurred */
1217ad540459SPierre 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");
121876bd3646SJed Brown     }
12199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(sf->nleaves, &remote));
122095fce210SBarry Smith     for (i = 0; i < sf->nleaves; i++) {
122195fce210SBarry Smith       remote[i].rank  = sf->remote[i].rank;
122238e7336fSToby Isaac       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
122395fce210SBarry Smith     }
12249566063dSJacob Faibussowitsch     PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, &sf->multi));
1225013b3241SStefano Zampini     sf->multi->multi = sf->multi;
12269566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(sf->multi, inoffset[sf->nroots], sf->nleaves, sf->mine, PETSC_COPY_VALUES, remote, PETSC_OWN_POINTER));
122795fce210SBarry Smith     if (sf->rankorder) { /* Sort the ranks */
122895fce210SBarry Smith       PetscMPIInt  rank;
122995fce210SBarry Smith       PetscInt    *inranks, *newoffset, *outranks, *newoutoffset, *tmpoffset, maxdegree;
123095fce210SBarry Smith       PetscSFNode *newremote;
12319566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
123295fce210SBarry Smith       for (i = 0, maxdegree = 0; i < sf->nroots; i++) maxdegree = PetscMax(maxdegree, indegree[i]);
12339566063dSJacob Faibussowitsch       PetscCall(PetscMalloc5(sf->multi->nroots, &inranks, sf->multi->nroots, &newoffset, maxlocal, &outranks, maxlocal, &newoutoffset, maxdegree, &tmpoffset));
12349837ea96SMatthew G. Knepley       for (i = 0; i < maxlocal; i++) outranks[i] = rank;
12359566063dSJacob Faibussowitsch       PetscCall(PetscSFReduceBegin(sf->multi, MPIU_INT, outranks, inranks, MPI_REPLACE));
12369566063dSJacob Faibussowitsch       PetscCall(PetscSFReduceEnd(sf->multi, MPIU_INT, outranks, inranks, MPI_REPLACE));
123795fce210SBarry Smith       /* Sort the incoming ranks at each vertex, build the inverse map */
123895fce210SBarry Smith       for (i = 0; i < sf->nroots; i++) {
123995fce210SBarry Smith         PetscInt j;
124095fce210SBarry Smith         for (j = 0; j < indegree[i]; j++) tmpoffset[j] = j;
12418e3a54c0SPierre Jolivet         PetscCall(PetscSortIntWithArray(indegree[i], PetscSafePointerPlusOffset(inranks, inoffset[i]), tmpoffset));
124295fce210SBarry Smith         for (j = 0; j < indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
124395fce210SBarry Smith       }
12449566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(sf->multi, MPIU_INT, newoffset, newoutoffset, MPI_REPLACE));
12459566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(sf->multi, MPIU_INT, newoffset, newoutoffset, MPI_REPLACE));
12469566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(sf->nleaves, &newremote));
124795fce210SBarry Smith       for (i = 0; i < sf->nleaves; i++) {
124895fce210SBarry Smith         newremote[i].rank  = sf->remote[i].rank;
124901365b40SToby Isaac         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
125095fce210SBarry Smith       }
12519566063dSJacob Faibussowitsch       PetscCall(PetscSFSetGraph(sf->multi, inoffset[sf->nroots], sf->nleaves, sf->mine, PETSC_COPY_VALUES, newremote, PETSC_OWN_POINTER));
12529566063dSJacob Faibussowitsch       PetscCall(PetscFree5(inranks, newoffset, outranks, newoutoffset, tmpoffset));
125395fce210SBarry Smith     }
12549566063dSJacob Faibussowitsch     PetscCall(PetscFree3(inoffset, outones, outoffset));
125595fce210SBarry Smith   }
125695fce210SBarry Smith   *multi = sf->multi;
12573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
125895fce210SBarry Smith }
125995fce210SBarry Smith 
126095fce210SBarry Smith /*@C
126120662ed9SBarry Smith   PetscSFCreateEmbeddedRootSF - removes edges from all but the selected roots of a `PetscSF`, does not remap indices
126295fce210SBarry Smith 
126395fce210SBarry Smith   Collective
126495fce210SBarry Smith 
12654165533cSJose E. Roman   Input Parameters:
126695fce210SBarry Smith + sf        - original star forest
1267ba2a7774SJunchao Zhang . nselected - number of selected roots on this process
1268ba2a7774SJunchao Zhang - selected  - indices of the selected roots on this process
126995fce210SBarry Smith 
12704165533cSJose E. Roman   Output Parameter:
1271cd620004SJunchao Zhang . esf - new star forest
127295fce210SBarry Smith 
127395fce210SBarry Smith   Level: advanced
127495fce210SBarry Smith 
127595fce210SBarry Smith   Note:
1276cab54364SBarry Smith   To use the new `PetscSF`, it may be necessary to know the indices of the leaves that are still participating. This can
127795fce210SBarry Smith   be done by calling PetscSFGetGraph().
127895fce210SBarry Smith 
1279cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
128095fce210SBarry Smith @*/
1281d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateEmbeddedRootSF(PetscSF sf, PetscInt nselected, const PetscInt *selected, PetscSF *esf)
1282d71ae5a4SJacob Faibussowitsch {
1283cd620004SJunchao Zhang   PetscInt           i, j, n, nroots, nleaves, esf_nleaves, *new_ilocal, minleaf, maxleaf, maxlocal;
1284cd620004SJunchao Zhang   const PetscInt    *ilocal;
1285cd620004SJunchao Zhang   signed char       *rootdata, *leafdata, *leafmem;
1286ba2a7774SJunchao Zhang   const PetscSFNode *iremote;
1287f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1288f659e5c7SJunchao Zhang   MPI_Comm           comm;
128995fce210SBarry Smith 
129095fce210SBarry Smith   PetscFunctionBegin;
129195fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
129229046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf, 1);
12934f572ea9SToby Isaac   if (nselected) PetscAssertPointer(selected, 3);
12944f572ea9SToby Isaac   PetscAssertPointer(esf, 4);
12950511a646SMatthew G. Knepley 
12969566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
12979566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_EmbedSF, sf, 0, 0, 0));
12989566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
12999566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, &iremote));
1300cd620004SJunchao Zhang 
130176bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) { /* Error out if selected[] has dups or out of range indices */
1302cd620004SJunchao Zhang     PetscBool dups;
13039566063dSJacob Faibussowitsch     PetscCall(PetscCheckDupsInt(nselected, selected, &dups));
130428b400f6SJacob Faibussowitsch     PetscCheck(!dups, comm, PETSC_ERR_ARG_WRONG, "selected[] has dups");
1305511e6246SStefano Zampini     for (i = 0; i < nselected; i++) PetscCheck(selected[i] >= 0 && selected[i] < nroots, comm, PETSC_ERR_ARG_OUTOFRANGE, "selected root index %" PetscInt_FMT " is out of [0,%" PetscInt_FMT ")", selected[i], nroots);
1306cd620004SJunchao Zhang   }
1307f659e5c7SJunchao Zhang 
1308dbbe0bcdSBarry Smith   if (sf->ops->CreateEmbeddedRootSF) PetscUseTypeMethod(sf, CreateEmbeddedRootSF, nselected, selected, esf);
1309dbbe0bcdSBarry Smith   else {
1310cd620004SJunchao Zhang     /* A generic version of creating embedded sf */
13119566063dSJacob Faibussowitsch     PetscCall(PetscSFGetLeafRange(sf, &minleaf, &maxleaf));
1312cd620004SJunchao Zhang     maxlocal = maxleaf - minleaf + 1;
13139566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(nroots, &rootdata, maxlocal, &leafmem));
13148e3a54c0SPierre Jolivet     leafdata = PetscSafePointerPlusOffset(leafmem, -minleaf);
1315cd620004SJunchao Zhang     /* Tag selected roots and bcast to leaves */
1316cd620004SJunchao Zhang     for (i = 0; i < nselected; i++) rootdata[selected[i]] = 1;
13179566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPI_SIGNED_CHAR, rootdata, leafdata, MPI_REPLACE));
13189566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPI_SIGNED_CHAR, rootdata, leafdata, MPI_REPLACE));
1319ba2a7774SJunchao Zhang 
1320cd620004SJunchao Zhang     /* Build esf with leaves that are still connected */
1321cd620004SJunchao Zhang     esf_nleaves = 0;
1322cd620004SJunchao Zhang     for (i = 0; i < nleaves; i++) {
1323cd620004SJunchao Zhang       j = ilocal ? ilocal[i] : i;
1324cd620004SJunchao Zhang       /* esf_nleaves += leafdata[j] should work in theory, but failed with SFWindow bugs
1325cd620004SJunchao Zhang          with PetscSFBcast. See https://gitlab.com/petsc/petsc/issues/555
1326cd620004SJunchao Zhang       */
1327cd620004SJunchao Zhang       esf_nleaves += (leafdata[j] ? 1 : 0);
1328cd620004SJunchao Zhang     }
13299566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(esf_nleaves, &new_ilocal));
13309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(esf_nleaves, &new_iremote));
1331cd620004SJunchao Zhang     for (i = n = 0; i < nleaves; i++) {
1332cd620004SJunchao Zhang       j = ilocal ? ilocal[i] : i;
1333cd620004SJunchao Zhang       if (leafdata[j]) {
1334cd620004SJunchao Zhang         new_ilocal[n]        = j;
1335cd620004SJunchao Zhang         new_iremote[n].rank  = iremote[i].rank;
1336cd620004SJunchao Zhang         new_iremote[n].index = iremote[i].index;
1337fc1ede2bSMatthew G. Knepley         ++n;
133895fce210SBarry Smith       }
133995fce210SBarry Smith     }
13409566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(comm, esf));
13419566063dSJacob Faibussowitsch     PetscCall(PetscSFSetFromOptions(*esf));
13429566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*esf, nroots, esf_nleaves, new_ilocal, PETSC_OWN_POINTER, new_iremote, PETSC_OWN_POINTER));
13439566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafmem));
1344f659e5c7SJunchao Zhang   }
13459566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_EmbedSF, sf, 0, 0, 0));
13463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
134795fce210SBarry Smith }
134895fce210SBarry Smith 
13492f5fb4c2SMatthew G. Knepley /*@C
135020662ed9SBarry Smith   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves of a `PetscSF`, does not remap indices
13512f5fb4c2SMatthew G. Knepley 
13522f5fb4c2SMatthew G. Knepley   Collective
13532f5fb4c2SMatthew G. Knepley 
13544165533cSJose E. Roman   Input Parameters:
13552f5fb4c2SMatthew G. Knepley + sf        - original star forest
1356f659e5c7SJunchao Zhang . nselected - number of selected leaves on this process
1357f659e5c7SJunchao Zhang - selected  - indices of the selected leaves on this process
13582f5fb4c2SMatthew G. Knepley 
13594165533cSJose E. Roman   Output Parameter:
13602f5fb4c2SMatthew G. Knepley . newsf - new star forest
13612f5fb4c2SMatthew G. Knepley 
13622f5fb4c2SMatthew G. Knepley   Level: advanced
13632f5fb4c2SMatthew G. Knepley 
1364cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreateEmbeddedRootSF()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
13652f5fb4c2SMatthew G. Knepley @*/
1366d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nselected, const PetscInt *selected, PetscSF *newsf)
1367d71ae5a4SJacob Faibussowitsch {
1368f659e5c7SJunchao Zhang   const PetscSFNode *iremote;
1369f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1370f659e5c7SJunchao Zhang   const PetscInt    *ilocal;
1371f659e5c7SJunchao Zhang   PetscInt           i, nroots, *leaves, *new_ilocal;
1372f659e5c7SJunchao Zhang   MPI_Comm           comm;
13732f5fb4c2SMatthew G. Knepley 
13742f5fb4c2SMatthew G. Knepley   PetscFunctionBegin;
13752f5fb4c2SMatthew G. Knepley   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
137629046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf, 1);
13774f572ea9SToby Isaac   if (nselected) PetscAssertPointer(selected, 3);
13784f572ea9SToby Isaac   PetscAssertPointer(newsf, 4);
13792f5fb4c2SMatthew G. Knepley 
1380f659e5c7SJunchao Zhang   /* Uniq selected[] and put results in leaves[] */
13819566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
13829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nselected, &leaves));
13839566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(leaves, selected, nselected));
13849566063dSJacob Faibussowitsch   PetscCall(PetscSortedRemoveDupsInt(&nselected, leaves));
138508401ef6SPierre 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);
1386f659e5c7SJunchao Zhang 
1387f659e5c7SJunchao Zhang   /* Optimize the routine only when sf is setup and hence we can reuse sf's communication pattern */
1388dbbe0bcdSBarry Smith   if (sf->setupcalled && sf->ops->CreateEmbeddedLeafSF) PetscUseTypeMethod(sf, CreateEmbeddedLeafSF, nselected, leaves, newsf);
1389dbbe0bcdSBarry Smith   else {
13909566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sf, &nroots, NULL, &ilocal, &iremote));
13919566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nselected, &new_ilocal));
13929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nselected, &new_iremote));
1393f659e5c7SJunchao Zhang     for (i = 0; i < nselected; ++i) {
1394f659e5c7SJunchao Zhang       const PetscInt l     = leaves[i];
1395f659e5c7SJunchao Zhang       new_ilocal[i]        = ilocal ? ilocal[l] : l;
1396f659e5c7SJunchao Zhang       new_iremote[i].rank  = iremote[l].rank;
1397f659e5c7SJunchao Zhang       new_iremote[i].index = iremote[l].index;
13982f5fb4c2SMatthew G. Knepley     }
13999566063dSJacob Faibussowitsch     PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_CONFONLY, newsf));
14009566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*newsf, nroots, nselected, new_ilocal, PETSC_OWN_POINTER, new_iremote, PETSC_OWN_POINTER));
1401f659e5c7SJunchao Zhang   }
14029566063dSJacob Faibussowitsch   PetscCall(PetscFree(leaves));
14033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14042f5fb4c2SMatthew G. Knepley }
14052f5fb4c2SMatthew G. Knepley 
140695fce210SBarry Smith /*@C
1407cab54364SBarry Smith   PetscSFBcastBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to `PetscSFBcastEnd()`
14083482bfa8SJunchao Zhang 
1409c3339decSBarry Smith   Collective
14103482bfa8SJunchao Zhang 
14114165533cSJose E. Roman   Input Parameters:
14123482bfa8SJunchao Zhang + sf       - star forest on which to communicate
14133482bfa8SJunchao Zhang . unit     - data type associated with each node
14143482bfa8SJunchao Zhang . rootdata - buffer to broadcast
14153482bfa8SJunchao Zhang - op       - operation to use for reduction
14163482bfa8SJunchao Zhang 
14174165533cSJose E. Roman   Output Parameter:
14183482bfa8SJunchao Zhang . leafdata - buffer to be reduced with values from each leaf's respective root
14193482bfa8SJunchao Zhang 
14203482bfa8SJunchao Zhang   Level: intermediate
14213482bfa8SJunchao Zhang 
142220662ed9SBarry Smith   Note:
142320662ed9SBarry Smith   When PETSc is configured with device support, it will use its own mechanism to figure out whether the given data pointers
1424da81f932SPierre Jolivet   are host pointers or device pointers, which may incur a noticeable cost. If you already knew the info, you should
1425cab54364SBarry Smith   use `PetscSFBcastWithMemTypeBegin()` instead.
1426cab54364SBarry Smith 
1427cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFBcastEnd()`, `PetscSFBcastWithMemTypeBegin()`
14283482bfa8SJunchao Zhang @*/
1429d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFBcastBegin(PetscSF sf, MPI_Datatype unit, const void *rootdata, void *leafdata, MPI_Op op)
1430d71ae5a4SJacob Faibussowitsch {
1431eb02082bSJunchao Zhang   PetscMemType rootmtype, leafmtype;
14323482bfa8SJunchao Zhang 
14333482bfa8SJunchao Zhang   PetscFunctionBegin;
14343482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
14359566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
14369566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_BcastBegin, sf, 0, 0, 0));
14379566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(rootdata, &rootmtype));
14389566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafdata, &leafmtype));
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));
14413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14423482bfa8SJunchao Zhang }
14433482bfa8SJunchao Zhang 
14443482bfa8SJunchao Zhang /*@C
144520662ed9SBarry Smith   PetscSFBcastWithMemTypeBegin - begin pointwise broadcast with root value being reduced to leaf value with explicit memory types, to be concluded with call
144620662ed9SBarry Smith   to `PetscSFBcastEnd()`
1447d0295fc0SJunchao Zhang 
1448c3339decSBarry Smith   Collective
1449d0295fc0SJunchao Zhang 
14504165533cSJose E. Roman   Input Parameters:
1451d0295fc0SJunchao Zhang + sf        - star forest on which to communicate
1452d0295fc0SJunchao Zhang . unit      - data type associated with each node
1453d0295fc0SJunchao Zhang . rootmtype - memory type of rootdata
1454d0295fc0SJunchao Zhang . rootdata  - buffer to broadcast
1455d0295fc0SJunchao Zhang . leafmtype - memory type of leafdata
1456d0295fc0SJunchao Zhang - op        - operation to use for reduction
1457d0295fc0SJunchao Zhang 
14584165533cSJose E. Roman   Output Parameter:
1459d0295fc0SJunchao Zhang . leafdata - buffer to be reduced with values from each leaf's respective root
1460d0295fc0SJunchao Zhang 
1461d0295fc0SJunchao Zhang   Level: intermediate
1462d0295fc0SJunchao Zhang 
1463cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFBcastEnd()`, `PetscSFBcastBegin()`
1464d0295fc0SJunchao Zhang @*/
1465d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFBcastWithMemTypeBegin(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, const void *rootdata, PetscMemType leafmtype, void *leafdata, MPI_Op op)
1466d71ae5a4SJacob Faibussowitsch {
1467d0295fc0SJunchao Zhang   PetscFunctionBegin;
1468d0295fc0SJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
14699566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
14709566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_BcastBegin, sf, 0, 0, 0));
1471dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, BcastBegin, unit, rootmtype, rootdata, leafmtype, leafdata, op);
14729566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_BcastBegin, sf, 0, 0, 0));
14733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1474d0295fc0SJunchao Zhang }
1475d0295fc0SJunchao Zhang 
1476d0295fc0SJunchao Zhang /*@C
147720662ed9SBarry Smith   PetscSFBcastEnd - end a broadcast and reduce operation started with `PetscSFBcastBegin()` or `PetscSFBcastWithMemTypeBegin()`
14783482bfa8SJunchao Zhang 
14793482bfa8SJunchao Zhang   Collective
14803482bfa8SJunchao Zhang 
14814165533cSJose E. Roman   Input Parameters:
14823482bfa8SJunchao Zhang + sf       - star forest
14833482bfa8SJunchao Zhang . unit     - data type
14843482bfa8SJunchao Zhang . rootdata - buffer to broadcast
14853482bfa8SJunchao Zhang - op       - operation to use for reduction
14863482bfa8SJunchao Zhang 
14874165533cSJose E. Roman   Output Parameter:
14883482bfa8SJunchao Zhang . leafdata - buffer to be reduced with values from each leaf's respective root
14893482bfa8SJunchao Zhang 
14903482bfa8SJunchao Zhang   Level: intermediate
14913482bfa8SJunchao Zhang 
1492cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFReduceEnd()`
14933482bfa8SJunchao Zhang @*/
1494d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFBcastEnd(PetscSF sf, MPI_Datatype unit, const void *rootdata, void *leafdata, MPI_Op op)
1495d71ae5a4SJacob Faibussowitsch {
14963482bfa8SJunchao Zhang   PetscFunctionBegin;
14973482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
14989566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_BcastEnd, sf, 0, 0, 0));
1499dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, BcastEnd, unit, rootdata, leafdata, op);
15009566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_BcastEnd, sf, 0, 0, 0));
15013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15023482bfa8SJunchao Zhang }
15033482bfa8SJunchao Zhang 
15043482bfa8SJunchao Zhang /*@C
1505cab54364SBarry Smith   PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to `PetscSFReduceEnd()`
150695fce210SBarry Smith 
150795fce210SBarry Smith   Collective
150895fce210SBarry Smith 
15094165533cSJose E. Roman   Input Parameters:
151095fce210SBarry Smith + sf       - star forest
151195fce210SBarry Smith . unit     - data type
151295fce210SBarry Smith . leafdata - values to reduce
151395fce210SBarry Smith - op       - reduction operation
151495fce210SBarry Smith 
15154165533cSJose E. Roman   Output Parameter:
151695fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root
151795fce210SBarry Smith 
151895fce210SBarry Smith   Level: intermediate
151995fce210SBarry Smith 
152020662ed9SBarry Smith   Note:
152120662ed9SBarry Smith   When PETSc is configured with device support, it will use its own mechanism to figure out whether the given data pointers
1522da81f932SPierre Jolivet   are host pointers or device pointers, which may incur a noticeable cost. If you already knew the info, you should
1523cab54364SBarry Smith   use `PetscSFReduceWithMemTypeBegin()` instead.
1524d0295fc0SJunchao Zhang 
152520662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFBcastBegin()`, `PetscSFReduceWithMemTypeBegin()`, `PetscSFReduceEnd()`
152695fce210SBarry Smith @*/
1527d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFReduceBegin(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *rootdata, MPI_Op op)
1528d71ae5a4SJacob Faibussowitsch {
1529eb02082bSJunchao Zhang   PetscMemType rootmtype, leafmtype;
153095fce210SBarry Smith 
153195fce210SBarry Smith   PetscFunctionBegin;
153295fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
15339566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
15349566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_ReduceBegin, sf, 0, 0, 0));
15359566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(rootdata, &rootmtype));
15369566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafdata, &leafmtype));
1537*f4f49eeaSPierre Jolivet   PetscCall(sf->ops->ReduceBegin(sf, unit, leafmtype, leafdata, rootmtype, rootdata, op));
15389566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_ReduceBegin, sf, 0, 0, 0));
15393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
154095fce210SBarry Smith }
154195fce210SBarry Smith 
154295fce210SBarry Smith /*@C
1543cab54364SBarry Smith   PetscSFReduceWithMemTypeBegin - begin reduction of leafdata into rootdata with explicit memory types, to be completed with call to `PetscSFReduceEnd()`
1544d0295fc0SJunchao Zhang 
1545d0295fc0SJunchao Zhang   Collective
1546d0295fc0SJunchao Zhang 
15474165533cSJose E. Roman   Input Parameters:
1548d0295fc0SJunchao Zhang + sf        - star forest
1549d0295fc0SJunchao Zhang . unit      - data type
1550d0295fc0SJunchao Zhang . leafmtype - memory type of leafdata
1551d0295fc0SJunchao Zhang . leafdata  - values to reduce
1552d0295fc0SJunchao Zhang . rootmtype - memory type of rootdata
1553d0295fc0SJunchao Zhang - op        - reduction operation
1554d0295fc0SJunchao Zhang 
15554165533cSJose E. Roman   Output Parameter:
1556d0295fc0SJunchao Zhang . rootdata - result of reduction of values from all leaves of each root
1557d0295fc0SJunchao Zhang 
1558d0295fc0SJunchao Zhang   Level: intermediate
1559d0295fc0SJunchao Zhang 
156020662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFBcastBegin()`, `PetscSFReduceBegin()`, `PetscSFReduceEnd()`
1561d0295fc0SJunchao Zhang @*/
1562d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFReduceWithMemTypeBegin(PetscSF sf, MPI_Datatype unit, PetscMemType leafmtype, const void *leafdata, PetscMemType rootmtype, void *rootdata, MPI_Op op)
1563d71ae5a4SJacob Faibussowitsch {
1564d0295fc0SJunchao Zhang   PetscFunctionBegin;
1565d0295fc0SJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
15669566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
15679566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_ReduceBegin, sf, 0, 0, 0));
1568*f4f49eeaSPierre Jolivet   PetscCall(sf->ops->ReduceBegin(sf, unit, leafmtype, leafdata, rootmtype, rootdata, op));
15699566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_ReduceBegin, sf, 0, 0, 0));
15703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1571d0295fc0SJunchao Zhang }
1572d0295fc0SJunchao Zhang 
1573d0295fc0SJunchao Zhang /*@C
157420662ed9SBarry Smith   PetscSFReduceEnd - end a reduction operation started with `PetscSFReduceBegin()` or `PetscSFReduceWithMemTypeBegin()`
157595fce210SBarry Smith 
157695fce210SBarry Smith   Collective
157795fce210SBarry Smith 
15784165533cSJose E. Roman   Input Parameters:
157995fce210SBarry Smith + sf       - star forest
158095fce210SBarry Smith . unit     - data type
158195fce210SBarry Smith . leafdata - values to reduce
158295fce210SBarry Smith - op       - reduction operation
158395fce210SBarry Smith 
15844165533cSJose E. Roman   Output Parameter:
158595fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root
158695fce210SBarry Smith 
158795fce210SBarry Smith   Level: intermediate
158895fce210SBarry Smith 
158920662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFBcastEnd()`, `PetscSFReduceBegin()`, `PetscSFReduceWithMemTypeBegin()`
159095fce210SBarry Smith @*/
1591d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFReduceEnd(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *rootdata, MPI_Op op)
1592d71ae5a4SJacob Faibussowitsch {
159395fce210SBarry Smith   PetscFunctionBegin;
159495fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
15959566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_ReduceEnd, sf, 0, 0, 0));
1596dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, ReduceEnd, unit, leafdata, rootdata, op);
15979566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_ReduceEnd, sf, 0, 0, 0));
15983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
159995fce210SBarry Smith }
160095fce210SBarry Smith 
160195fce210SBarry Smith /*@C
1602cab54364SBarry Smith   PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value,
1603cab54364SBarry Smith   to be completed with `PetscSFFetchAndOpEnd()`
1604a1729e3fSJunchao Zhang 
1605a1729e3fSJunchao Zhang   Collective
1606a1729e3fSJunchao Zhang 
16074165533cSJose E. Roman   Input Parameters:
1608a1729e3fSJunchao Zhang + sf       - star forest
1609a1729e3fSJunchao Zhang . unit     - data type
1610a1729e3fSJunchao Zhang . leafdata - leaf values to use in reduction
1611a1729e3fSJunchao Zhang - op       - operation to use for reduction
1612a1729e3fSJunchao Zhang 
16134165533cSJose E. Roman   Output Parameters:
1614a1729e3fSJunchao Zhang + rootdata   - root values to be updated, input state is seen by first process to perform an update
1615a1729e3fSJunchao Zhang - leafupdate - state at each leaf's respective root immediately prior to my atomic update
1616a1729e3fSJunchao Zhang 
1617a1729e3fSJunchao Zhang   Level: advanced
1618a1729e3fSJunchao Zhang 
1619a1729e3fSJunchao Zhang   Note:
1620a1729e3fSJunchao Zhang   The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
1621a1729e3fSJunchao Zhang   might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
1622a1729e3fSJunchao Zhang   not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
1623a1729e3fSJunchao Zhang   integers.
1624a1729e3fSJunchao Zhang 
1625cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeBegin()`, `PetscSFReduceBegin()`, `PetscSFSetGraph()`
1626a1729e3fSJunchao Zhang @*/
1627d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf, MPI_Datatype unit, void *rootdata, const void *leafdata, void *leafupdate, MPI_Op op)
1628d71ae5a4SJacob Faibussowitsch {
1629eb02082bSJunchao Zhang   PetscMemType rootmtype, leafmtype, leafupdatemtype;
1630a1729e3fSJunchao Zhang 
1631a1729e3fSJunchao Zhang   PetscFunctionBegin;
1632a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
16339566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
16349566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_FetchAndOpBegin, sf, 0, 0, 0));
16359566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(rootdata, &rootmtype));
16369566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafdata, &leafmtype));
16379566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafupdate, &leafupdatemtype));
163808401ef6SPierre Jolivet   PetscCheck(leafmtype == leafupdatemtype, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for leafdata and leafupdate in different memory types");
1639dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, FetchAndOpBegin, unit, rootmtype, rootdata, leafmtype, leafdata, leafupdate, op);
16409566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_FetchAndOpBegin, sf, 0, 0, 0));
16413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1642a1729e3fSJunchao Zhang }
1643a1729e3fSJunchao Zhang 
1644a1729e3fSJunchao Zhang /*@C
1645cab54364SBarry Smith   PetscSFFetchAndOpWithMemTypeBegin - begin operation with explicit memory types that fetches values from root and updates atomically by
1646cab54364SBarry Smith   applying operation using my leaf value, to be completed with `PetscSFFetchAndOpEnd()`
1647d3b3e55cSJunchao Zhang 
1648d3b3e55cSJunchao Zhang   Collective
1649d3b3e55cSJunchao Zhang 
1650d3b3e55cSJunchao Zhang   Input Parameters:
1651d3b3e55cSJunchao Zhang + sf              - star forest
1652d3b3e55cSJunchao Zhang . unit            - data type
1653d3b3e55cSJunchao Zhang . rootmtype       - memory type of rootdata
1654d3b3e55cSJunchao Zhang . leafmtype       - memory type of leafdata
1655d3b3e55cSJunchao Zhang . leafdata        - leaf values to use in reduction
1656d3b3e55cSJunchao Zhang . leafupdatemtype - memory type of leafupdate
1657d3b3e55cSJunchao Zhang - op              - operation to use for reduction
1658d3b3e55cSJunchao Zhang 
1659d3b3e55cSJunchao Zhang   Output Parameters:
1660d3b3e55cSJunchao Zhang + rootdata   - root values to be updated, input state is seen by first process to perform an update
1661d3b3e55cSJunchao Zhang - leafupdate - state at each leaf's respective root immediately prior to my atomic update
1662d3b3e55cSJunchao Zhang 
1663d3b3e55cSJunchao Zhang   Level: advanced
1664d3b3e55cSJunchao Zhang 
1665cab54364SBarry Smith   Note:
1666cab54364SBarry Smith   See `PetscSFFetchAndOpBegin()` for more details.
1667d3b3e55cSJunchao Zhang 
166820662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFFetchAndOpBegin()`, `PetscSFComputeDegreeBegin()`, `PetscSFReduceBegin()`, `PetscSFSetGraph()`, `PetscSFFetchAndOpEnd()`
1669d3b3e55cSJunchao Zhang @*/
1670d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFFetchAndOpWithMemTypeBegin(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, void *rootdata, PetscMemType leafmtype, const void *leafdata, PetscMemType leafupdatemtype, void *leafupdate, MPI_Op op)
1671d71ae5a4SJacob Faibussowitsch {
1672d3b3e55cSJunchao Zhang   PetscFunctionBegin;
1673d3b3e55cSJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
16749566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
16759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_FetchAndOpBegin, sf, 0, 0, 0));
167608401ef6SPierre Jolivet   PetscCheck(leafmtype == leafupdatemtype, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for leafdata and leafupdate in different memory types");
1677dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, FetchAndOpBegin, unit, rootmtype, rootdata, leafmtype, leafdata, leafupdate, op);
16789566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_FetchAndOpBegin, sf, 0, 0, 0));
16793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1680d3b3e55cSJunchao Zhang }
1681d3b3e55cSJunchao Zhang 
1682d3b3e55cSJunchao Zhang /*@C
168320662ed9SBarry Smith   PetscSFFetchAndOpEnd - end operation started in matching call to `PetscSFFetchAndOpBegin()` or `PetscSFFetchAndOpWithMemTypeBegin()`
168420662ed9SBarry Smith   to fetch values from roots and update atomically by applying operation using my leaf value
1685a1729e3fSJunchao Zhang 
1686a1729e3fSJunchao Zhang   Collective
1687a1729e3fSJunchao Zhang 
16884165533cSJose E. Roman   Input Parameters:
1689a1729e3fSJunchao Zhang + sf       - star forest
1690a1729e3fSJunchao Zhang . unit     - data type
1691a1729e3fSJunchao Zhang . leafdata - leaf values to use in reduction
1692a1729e3fSJunchao Zhang - op       - operation to use for reduction
1693a1729e3fSJunchao Zhang 
16944165533cSJose E. Roman   Output Parameters:
1695a1729e3fSJunchao Zhang + rootdata   - root values to be updated, input state is seen by first process to perform an update
1696a1729e3fSJunchao Zhang - leafupdate - state at each leaf's respective root immediately prior to my atomic update
1697a1729e3fSJunchao Zhang 
1698a1729e3fSJunchao Zhang   Level: advanced
1699a1729e3fSJunchao Zhang 
170020662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeEnd()`, `PetscSFReduceEnd()`, `PetscSFSetGraph()`, `PetscSFFetchAndOpBegin()`, `PetscSFFetchAndOpWithMemTypeBegin()`
1701a1729e3fSJunchao Zhang @*/
1702d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf, MPI_Datatype unit, void *rootdata, const void *leafdata, void *leafupdate, MPI_Op op)
1703d71ae5a4SJacob Faibussowitsch {
1704a1729e3fSJunchao Zhang   PetscFunctionBegin;
1705a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
17069566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_FetchAndOpEnd, sf, 0, 0, 0));
1707dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, FetchAndOpEnd, unit, rootdata, leafdata, leafupdate, op);
17089566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_FetchAndOpEnd, sf, 0, 0, 0));
17093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1710a1729e3fSJunchao Zhang }
1711a1729e3fSJunchao Zhang 
1712a1729e3fSJunchao Zhang /*@C
1713cab54364SBarry Smith   PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with `PetscSFComputeDegreeEnd()`
171495fce210SBarry Smith 
171595fce210SBarry Smith   Collective
171695fce210SBarry Smith 
17174165533cSJose E. Roman   Input Parameter:
171895fce210SBarry Smith . sf - star forest
171995fce210SBarry Smith 
17204165533cSJose E. Roman   Output Parameter:
172195fce210SBarry Smith . degree - degree of each root vertex
172295fce210SBarry Smith 
172395fce210SBarry Smith   Level: advanced
172495fce210SBarry Smith 
1725cab54364SBarry Smith   Note:
172620662ed9SBarry Smith   The returned array is owned by `PetscSF` and automatically freed by `PetscSFDestroy()`. Hence there is no need to call `PetscFree()` on it.
1727ffe67aa5SVáclav Hapla 
1728cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGatherBegin()`, `PetscSFComputeDegreeEnd()`
172995fce210SBarry Smith @*/
1730d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf, const PetscInt **degree)
1731d71ae5a4SJacob Faibussowitsch {
173295fce210SBarry Smith   PetscFunctionBegin;
173395fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
173495fce210SBarry Smith   PetscSFCheckGraphSet(sf, 1);
17354f572ea9SToby Isaac   PetscAssertPointer(degree, 2);
1736803bd9e8SMatthew G. Knepley   if (!sf->degreeknown) {
17375b0d146aSStefano Zampini     PetscInt i, nroots = sf->nroots, maxlocal;
173828b400f6SJacob Faibussowitsch     PetscCheck(!sf->degree, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
17395b0d146aSStefano Zampini     maxlocal = sf->maxleaf - sf->minleaf + 1;
17409566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nroots, &sf->degree));
17419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PetscMax(maxlocal, 1), &sf->degreetmp)); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
174229046d53SLisandro Dalcin     for (i = 0; i < nroots; i++) sf->degree[i] = 0;
17439837ea96SMatthew G. Knepley     for (i = 0; i < maxlocal; i++) sf->degreetmp[i] = 1;
17449566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf, MPIU_INT, sf->degreetmp - sf->minleaf, sf->degree, MPI_SUM));
174595fce210SBarry Smith   }
174695fce210SBarry Smith   *degree = NULL;
17473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
174895fce210SBarry Smith }
174995fce210SBarry Smith 
175095fce210SBarry Smith /*@C
1751cab54364SBarry Smith   PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with `PetscSFComputeDegreeBegin()`
175295fce210SBarry Smith 
175395fce210SBarry Smith   Collective
175495fce210SBarry Smith 
17554165533cSJose E. Roman   Input Parameter:
175695fce210SBarry Smith . sf - star forest
175795fce210SBarry Smith 
17584165533cSJose E. Roman   Output Parameter:
175995fce210SBarry Smith . degree - degree of each root vertex
176095fce210SBarry Smith 
176195fce210SBarry Smith   Level: developer
176295fce210SBarry Smith 
1763cab54364SBarry Smith   Note:
176420662ed9SBarry Smith   The returned array is owned by `PetscSF` and automatically freed by `PetscSFDestroy()`. Hence there is no need to call `PetscFree()` on it.
1765ffe67aa5SVáclav Hapla 
1766cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGatherBegin()`, `PetscSFComputeDegreeBegin()`
176795fce210SBarry Smith @*/
1768d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf, const PetscInt **degree)
1769d71ae5a4SJacob Faibussowitsch {
177095fce210SBarry Smith   PetscFunctionBegin;
177195fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
177295fce210SBarry Smith   PetscSFCheckGraphSet(sf, 1);
17734f572ea9SToby Isaac   PetscAssertPointer(degree, 2);
177495fce210SBarry Smith   if (!sf->degreeknown) {
177528b400f6SJacob Faibussowitsch     PetscCheck(sf->degreetmp, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
17769566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf, MPIU_INT, sf->degreetmp - sf->minleaf, sf->degree, MPI_SUM));
17779566063dSJacob Faibussowitsch     PetscCall(PetscFree(sf->degreetmp));
177895fce210SBarry Smith     sf->degreeknown = PETSC_TRUE;
177995fce210SBarry Smith   }
178095fce210SBarry Smith   *degree = sf->degree;
17813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
178295fce210SBarry Smith }
178395fce210SBarry Smith 
1784673100f5SVaclav Hapla /*@C
178520662ed9SBarry Smith   PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-`PetscSF` returned by `PetscSFGetMultiSF()`).
178666dfcd1aSVaclav Hapla   Each multi-root is assigned index of the corresponding original root.
1787673100f5SVaclav Hapla 
1788673100f5SVaclav Hapla   Collective
1789673100f5SVaclav Hapla 
17904165533cSJose E. Roman   Input Parameters:
1791673100f5SVaclav Hapla + sf     - star forest
1792cab54364SBarry Smith - degree - degree of each root vertex, computed with `PetscSFComputeDegreeBegin()`/`PetscSFComputeDegreeEnd()`
1793673100f5SVaclav Hapla 
17944165533cSJose E. Roman   Output Parameters:
179520662ed9SBarry Smith + nMultiRoots             - (optional) number of multi-roots (roots of multi-`PetscSF`)
179620662ed9SBarry Smith - multiRootsOrigNumbering - original indices of multi-roots; length of this array is `nMultiRoots`
1797673100f5SVaclav Hapla 
1798673100f5SVaclav Hapla   Level: developer
1799673100f5SVaclav Hapla 
1800cab54364SBarry Smith   Note:
180120662ed9SBarry Smith   The returned array `multiRootsOrigNumbering` is newly allocated and should be destroyed with `PetscFree()` when no longer needed.
1802ffe67aa5SVáclav Hapla 
1803cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeBegin()`, `PetscSFComputeDegreeEnd()`, `PetscSFGetMultiSF()`
1804673100f5SVaclav Hapla @*/
1805d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[])
1806d71ae5a4SJacob Faibussowitsch {
1807673100f5SVaclav Hapla   PetscSF  msf;
1808673100f5SVaclav Hapla   PetscInt i, j, k, nroots, nmroots;
1809673100f5SVaclav Hapla 
1810673100f5SVaclav Hapla   PetscFunctionBegin;
1811673100f5SVaclav Hapla   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
18129566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
18134f572ea9SToby Isaac   if (nroots) PetscAssertPointer(degree, 2);
18144f572ea9SToby Isaac   if (nMultiRoots) PetscAssertPointer(nMultiRoots, 3);
18154f572ea9SToby Isaac   PetscAssertPointer(multiRootsOrigNumbering, 4);
18169566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf, &msf));
18179566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL));
18189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nmroots, multiRootsOrigNumbering));
1819673100f5SVaclav Hapla   for (i = 0, j = 0, k = 0; i < nroots; i++) {
1820673100f5SVaclav Hapla     if (!degree[i]) continue;
1821ad540459SPierre Jolivet     for (j = 0; j < degree[i]; j++, k++) (*multiRootsOrigNumbering)[k] = i;
1822673100f5SVaclav Hapla   }
182308401ef6SPierre Jolivet   PetscCheck(k == nmroots, PETSC_COMM_SELF, PETSC_ERR_PLIB, "sanity check fail");
182466dfcd1aSVaclav Hapla   if (nMultiRoots) *nMultiRoots = nmroots;
18253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1826673100f5SVaclav Hapla }
1827673100f5SVaclav Hapla 
182895fce210SBarry Smith /*@C
1829cab54364SBarry Smith   PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with `PetscSFGatherEnd()`
183095fce210SBarry Smith 
183195fce210SBarry Smith   Collective
183295fce210SBarry Smith 
18334165533cSJose E. Roman   Input Parameters:
183495fce210SBarry Smith + sf       - star forest
183595fce210SBarry Smith . unit     - data type
183695fce210SBarry Smith - leafdata - leaf data to gather to roots
183795fce210SBarry Smith 
18384165533cSJose E. Roman   Output Parameter:
183995fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree
184095fce210SBarry Smith 
184195fce210SBarry Smith   Level: intermediate
184295fce210SBarry Smith 
1843cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeBegin()`, `PetscSFScatterBegin()`
184495fce210SBarry Smith @*/
1845d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGatherBegin(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *multirootdata)
1846d71ae5a4SJacob Faibussowitsch {
1847a5526d50SJunchao Zhang   PetscSF multi = NULL;
184895fce210SBarry Smith 
184995fce210SBarry Smith   PetscFunctionBegin;
185095fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
18519566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
18529566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf, &multi));
18539566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(multi, unit, leafdata, multirootdata, MPI_REPLACE));
18543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
185595fce210SBarry Smith }
185695fce210SBarry Smith 
185795fce210SBarry Smith /*@C
1858cab54364SBarry Smith   PetscSFGatherEnd - ends pointwise gather operation that was started with `PetscSFGatherBegin()`
185995fce210SBarry Smith 
186095fce210SBarry Smith   Collective
186195fce210SBarry Smith 
18624165533cSJose E. Roman   Input Parameters:
186395fce210SBarry Smith + sf       - star forest
186495fce210SBarry Smith . unit     - data type
186595fce210SBarry Smith - leafdata - leaf data to gather to roots
186695fce210SBarry Smith 
18674165533cSJose E. Roman   Output Parameter:
186895fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree
186995fce210SBarry Smith 
187095fce210SBarry Smith   Level: intermediate
187195fce210SBarry Smith 
1872cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeEnd()`, `PetscSFScatterEnd()`
187395fce210SBarry Smith @*/
1874d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGatherEnd(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *multirootdata)
1875d71ae5a4SJacob Faibussowitsch {
1876a5526d50SJunchao Zhang   PetscSF multi = NULL;
187795fce210SBarry Smith 
187895fce210SBarry Smith   PetscFunctionBegin;
187995fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
18809566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf, &multi));
18819566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(multi, unit, leafdata, multirootdata, MPI_REPLACE));
18823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
188395fce210SBarry Smith }
188495fce210SBarry Smith 
188595fce210SBarry Smith /*@C
1886cab54364SBarry Smith   PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with `PetscSFScatterEnd()`
188795fce210SBarry Smith 
188895fce210SBarry Smith   Collective
188995fce210SBarry Smith 
18904165533cSJose E. Roman   Input Parameters:
189195fce210SBarry Smith + sf            - star forest
189295fce210SBarry Smith . unit          - data type
189395fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf
189495fce210SBarry Smith 
18954165533cSJose E. Roman   Output Parameter:
189695fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root
189795fce210SBarry Smith 
189895fce210SBarry Smith   Level: intermediate
189995fce210SBarry Smith 
190020662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeBegin()`, `PetscSFScatterEnd()`
190195fce210SBarry Smith @*/
1902d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFScatterBegin(PetscSF sf, MPI_Datatype unit, const void *multirootdata, void *leafdata)
1903d71ae5a4SJacob Faibussowitsch {
1904a5526d50SJunchao Zhang   PetscSF multi = NULL;
190595fce210SBarry Smith 
190695fce210SBarry Smith   PetscFunctionBegin;
190795fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
19089566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
19099566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf, &multi));
19109566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(multi, unit, multirootdata, leafdata, MPI_REPLACE));
19113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
191295fce210SBarry Smith }
191395fce210SBarry Smith 
191495fce210SBarry Smith /*@C
1915cab54364SBarry Smith   PetscSFScatterEnd - ends pointwise scatter operation that was started with `PetscSFScatterBegin()`
191695fce210SBarry Smith 
191795fce210SBarry Smith   Collective
191895fce210SBarry Smith 
19194165533cSJose E. Roman   Input Parameters:
192095fce210SBarry Smith + sf            - star forest
192195fce210SBarry Smith . unit          - data type
192295fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf
192395fce210SBarry Smith 
19244165533cSJose E. Roman   Output Parameter:
192595fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root
192695fce210SBarry Smith 
192795fce210SBarry Smith   Level: intermediate
192895fce210SBarry Smith 
192920662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeEnd()`, `PetscSFScatterBegin()`
193095fce210SBarry Smith @*/
1931d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFScatterEnd(PetscSF sf, MPI_Datatype unit, const void *multirootdata, void *leafdata)
1932d71ae5a4SJacob Faibussowitsch {
1933a5526d50SJunchao Zhang   PetscSF multi = NULL;
193495fce210SBarry Smith 
193595fce210SBarry Smith   PetscFunctionBegin;
193695fce210SBarry Smith   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
19379566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf, &multi));
19389566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(multi, unit, multirootdata, leafdata, MPI_REPLACE));
19393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
194095fce210SBarry Smith }
1941a7b3aa13SAta Mesgarnejad 
1942d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFCheckLeavesUnique_Private(PetscSF sf)
1943d71ae5a4SJacob Faibussowitsch {
1944a072220fSLawrence Mitchell   PetscInt        i, n, nleaves;
1945a072220fSLawrence Mitchell   const PetscInt *ilocal = NULL;
1946a072220fSLawrence Mitchell   PetscHSetI      seen;
1947a072220fSLawrence Mitchell 
1948a072220fSLawrence Mitchell   PetscFunctionBegin;
1949b458e8f1SJose E. Roman   if (PetscDefined(USE_DEBUG)) {
19509566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, &ilocal, NULL));
19519566063dSJacob Faibussowitsch     PetscCall(PetscHSetICreate(&seen));
1952a072220fSLawrence Mitchell     for (i = 0; i < nleaves; i++) {
1953a072220fSLawrence Mitchell       const PetscInt leaf = ilocal ? ilocal[i] : i;
19549566063dSJacob Faibussowitsch       PetscCall(PetscHSetIAdd(seen, leaf));
1955a072220fSLawrence Mitchell     }
19569566063dSJacob Faibussowitsch     PetscCall(PetscHSetIGetSize(seen, &n));
195708401ef6SPierre Jolivet     PetscCheck(n == nleaves, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Provided leaves have repeated values: all leaves must be unique");
19589566063dSJacob Faibussowitsch     PetscCall(PetscHSetIDestroy(&seen));
1959b458e8f1SJose E. Roman   }
19603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1961a072220fSLawrence Mitchell }
196254729392SStefano Zampini 
1963a7b3aa13SAta Mesgarnejad /*@
1964cab54364SBarry Smith   PetscSFCompose - Compose a new `PetscSF` by putting the second `PetscSF` under the first one in a top (roots) down (leaves) view
1965a7b3aa13SAta Mesgarnejad 
1966a7b3aa13SAta Mesgarnejad   Input Parameters:
1967cab54364SBarry Smith + sfA - The first `PetscSF`
1968cab54364SBarry Smith - sfB - The second `PetscSF`
1969a7b3aa13SAta Mesgarnejad 
19702fe279fdSBarry Smith   Output Parameter:
1971cab54364SBarry Smith . sfBA - The composite `PetscSF`
1972a7b3aa13SAta Mesgarnejad 
1973a7b3aa13SAta Mesgarnejad   Level: developer
1974a7b3aa13SAta Mesgarnejad 
1975a072220fSLawrence Mitchell   Notes:
1976cab54364SBarry Smith   Currently, the two `PetscSF`s must be defined on congruent communicators and they must be true star
197754729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots.
197854729392SStefano Zampini 
197920662ed9SBarry Smith   `sfA`'s leaf space and `sfB`'s root space might be partially overlapped. The composition builds
198020662ed9SBarry Smith   a graph with `sfA`'s roots and `sfB`'s leaves only when there is a path between them. Unconnected
198120662ed9SBarry Smith   nodes (roots or leaves) are not in `sfBA`. Doing a `PetscSFBcastBegin()`/`PetscSFBcastEnd()` on the new `PetscSF` is equivalent to doing a
198220662ed9SBarry Smith   `PetscSFBcastBegin()`/`PetscSFBcastEnd()` on `sfA`, then a `PetscSFBcastBegin()`/`PetscSFBcastEnd()` on `sfB`, on connected nodes.
1983a072220fSLawrence Mitchell 
1984db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFComposeInverse()`, `PetscSFGetGraph()`, `PetscSFSetGraph()`
1985a7b3aa13SAta Mesgarnejad @*/
1986d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA)
1987d71ae5a4SJacob Faibussowitsch {
1988a7b3aa13SAta Mesgarnejad   const PetscSFNode *remotePointsA, *remotePointsB;
1989d41018fbSJunchao Zhang   PetscSFNode       *remotePointsBA = NULL, *reorderedRemotePointsA = NULL, *leafdataB;
199054729392SStefano Zampini   const PetscInt    *localPointsA, *localPointsB;
199154729392SStefano Zampini   PetscInt          *localPointsBA;
199254729392SStefano Zampini   PetscInt           i, numRootsA, numLeavesA, numRootsB, numLeavesB, minleaf, maxleaf, numLeavesBA;
199354729392SStefano Zampini   PetscBool          denseB;
1994a7b3aa13SAta Mesgarnejad 
1995a7b3aa13SAta Mesgarnejad   PetscFunctionBegin;
1996a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1);
199729046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfA, 1);
199829046d53SLisandro Dalcin   PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 2);
199929046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfB, 2);
200054729392SStefano Zampini   PetscCheckSameComm(sfA, 1, sfB, 2);
20014f572ea9SToby Isaac   PetscAssertPointer(sfBA, 3);
20029566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckLeavesUnique_Private(sfA));
20039566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckLeavesUnique_Private(sfB));
200454729392SStefano Zampini 
20059566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA));
20069566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB));
200720662ed9SBarry Smith   /* Make sure that PetscSFBcast{Begin, End}(sfB, ...) works with root data of size
200820662ed9SBarry Smith      numRootsB; otherwise, garbage will be broadcasted.
200920662ed9SBarry Smith      Example (comm size = 1):
201020662ed9SBarry Smith      sfA: 0 <- (0, 0)
201120662ed9SBarry Smith      sfB: 100 <- (0, 0)
201220662ed9SBarry Smith           101 <- (0, 1)
201320662ed9SBarry Smith      Here, we have remotePointsA = [(0, 0)], but for remotePointsA to be a valid tartget
201420662ed9SBarry Smith      of sfB, it has to be recasted as [(0, 0), (-1, -1)] so that points 100 and 101 would
201520662ed9SBarry Smith      receive (0, 0) and (-1, -1), respectively, when PetscSFBcast(sfB, ...) is called on
201620662ed9SBarry Smith      remotePointsA; if not recasted, point 101 would receive a garbage value.             */
20179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numRootsB, &reorderedRemotePointsA));
201854729392SStefano Zampini   for (i = 0; i < numRootsB; i++) {
201954729392SStefano Zampini     reorderedRemotePointsA[i].rank  = -1;
202054729392SStefano Zampini     reorderedRemotePointsA[i].index = -1;
202154729392SStefano Zampini   }
202254729392SStefano Zampini   for (i = 0; i < numLeavesA; i++) {
20230ea77edaSksagiyam     PetscInt localp = localPointsA ? localPointsA[i] : i;
20240ea77edaSksagiyam 
20250ea77edaSksagiyam     if (localp >= numRootsB) continue;
20260ea77edaSksagiyam     reorderedRemotePointsA[localp] = remotePointsA[i];
202754729392SStefano Zampini   }
2028d41018fbSJunchao Zhang   remotePointsA = reorderedRemotePointsA;
20299566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(sfB, &minleaf, &maxleaf));
20309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxleaf - minleaf + 1, &leafdataB));
20310ea77edaSksagiyam   for (i = 0; i < maxleaf - minleaf + 1; i++) {
20320ea77edaSksagiyam     leafdataB[i].rank  = -1;
20330ea77edaSksagiyam     leafdataB[i].index = -1;
20340ea77edaSksagiyam   }
20358e3a54c0SPierre Jolivet   PetscCall(PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, PetscSafePointerPlusOffset(leafdataB, -minleaf), MPI_REPLACE));
20368e3a54c0SPierre Jolivet   PetscCall(PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, PetscSafePointerPlusOffset(leafdataB, -minleaf), MPI_REPLACE));
20379566063dSJacob Faibussowitsch   PetscCall(PetscFree(reorderedRemotePointsA));
2038d41018fbSJunchao Zhang 
203954729392SStefano Zampini   denseB = (PetscBool)!localPointsB;
204054729392SStefano Zampini   for (i = 0, numLeavesBA = 0; i < numLeavesB; i++) {
204154729392SStefano Zampini     if (leafdataB[localPointsB ? localPointsB[i] - minleaf : i].rank == -1) denseB = PETSC_FALSE;
204254729392SStefano Zampini     else numLeavesBA++;
204354729392SStefano Zampini   }
204454729392SStefano Zampini   if (denseB) {
2045d41018fbSJunchao Zhang     localPointsBA  = NULL;
2046d41018fbSJunchao Zhang     remotePointsBA = leafdataB;
2047d41018fbSJunchao Zhang   } else {
20489566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numLeavesBA, &localPointsBA));
20499566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numLeavesBA, &remotePointsBA));
205054729392SStefano Zampini     for (i = 0, numLeavesBA = 0; i < numLeavesB; i++) {
205154729392SStefano Zampini       const PetscInt l = localPointsB ? localPointsB[i] : i;
205254729392SStefano Zampini 
205354729392SStefano Zampini       if (leafdataB[l - minleaf].rank == -1) continue;
205454729392SStefano Zampini       remotePointsBA[numLeavesBA] = leafdataB[l - minleaf];
205554729392SStefano Zampini       localPointsBA[numLeavesBA]  = l;
205654729392SStefano Zampini       numLeavesBA++;
205754729392SStefano Zampini     }
20589566063dSJacob Faibussowitsch     PetscCall(PetscFree(leafdataB));
2059d41018fbSJunchao Zhang   }
20609566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfA), sfBA));
20619566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(*sfBA));
20629566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sfBA, numRootsA, numLeavesBA, localPointsBA, PETSC_OWN_POINTER, remotePointsBA, PETSC_OWN_POINTER));
20633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2064a7b3aa13SAta Mesgarnejad }
20651c6ba672SJunchao Zhang 
206604c0ada0SJunchao Zhang /*@
2067cab54364SBarry Smith   PetscSFComposeInverse - Compose a new `PetscSF` by putting the inverse of the second `PetscSF` under the first one
206804c0ada0SJunchao Zhang 
206904c0ada0SJunchao Zhang   Input Parameters:
2070cab54364SBarry Smith + sfA - The first `PetscSF`
2071cab54364SBarry Smith - sfB - The second `PetscSF`
207204c0ada0SJunchao Zhang 
20732fe279fdSBarry Smith   Output Parameter:
2074cab54364SBarry Smith . sfBA - The composite `PetscSF`.
207504c0ada0SJunchao Zhang 
207604c0ada0SJunchao Zhang   Level: developer
207704c0ada0SJunchao Zhang 
207854729392SStefano Zampini   Notes:
207920662ed9SBarry Smith   Currently, the two `PetscSF`s must be defined on congruent communicators and they must be true star
208054729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots. Even more, all roots of the
208120662ed9SBarry Smith   second `PetscSF` must have a degree of 1, i.e., no roots have more than one leaf connected.
208254729392SStefano Zampini 
208320662ed9SBarry Smith   `sfA`'s leaf space and `sfB`'s leaf space might be partially overlapped. The composition builds
208420662ed9SBarry Smith   a graph with `sfA`'s roots and `sfB`'s roots only when there is a path between them. Unconnected
208520662ed9SBarry Smith   roots are not in `sfBA`. Doing a `PetscSFBcastBegin()`/`PetscSFBcastEnd()` on the new `PetscSF` is equivalent to doing a `PetscSFBcastBegin()`/`PetscSFBcastEnd()`
208620662ed9SBarry Smith   on `sfA`, then
208720662ed9SBarry Smith   a `PetscSFReduceBegin()`/`PetscSFReduceEnd()` on `sfB`, on connected roots.
208854729392SStefano Zampini 
2089db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFCompose()`, `PetscSFGetGraph()`, `PetscSFSetGraph()`, `PetscSFCreateInverseSF()`
209004c0ada0SJunchao Zhang @*/
2091d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFComposeInverse(PetscSF sfA, PetscSF sfB, PetscSF *sfBA)
2092d71ae5a4SJacob Faibussowitsch {
209304c0ada0SJunchao Zhang   const PetscSFNode *remotePointsA, *remotePointsB;
209404c0ada0SJunchao Zhang   PetscSFNode       *remotePointsBA;
209504c0ada0SJunchao Zhang   const PetscInt    *localPointsA, *localPointsB;
209654729392SStefano Zampini   PetscSFNode       *reorderedRemotePointsA = NULL;
209754729392SStefano Zampini   PetscInt           i, numRootsA, numLeavesA, numLeavesBA, numRootsB, numLeavesB, minleaf, maxleaf, *localPointsBA;
20985b0d146aSStefano Zampini   MPI_Op             op;
20995b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
21005b0d146aSStefano Zampini   PetscBool iswin;
21015b0d146aSStefano Zampini #endif
210204c0ada0SJunchao Zhang 
210304c0ada0SJunchao Zhang   PetscFunctionBegin;
210404c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1);
210504c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfA, 1);
210604c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 2);
210704c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfB, 2);
210854729392SStefano Zampini   PetscCheckSameComm(sfA, 1, sfB, 2);
21094f572ea9SToby Isaac   PetscAssertPointer(sfBA, 3);
21109566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckLeavesUnique_Private(sfA));
21119566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckLeavesUnique_Private(sfB));
211254729392SStefano Zampini 
21139566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA));
21149566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB));
21155b0d146aSStefano Zampini 
21165b0d146aSStefano Zampini   /* TODO: Check roots of sfB have degree of 1 */
21175b0d146aSStefano Zampini   /* Once we implement it, we can replace the MPI_MAXLOC
211883df288dSJunchao Zhang      with MPI_REPLACE. In that case, MPI_MAXLOC and MPI_REPLACE have the same effect.
21195b0d146aSStefano Zampini      We use MPI_MAXLOC only to have a deterministic output from this routine if
21205b0d146aSStefano Zampini      the root condition is not meet.
21215b0d146aSStefano Zampini    */
21225b0d146aSStefano Zampini   op = MPI_MAXLOC;
21235b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
21245b0d146aSStefano Zampini   /* we accept a non-deterministic output (if any) with PETSCSFWINDOW, since MPI_MAXLOC cannot operate on MPIU_2INT with MPI_Accumulate */
21259566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)sfB, PETSCSFWINDOW, &iswin));
212683df288dSJunchao Zhang   if (iswin) op = MPI_REPLACE;
21275b0d146aSStefano Zampini #endif
21285b0d146aSStefano Zampini 
21299566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(sfB, &minleaf, &maxleaf));
21309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxleaf - minleaf + 1, &reorderedRemotePointsA));
213154729392SStefano Zampini   for (i = 0; i < maxleaf - minleaf + 1; i++) {
213254729392SStefano Zampini     reorderedRemotePointsA[i].rank  = -1;
213354729392SStefano Zampini     reorderedRemotePointsA[i].index = -1;
213454729392SStefano Zampini   }
213554729392SStefano Zampini   if (localPointsA) {
213654729392SStefano Zampini     for (i = 0; i < numLeavesA; i++) {
213754729392SStefano Zampini       if (localPointsA[i] > maxleaf || localPointsA[i] < minleaf) continue;
213854729392SStefano Zampini       reorderedRemotePointsA[localPointsA[i] - minleaf] = remotePointsA[i];
213954729392SStefano Zampini     }
214054729392SStefano Zampini   } else {
214154729392SStefano Zampini     for (i = 0; i < numLeavesA; i++) {
214254729392SStefano Zampini       if (i > maxleaf || i < minleaf) continue;
214354729392SStefano Zampini       reorderedRemotePointsA[i - minleaf] = remotePointsA[i];
214454729392SStefano Zampini     }
214554729392SStefano Zampini   }
214654729392SStefano Zampini 
21479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numRootsB, &localPointsBA));
21489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numRootsB, &remotePointsBA));
214954729392SStefano Zampini   for (i = 0; i < numRootsB; i++) {
215054729392SStefano Zampini     remotePointsBA[i].rank  = -1;
215154729392SStefano Zampini     remotePointsBA[i].index = -1;
215254729392SStefano Zampini   }
215354729392SStefano Zampini 
21548e3a54c0SPierre Jolivet   PetscCall(PetscSFReduceBegin(sfB, MPIU_2INT, PetscSafePointerPlusOffset(reorderedRemotePointsA, -minleaf), remotePointsBA, op));
21558e3a54c0SPierre Jolivet   PetscCall(PetscSFReduceEnd(sfB, MPIU_2INT, PetscSafePointerPlusOffset(reorderedRemotePointsA, -minleaf), remotePointsBA, op));
21569566063dSJacob Faibussowitsch   PetscCall(PetscFree(reorderedRemotePointsA));
215754729392SStefano Zampini   for (i = 0, numLeavesBA = 0; i < numRootsB; i++) {
215854729392SStefano Zampini     if (remotePointsBA[i].rank == -1) continue;
215954729392SStefano Zampini     remotePointsBA[numLeavesBA].rank  = remotePointsBA[i].rank;
216054729392SStefano Zampini     remotePointsBA[numLeavesBA].index = remotePointsBA[i].index;
216154729392SStefano Zampini     localPointsBA[numLeavesBA]        = i;
216254729392SStefano Zampini     numLeavesBA++;
216354729392SStefano Zampini   }
21649566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfA), sfBA));
21659566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(*sfBA));
21669566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sfBA, numRootsA, numLeavesBA, localPointsBA, PETSC_OWN_POINTER, remotePointsBA, PETSC_OWN_POINTER));
21673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
216804c0ada0SJunchao Zhang }
216904c0ada0SJunchao Zhang 
21701c6ba672SJunchao Zhang /*
2171cab54364SBarry Smith   PetscSFCreateLocalSF_Private - Creates a local `PetscSF` that only has intra-process edges of the global `PetscSF`
21721c6ba672SJunchao Zhang 
21732fe279fdSBarry Smith   Input Parameter:
2174cab54364SBarry Smith . sf - The global `PetscSF`
21751c6ba672SJunchao Zhang 
21762fe279fdSBarry Smith   Output Parameter:
2177cab54364SBarry Smith . out - The local `PetscSF`
2178cab54364SBarry Smith 
2179cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`
21801c6ba672SJunchao Zhang  */
2181d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateLocalSF_Private(PetscSF sf, PetscSF *out)
2182d71ae5a4SJacob Faibussowitsch {
21831c6ba672SJunchao Zhang   MPI_Comm           comm;
21841c6ba672SJunchao Zhang   PetscMPIInt        myrank;
21851c6ba672SJunchao Zhang   const PetscInt    *ilocal;
21861c6ba672SJunchao Zhang   const PetscSFNode *iremote;
21871c6ba672SJunchao Zhang   PetscInt           i, j, nroots, nleaves, lnleaves, *lilocal;
21881c6ba672SJunchao Zhang   PetscSFNode       *liremote;
21891c6ba672SJunchao Zhang   PetscSF            lsf;
21901c6ba672SJunchao Zhang 
21911c6ba672SJunchao Zhang   PetscFunctionBegin;
21921c6ba672SJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
2193dbbe0bcdSBarry Smith   if (sf->ops->CreateLocalSF) PetscUseTypeMethod(sf, CreateLocalSF, out);
2194dbbe0bcdSBarry Smith   else {
21951c6ba672SJunchao Zhang     /* Could use PetscSFCreateEmbeddedLeafSF, but since we know the comm is PETSC_COMM_SELF, we can make it fast */
21969566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
21979566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm, &myrank));
21981c6ba672SJunchao Zhang 
21991c6ba672SJunchao Zhang     /* Find out local edges and build a local SF */
22009566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, &iremote));
22019371c9d4SSatish Balay     for (i = lnleaves = 0; i < nleaves; i++) {
22029371c9d4SSatish Balay       if (iremote[i].rank == (PetscInt)myrank) lnleaves++;
22039371c9d4SSatish Balay     }
22049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lnleaves, &lilocal));
22059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lnleaves, &liremote));
22061c6ba672SJunchao Zhang 
22071c6ba672SJunchao Zhang     for (i = j = 0; i < nleaves; i++) {
22081c6ba672SJunchao Zhang       if (iremote[i].rank == (PetscInt)myrank) {
22091c6ba672SJunchao Zhang         lilocal[j]        = ilocal ? ilocal[i] : i; /* ilocal=NULL for contiguous storage */
22101c6ba672SJunchao Zhang         liremote[j].rank  = 0;                      /* rank in PETSC_COMM_SELF */
22111c6ba672SJunchao Zhang         liremote[j].index = iremote[i].index;
22121c6ba672SJunchao Zhang         j++;
22131c6ba672SJunchao Zhang       }
22141c6ba672SJunchao Zhang     }
22159566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PETSC_COMM_SELF, &lsf));
22169566063dSJacob Faibussowitsch     PetscCall(PetscSFSetFromOptions(lsf));
22179566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(lsf, nroots, lnleaves, lilocal, PETSC_OWN_POINTER, liremote, PETSC_OWN_POINTER));
22189566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(lsf));
22191c6ba672SJunchao Zhang     *out = lsf;
22201c6ba672SJunchao Zhang   }
22213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22221c6ba672SJunchao Zhang }
2223dd5b3ca6SJunchao Zhang 
2224dd5b3ca6SJunchao Zhang /* Similar to PetscSFBcast, but only Bcast to leaves on rank 0 */
2225d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFBcastToZero_Private(PetscSF sf, MPI_Datatype unit, const void *rootdata, void *leafdata)
2226d71ae5a4SJacob Faibussowitsch {
2227eb02082bSJunchao Zhang   PetscMemType rootmtype, leafmtype;
2228dd5b3ca6SJunchao Zhang 
2229dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
2230dd5b3ca6SJunchao Zhang   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
22319566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
22329566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_BcastBegin, sf, 0, 0, 0));
22339566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(rootdata, &rootmtype));
22349566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafdata, &leafmtype));
2235dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf, BcastToZero, unit, rootmtype, rootdata, leafmtype, leafdata);
22369566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_BcastBegin, sf, 0, 0, 0));
22373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2238dd5b3ca6SJunchao Zhang }
2239dd5b3ca6SJunchao Zhang 
2240157edd7aSVaclav Hapla /*@
2241cab54364SBarry Smith   PetscSFConcatenate - concatenate multiple `PetscSF` into one
2242157edd7aSVaclav Hapla 
2243157edd7aSVaclav Hapla   Input Parameters:
2244157edd7aSVaclav Hapla + comm        - the communicator
2245cab54364SBarry Smith . nsfs        - the number of input `PetscSF`
2246cab54364SBarry Smith . sfs         - the array of input `PetscSF`
22471f40158dSVaclav Hapla . rootMode    - the root mode specifying how roots are handled
224820662ed9SBarry Smith - leafOffsets - the array of local leaf offsets, one for each input `PetscSF`, or `NULL` for contiguous storage
2249157edd7aSVaclav Hapla 
22502fe279fdSBarry Smith   Output Parameter:
2251cab54364SBarry Smith . newsf - The resulting `PetscSF`
2252157edd7aSVaclav Hapla 
22531f40158dSVaclav Hapla   Level: advanced
2254157edd7aSVaclav Hapla 
2255157edd7aSVaclav Hapla   Notes:
225620662ed9SBarry Smith   The communicator of all `PetscSF`s in `sfs` must be comm.
2257157edd7aSVaclav Hapla 
225820662ed9SBarry Smith   Leaves are always concatenated locally, keeping them ordered by the input `PetscSF` index and original local order.
225920662ed9SBarry Smith 
226020662ed9SBarry Smith   The offsets in `leafOffsets` are added to the original leaf indices.
226120662ed9SBarry Smith 
226220662ed9SBarry Smith   If all input SFs use contiguous leaf storage (`ilocal` = `NULL`), `leafOffsets` can be passed as `NULL` as well.
226320662ed9SBarry Smith   In this case, `NULL` is also passed as `ilocal` to the resulting `PetscSF`.
226420662ed9SBarry Smith 
226520662ed9SBarry Smith   If any input `PetscSF` has non-null `ilocal`, `leafOffsets` is needed to distinguish leaves from different input `PetscSF`s.
2266157edd7aSVaclav Hapla   In this case, user is responsible to provide correct offsets so that the resulting leaves are unique (otherwise an error occurs).
2267157edd7aSVaclav Hapla 
226820662ed9SBarry Smith   All root modes retain the essential connectivity condition.
226920662ed9SBarry Smith   If two leaves of the same input `PetscSF` are connected (sharing the same root), they are also connected in the output `PetscSF`.
227020662ed9SBarry Smith   Parameter `rootMode` controls how the input root spaces are combined.
227120662ed9SBarry Smith   For `PETSCSF_CONCATENATE_ROOTMODE_SHARED`, the root space is considered the same for each input `PetscSF` (checked in debug mode)
227220662ed9SBarry Smith   and is also the same in the output `PetscSF`.
22731f40158dSVaclav Hapla   For `PETSCSF_CONCATENATE_ROOTMODE_LOCAL` and `PETSCSF_CONCATENATE_ROOTMODE_GLOBAL`, the input root spaces are taken as separate and joined.
22741f40158dSVaclav Hapla   `PETSCSF_CONCATENATE_ROOTMODE_LOCAL` joins the root spaces locally;
227520662ed9SBarry Smith   roots of sfs[0], sfs[1], sfs[2], ... are joined on each rank separately, ordered by input `PetscSF` and original local index, and renumbered contiguously.
22761f40158dSVaclav Hapla   `PETSCSF_CONCATENATE_ROOTMODE_GLOBAL` joins the root spaces globally;
22771593df67SStefano Zampini   roots of sfs[0], sfs[1], sfs[2], ... are joined globally, ordered by input `PetscSF` index and original global index, and renumbered contiguously;
22781f40158dSVaclav Hapla   the original root ranks are ignored.
22791f40158dSVaclav Hapla   For both `PETSCSF_CONCATENATE_ROOTMODE_LOCAL` and `PETSCSF_CONCATENATE_ROOTMODE_GLOBAL`,
228020662ed9SBarry Smith   the output `PetscSF`'s root layout is such that the local number of roots is a sum of the input `PetscSF`'s local numbers of roots on each rank
228120662ed9SBarry Smith   to keep the load balancing.
228220662ed9SBarry Smith   However, for `PETSCSF_CONCATENATE_ROOTMODE_GLOBAL`, roots can move to different ranks.
22831f40158dSVaclav Hapla 
22841f40158dSVaclav Hapla   Example:
22851f40158dSVaclav Hapla   We can use src/vec/is/sf/tests/ex18.c to compare the root modes. By running
228620662ed9SBarry Smith .vb
228720662ed9SBarry Smith   make -C $PETSC_DIR/src/vec/is/sf/tests ex18
228820662ed9SBarry Smith   for m in {local,global,shared}; do
228920662ed9SBarry Smith     mpirun -n 2 $PETSC_DIR/src/vec/is/sf/tests/ex18 -nsfs 2 -n 2 -root_mode $m -sf_view
229020662ed9SBarry Smith   done
229120662ed9SBarry Smith .ve
229220662ed9SBarry Smith   we generate two identical `PetscSF`s sf_0 and sf_1,
229320662ed9SBarry Smith .vb
229420662ed9SBarry Smith   PetscSF Object: sf_0 2 MPI processes
229520662ed9SBarry Smith     type: basic
229620662ed9SBarry Smith     rank #leaves #roots
229720662ed9SBarry Smith     [ 0]       4      2
229820662ed9SBarry Smith     [ 1]       4      2
229920662ed9SBarry Smith     leaves      roots       roots in global numbering
230020662ed9SBarry Smith     ( 0,  0) <- ( 0,  0)  =   0
230120662ed9SBarry Smith     ( 0,  1) <- ( 0,  1)  =   1
230220662ed9SBarry Smith     ( 0,  2) <- ( 1,  0)  =   2
230320662ed9SBarry Smith     ( 0,  3) <- ( 1,  1)  =   3
230420662ed9SBarry Smith     ( 1,  0) <- ( 0,  0)  =   0
230520662ed9SBarry Smith     ( 1,  1) <- ( 0,  1)  =   1
230620662ed9SBarry Smith     ( 1,  2) <- ( 1,  0)  =   2
230720662ed9SBarry Smith     ( 1,  3) <- ( 1,  1)  =   3
230820662ed9SBarry Smith .ve
2309e33f79d8SJacob Faibussowitsch   and pass them to `PetscSFConcatenate()` along with different choices of `rootMode`, yielding different result_sf\:
231020662ed9SBarry Smith .vb
231120662ed9SBarry Smith   rootMode = local:
231220662ed9SBarry Smith   PetscSF Object: result_sf 2 MPI processes
231320662ed9SBarry Smith     type: basic
231420662ed9SBarry Smith     rank #leaves #roots
231520662ed9SBarry Smith     [ 0]       8      4
231620662ed9SBarry Smith     [ 1]       8      4
231720662ed9SBarry Smith     leaves      roots       roots in global numbering
231820662ed9SBarry Smith     ( 0,  0) <- ( 0,  0)  =   0
231920662ed9SBarry Smith     ( 0,  1) <- ( 0,  1)  =   1
232020662ed9SBarry Smith     ( 0,  2) <- ( 1,  0)  =   4
232120662ed9SBarry Smith     ( 0,  3) <- ( 1,  1)  =   5
232220662ed9SBarry Smith     ( 0,  4) <- ( 0,  2)  =   2
232320662ed9SBarry Smith     ( 0,  5) <- ( 0,  3)  =   3
232420662ed9SBarry Smith     ( 0,  6) <- ( 1,  2)  =   6
232520662ed9SBarry Smith     ( 0,  7) <- ( 1,  3)  =   7
232620662ed9SBarry Smith     ( 1,  0) <- ( 0,  0)  =   0
232720662ed9SBarry Smith     ( 1,  1) <- ( 0,  1)  =   1
232820662ed9SBarry Smith     ( 1,  2) <- ( 1,  0)  =   4
232920662ed9SBarry Smith     ( 1,  3) <- ( 1,  1)  =   5
233020662ed9SBarry Smith     ( 1,  4) <- ( 0,  2)  =   2
233120662ed9SBarry Smith     ( 1,  5) <- ( 0,  3)  =   3
233220662ed9SBarry Smith     ( 1,  6) <- ( 1,  2)  =   6
233320662ed9SBarry Smith     ( 1,  7) <- ( 1,  3)  =   7
233420662ed9SBarry Smith 
233520662ed9SBarry Smith   rootMode = global:
233620662ed9SBarry Smith   PetscSF Object: result_sf 2 MPI processes
233720662ed9SBarry Smith     type: basic
233820662ed9SBarry Smith     rank #leaves #roots
233920662ed9SBarry Smith     [ 0]       8      4
234020662ed9SBarry Smith     [ 1]       8      4
234120662ed9SBarry Smith     leaves      roots       roots in global numbering
234220662ed9SBarry Smith     ( 0,  0) <- ( 0,  0)  =   0
234320662ed9SBarry Smith     ( 0,  1) <- ( 0,  1)  =   1
234420662ed9SBarry Smith     ( 0,  2) <- ( 0,  2)  =   2
234520662ed9SBarry Smith     ( 0,  3) <- ( 0,  3)  =   3
234620662ed9SBarry Smith     ( 0,  4) <- ( 1,  0)  =   4
234720662ed9SBarry Smith     ( 0,  5) <- ( 1,  1)  =   5
234820662ed9SBarry Smith     ( 0,  6) <- ( 1,  2)  =   6
234920662ed9SBarry Smith     ( 0,  7) <- ( 1,  3)  =   7
235020662ed9SBarry Smith     ( 1,  0) <- ( 0,  0)  =   0
235120662ed9SBarry Smith     ( 1,  1) <- ( 0,  1)  =   1
235220662ed9SBarry Smith     ( 1,  2) <- ( 0,  2)  =   2
235320662ed9SBarry Smith     ( 1,  3) <- ( 0,  3)  =   3
235420662ed9SBarry Smith     ( 1,  4) <- ( 1,  0)  =   4
235520662ed9SBarry Smith     ( 1,  5) <- ( 1,  1)  =   5
235620662ed9SBarry Smith     ( 1,  6) <- ( 1,  2)  =   6
235720662ed9SBarry Smith     ( 1,  7) <- ( 1,  3)  =   7
235820662ed9SBarry Smith 
235920662ed9SBarry Smith   rootMode = shared:
236020662ed9SBarry Smith   PetscSF Object: result_sf 2 MPI processes
236120662ed9SBarry Smith     type: basic
236220662ed9SBarry Smith     rank #leaves #roots
236320662ed9SBarry Smith     [ 0]       8      2
236420662ed9SBarry Smith     [ 1]       8      2
236520662ed9SBarry Smith     leaves      roots       roots in global numbering
236620662ed9SBarry Smith     ( 0,  0) <- ( 0,  0)  =   0
236720662ed9SBarry Smith     ( 0,  1) <- ( 0,  1)  =   1
236820662ed9SBarry Smith     ( 0,  2) <- ( 1,  0)  =   2
236920662ed9SBarry Smith     ( 0,  3) <- ( 1,  1)  =   3
237020662ed9SBarry Smith     ( 0,  4) <- ( 0,  0)  =   0
237120662ed9SBarry Smith     ( 0,  5) <- ( 0,  1)  =   1
237220662ed9SBarry Smith     ( 0,  6) <- ( 1,  0)  =   2
237320662ed9SBarry Smith     ( 0,  7) <- ( 1,  1)  =   3
237420662ed9SBarry Smith     ( 1,  0) <- ( 0,  0)  =   0
237520662ed9SBarry Smith     ( 1,  1) <- ( 0,  1)  =   1
237620662ed9SBarry Smith     ( 1,  2) <- ( 1,  0)  =   2
237720662ed9SBarry Smith     ( 1,  3) <- ( 1,  1)  =   3
237820662ed9SBarry Smith     ( 1,  4) <- ( 0,  0)  =   0
237920662ed9SBarry Smith     ( 1,  5) <- ( 0,  1)  =   1
238020662ed9SBarry Smith     ( 1,  6) <- ( 1,  0)  =   2
238120662ed9SBarry Smith     ( 1,  7) <- ( 1,  1)  =   3
238220662ed9SBarry Smith .ve
23831f40158dSVaclav Hapla 
23841f40158dSVaclav Hapla .seealso: `PetscSF`, `PetscSFCompose()`, `PetscSFGetGraph()`, `PetscSFSetGraph()`, `PetscSFConcatenateRootMode`
2385157edd7aSVaclav Hapla @*/
23861f40158dSVaclav Hapla PetscErrorCode PetscSFConcatenate(MPI_Comm comm, PetscInt nsfs, PetscSF sfs[], PetscSFConcatenateRootMode rootMode, PetscInt leafOffsets[], PetscSF *newsf)
2387d71ae5a4SJacob Faibussowitsch {
2388157edd7aSVaclav Hapla   PetscInt     i, s, nLeaves, nRoots;
2389157edd7aSVaclav Hapla   PetscInt    *leafArrayOffsets;
2390157edd7aSVaclav Hapla   PetscInt    *ilocal_new;
2391157edd7aSVaclav Hapla   PetscSFNode *iremote_new;
2392157edd7aSVaclav Hapla   PetscBool    all_ilocal_null = PETSC_FALSE;
23931f40158dSVaclav Hapla   PetscLayout  glayout         = NULL;
23941f40158dSVaclav Hapla   PetscInt    *gremote         = NULL;
23951f40158dSVaclav Hapla   PetscMPIInt  rank, size;
2396157edd7aSVaclav Hapla 
2397157edd7aSVaclav Hapla   PetscFunctionBegin;
239812f479c1SVaclav Hapla   if (PetscDefined(USE_DEBUG)) {
2399157edd7aSVaclav Hapla     PetscSF dummy; /* just to have a PetscObject on comm for input validation */
2400157edd7aSVaclav Hapla 
24019566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(comm, &dummy));
2402157edd7aSVaclav Hapla     PetscValidLogicalCollectiveInt(dummy, nsfs, 2);
24034f572ea9SToby Isaac     PetscAssertPointer(sfs, 3);
2404157edd7aSVaclav Hapla     for (i = 0; i < nsfs; i++) {
2405157edd7aSVaclav Hapla       PetscValidHeaderSpecific(sfs[i], PETSCSF_CLASSID, 3);
2406157edd7aSVaclav Hapla       PetscCheckSameComm(dummy, 1, sfs[i], 3);
2407157edd7aSVaclav Hapla     }
24081f40158dSVaclav Hapla     PetscValidLogicalCollectiveEnum(dummy, rootMode, 4);
24094f572ea9SToby Isaac     if (leafOffsets) PetscAssertPointer(leafOffsets, 5);
24104f572ea9SToby Isaac     PetscAssertPointer(newsf, 6);
24119566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&dummy));
2412157edd7aSVaclav Hapla   }
2413157edd7aSVaclav Hapla   if (!nsfs) {
24149566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(comm, newsf));
24159566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*newsf, 0, 0, NULL, PETSC_OWN_POINTER, NULL, PETSC_OWN_POINTER));
24163ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2417157edd7aSVaclav Hapla   }
24189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
24191f40158dSVaclav Hapla   PetscCallMPI(MPI_Comm_size(comm, &size));
2420157edd7aSVaclav Hapla 
24211f40158dSVaclav Hapla   /* Calculate leaf array offsets */
24229566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsfs + 1, &leafArrayOffsets));
2423157edd7aSVaclav Hapla   leafArrayOffsets[0] = 0;
2424157edd7aSVaclav Hapla   for (s = 0; s < nsfs; s++) {
2425157edd7aSVaclav Hapla     PetscInt nl;
2426157edd7aSVaclav Hapla 
24279566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfs[s], NULL, &nl, NULL, NULL));
2428157edd7aSVaclav Hapla     leafArrayOffsets[s + 1] = leafArrayOffsets[s] + nl;
2429157edd7aSVaclav Hapla   }
2430157edd7aSVaclav Hapla   nLeaves = leafArrayOffsets[nsfs];
2431157edd7aSVaclav Hapla 
24321f40158dSVaclav Hapla   /* Calculate number of roots */
24331f40158dSVaclav Hapla   switch (rootMode) {
24341f40158dSVaclav Hapla   case PETSCSF_CONCATENATE_ROOTMODE_SHARED: {
24351f40158dSVaclav Hapla     PetscCall(PetscSFGetGraph(sfs[0], &nRoots, NULL, NULL, NULL));
24361f40158dSVaclav Hapla     if (PetscDefined(USE_DEBUG)) {
24371f40158dSVaclav Hapla       for (s = 1; s < nsfs; s++) {
24381f40158dSVaclav Hapla         PetscInt nr;
24391f40158dSVaclav Hapla 
24401f40158dSVaclav Hapla         PetscCall(PetscSFGetGraph(sfs[s], &nr, NULL, NULL, NULL));
24411f40158dSVaclav Hapla         PetscCheck(nr == nRoots, comm, PETSC_ERR_ARG_SIZ, "rootMode = %s but sfs[%" PetscInt_FMT "] has a different number of roots (%" PetscInt_FMT ") than sfs[0] (%" PetscInt_FMT ")", PetscSFConcatenateRootModes[rootMode], s, nr, nRoots);
24421f40158dSVaclav Hapla       }
24431f40158dSVaclav Hapla     }
24441f40158dSVaclav Hapla   } break;
24451f40158dSVaclav Hapla   case PETSCSF_CONCATENATE_ROOTMODE_GLOBAL: {
24461f40158dSVaclav Hapla     /* Calculate also global layout in this case */
24471f40158dSVaclav Hapla     PetscInt    *nls;
24481f40158dSVaclav Hapla     PetscLayout *lts;
24491f40158dSVaclav Hapla     PetscInt   **inds;
24501f40158dSVaclav Hapla     PetscInt     j;
24511f40158dSVaclav Hapla     PetscInt     rootOffset = 0;
24521f40158dSVaclav Hapla 
24531f40158dSVaclav Hapla     PetscCall(PetscCalloc3(nsfs, &lts, nsfs, &nls, nsfs, &inds));
24541f40158dSVaclav Hapla     PetscCall(PetscLayoutCreate(comm, &glayout));
24551f40158dSVaclav Hapla     glayout->bs = 1;
24561f40158dSVaclav Hapla     glayout->n  = 0;
24571f40158dSVaclav Hapla     glayout->N  = 0;
24581f40158dSVaclav Hapla     for (s = 0; s < nsfs; s++) {
24591f40158dSVaclav Hapla       PetscCall(PetscSFGetGraphLayout(sfs[s], &lts[s], &nls[s], NULL, &inds[s]));
24601f40158dSVaclav Hapla       glayout->n += lts[s]->n;
24611f40158dSVaclav Hapla       glayout->N += lts[s]->N;
24621f40158dSVaclav Hapla     }
24631f40158dSVaclav Hapla     PetscCall(PetscLayoutSetUp(glayout));
24641f40158dSVaclav Hapla     PetscCall(PetscMalloc1(nLeaves, &gremote));
24651f40158dSVaclav Hapla     for (s = 0, j = 0; s < nsfs; s++) {
24661f40158dSVaclav Hapla       for (i = 0; i < nls[s]; i++, j++) gremote[j] = inds[s][i] + rootOffset;
24671f40158dSVaclav Hapla       rootOffset += lts[s]->N;
24681f40158dSVaclav Hapla       PetscCall(PetscLayoutDestroy(&lts[s]));
24691f40158dSVaclav Hapla       PetscCall(PetscFree(inds[s]));
24701f40158dSVaclav Hapla     }
24711f40158dSVaclav Hapla     PetscCall(PetscFree3(lts, nls, inds));
24721f40158dSVaclav Hapla     nRoots = glayout->N;
24731f40158dSVaclav Hapla   } break;
24741f40158dSVaclav Hapla   case PETSCSF_CONCATENATE_ROOTMODE_LOCAL:
24751f40158dSVaclav Hapla     /* nRoots calculated later in this case */
24761f40158dSVaclav Hapla     break;
24771f40158dSVaclav Hapla   default:
24781f40158dSVaclav Hapla     SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid PetscSFConcatenateRootMode %d", rootMode);
24791f40158dSVaclav Hapla   }
24801f40158dSVaclav Hapla 
2481157edd7aSVaclav Hapla   if (!leafOffsets) {
2482157edd7aSVaclav Hapla     all_ilocal_null = PETSC_TRUE;
2483157edd7aSVaclav Hapla     for (s = 0; s < nsfs; s++) {
2484157edd7aSVaclav Hapla       const PetscInt *ilocal;
2485157edd7aSVaclav Hapla 
24869566063dSJacob Faibussowitsch       PetscCall(PetscSFGetGraph(sfs[s], NULL, NULL, &ilocal, NULL));
2487157edd7aSVaclav Hapla       if (ilocal) {
2488157edd7aSVaclav Hapla         all_ilocal_null = PETSC_FALSE;
2489157edd7aSVaclav Hapla         break;
2490157edd7aSVaclav Hapla       }
2491157edd7aSVaclav Hapla     }
2492157edd7aSVaclav 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");
2493157edd7aSVaclav Hapla   }
2494157edd7aSVaclav Hapla 
2495157edd7aSVaclav Hapla   /* Renumber and concatenate local leaves */
2496157edd7aSVaclav Hapla   ilocal_new = NULL;
2497157edd7aSVaclav Hapla   if (!all_ilocal_null) {
24989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nLeaves, &ilocal_new));
2499157edd7aSVaclav Hapla     for (i = 0; i < nLeaves; i++) ilocal_new[i] = -1;
2500157edd7aSVaclav Hapla     for (s = 0; s < nsfs; s++) {
2501157edd7aSVaclav Hapla       const PetscInt *ilocal;
25028e3a54c0SPierre Jolivet       PetscInt       *ilocal_l = PetscSafePointerPlusOffset(ilocal_new, leafArrayOffsets[s]);
2503157edd7aSVaclav Hapla       PetscInt        i, nleaves_l;
2504157edd7aSVaclav Hapla 
25059566063dSJacob Faibussowitsch       PetscCall(PetscSFGetGraph(sfs[s], NULL, &nleaves_l, &ilocal, NULL));
2506157edd7aSVaclav Hapla       for (i = 0; i < nleaves_l; i++) ilocal_l[i] = (ilocal ? ilocal[i] : i) + leafOffsets[s];
2507157edd7aSVaclav Hapla     }
2508157edd7aSVaclav Hapla   }
2509157edd7aSVaclav Hapla 
2510157edd7aSVaclav Hapla   /* Renumber and concatenate remote roots */
25111f40158dSVaclav Hapla   if (rootMode == PETSCSF_CONCATENATE_ROOTMODE_LOCAL || rootMode == PETSCSF_CONCATENATE_ROOTMODE_SHARED) {
25121f40158dSVaclav Hapla     PetscInt rootOffset = 0;
25131f40158dSVaclav Hapla 
25149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nLeaves, &iremote_new));
2515157edd7aSVaclav Hapla     for (i = 0; i < nLeaves; i++) {
2516157edd7aSVaclav Hapla       iremote_new[i].rank  = -1;
2517157edd7aSVaclav Hapla       iremote_new[i].index = -1;
2518157edd7aSVaclav Hapla     }
2519157edd7aSVaclav Hapla     for (s = 0; s < nsfs; s++) {
2520157edd7aSVaclav Hapla       PetscInt           i, nl, nr;
2521157edd7aSVaclav Hapla       PetscSF            tmp_sf;
2522157edd7aSVaclav Hapla       const PetscSFNode *iremote;
2523157edd7aSVaclav Hapla       PetscSFNode       *tmp_rootdata;
25248e3a54c0SPierre Jolivet       PetscSFNode       *tmp_leafdata = PetscSafePointerPlusOffset(iremote_new, leafArrayOffsets[s]);
2525157edd7aSVaclav Hapla 
25269566063dSJacob Faibussowitsch       PetscCall(PetscSFGetGraph(sfs[s], &nr, &nl, NULL, &iremote));
25279566063dSJacob Faibussowitsch       PetscCall(PetscSFCreate(comm, &tmp_sf));
2528157edd7aSVaclav Hapla       /* create helper SF with contiguous leaves */
25299566063dSJacob Faibussowitsch       PetscCall(PetscSFSetGraph(tmp_sf, nr, nl, NULL, PETSC_USE_POINTER, (PetscSFNode *)iremote, PETSC_COPY_VALUES));
25309566063dSJacob Faibussowitsch       PetscCall(PetscSFSetUp(tmp_sf));
25319566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nr, &tmp_rootdata));
25321f40158dSVaclav Hapla       if (rootMode == PETSCSF_CONCATENATE_ROOTMODE_LOCAL) {
2533157edd7aSVaclav Hapla         for (i = 0; i < nr; i++) {
25341f40158dSVaclav Hapla           tmp_rootdata[i].index = i + rootOffset;
2535157edd7aSVaclav Hapla           tmp_rootdata[i].rank  = (PetscInt)rank;
2536157edd7aSVaclav Hapla         }
25371f40158dSVaclav Hapla         rootOffset += nr;
25381f40158dSVaclav Hapla       } else {
25391f40158dSVaclav Hapla         for (i = 0; i < nr; i++) {
25401f40158dSVaclav Hapla           tmp_rootdata[i].index = i;
25411f40158dSVaclav Hapla           tmp_rootdata[i].rank  = (PetscInt)rank;
25421f40158dSVaclav Hapla         }
25431f40158dSVaclav Hapla       }
25449566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(tmp_sf, MPIU_2INT, tmp_rootdata, tmp_leafdata, MPI_REPLACE));
25459566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(tmp_sf, MPIU_2INT, tmp_rootdata, tmp_leafdata, MPI_REPLACE));
25469566063dSJacob Faibussowitsch       PetscCall(PetscSFDestroy(&tmp_sf));
25479566063dSJacob Faibussowitsch       PetscCall(PetscFree(tmp_rootdata));
2548157edd7aSVaclav Hapla     }
2549aa624791SPierre Jolivet     if (rootMode == PETSCSF_CONCATENATE_ROOTMODE_LOCAL) nRoots = rootOffset; // else nRoots already calculated above
2550157edd7aSVaclav Hapla 
2551157edd7aSVaclav Hapla     /* Build the new SF */
25529566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(comm, newsf));
25539566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*newsf, nRoots, nLeaves, ilocal_new, PETSC_OWN_POINTER, iremote_new, PETSC_OWN_POINTER));
25541f40158dSVaclav Hapla   } else {
25551f40158dSVaclav Hapla     /* Build the new SF */
25561f40158dSVaclav Hapla     PetscCall(PetscSFCreate(comm, newsf));
25571f40158dSVaclav Hapla     PetscCall(PetscSFSetGraphLayout(*newsf, glayout, nLeaves, ilocal_new, PETSC_OWN_POINTER, gremote));
25581f40158dSVaclav Hapla   }
25599566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(*newsf));
25601f40158dSVaclav Hapla   PetscCall(PetscSFViewFromOptions(*newsf, NULL, "-sf_concat_view"));
25611f40158dSVaclav Hapla   PetscCall(PetscLayoutDestroy(&glayout));
25621f40158dSVaclav Hapla   PetscCall(PetscFree(gremote));
25639566063dSJacob Faibussowitsch   PetscCall(PetscFree(leafArrayOffsets));
25643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2565157edd7aSVaclav Hapla }
2566