xref: /petsc/src/vec/is/sf/interface/sf.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
1af0996ceSBarry Smith #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/
2c4e6a40aSLawrence Mitchell #include <petsc/private/hashseti.h>
353dd6d7dSJunchao Zhang #include <petsc/private/viewerimpl.h>
495fce210SBarry Smith #include <petscctable.h>
595fce210SBarry Smith 
67fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_CUDA)
77fd2d3dbSJunchao Zhang   #include <cuda_runtime.h>
87fd2d3dbSJunchao Zhang #endif
97fd2d3dbSJunchao Zhang 
107fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_HIP)
117fd2d3dbSJunchao Zhang   #include <hip/hip_runtime.h>
127fd2d3dbSJunchao Zhang #endif
137fd2d3dbSJunchao Zhang 
142abc8c78SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER)
152abc8c78SJacob Faibussowitsch void PetscSFCheckGraphSet(PetscSF,int);
162abc8c78SJacob Faibussowitsch #else
1795fce210SBarry Smith #if defined(PETSC_USE_DEBUG)
187a46b595SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) PetscCheck((sf)->graphset,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetGraph() or PetscSFSetGraphWithPattern() on argument %d \"%s\" before %s()",(arg),#sf,PETSC_FUNCTION_NAME);
1995fce210SBarry Smith #else
2095fce210SBarry Smith #  define PetscSFCheckGraphSet(sf,arg) do {} while (0)
2195fce210SBarry Smith #endif
222abc8c78SJacob Faibussowitsch #endif
2395fce210SBarry Smith 
244c8fdceaSLisandro Dalcin const char *const PetscSFDuplicateOptions[] = {"CONFONLY","RANKS","GRAPH","PetscSFDuplicateOption","PETSCSF_DUPLICATE_",NULL};
2595fce210SBarry Smith 
268af6ec1cSBarry Smith /*@
2795fce210SBarry Smith    PetscSFCreate - create a star forest communication context
2895fce210SBarry Smith 
29d083f849SBarry Smith    Collective
3095fce210SBarry Smith 
314165533cSJose E. Roman    Input Parameter:
3295fce210SBarry Smith .  comm - communicator on which the star forest will operate
3395fce210SBarry Smith 
344165533cSJose E. Roman    Output Parameter:
3595fce210SBarry Smith .  sf - new star forest context
3695fce210SBarry Smith 
37dd5b3ca6SJunchao Zhang    Options Database Keys:
38dd5b3ca6SJunchao Zhang +  -sf_type basic     -Use MPI persistent Isend/Irecv for communication (Default)
39dd5b3ca6SJunchao Zhang .  -sf_type window    -Use MPI-3 one-sided window for communication
40dd5b3ca6SJunchao Zhang -  -sf_type neighbor  -Use MPI-3 neighborhood collectives for communication
41dd5b3ca6SJunchao Zhang 
4295fce210SBarry Smith    Level: intermediate
4395fce210SBarry Smith 
44dd5b3ca6SJunchao Zhang    Notes:
45dd5b3ca6SJunchao Zhang    When one knows the communication graph is one of the predefined graph, such as MPI_Alltoall, MPI_Allgatherv,
46dd5b3ca6SJunchao Zhang    MPI_Gatherv, one can create a PetscSF and then set its graph with PetscSFSetGraphWithPattern(). These special
47dd5b3ca6SJunchao Zhang    SFs are optimized and they have better performance than general SFs.
48dd5b3ca6SJunchao Zhang 
49db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`, `PetscSFSetGraphWithPattern()`, `PetscSFDestroy()`
5095fce210SBarry Smith @*/
5195fce210SBarry Smith PetscErrorCode PetscSFCreate(MPI_Comm comm,PetscSF *sf)
5295fce210SBarry Smith {
5395fce210SBarry Smith   PetscSF        b;
5495fce210SBarry Smith 
5595fce210SBarry Smith   PetscFunctionBegin;
5695fce210SBarry Smith   PetscValidPointer(sf,2);
579566063dSJacob Faibussowitsch   PetscCall(PetscSFInitializePackage());
5895fce210SBarry Smith 
599566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(b,PETSCSF_CLASSID,"PetscSF","Star Forest","PetscSF",comm,PetscSFDestroy,PetscSFView));
6095fce210SBarry Smith 
6195fce210SBarry Smith   b->nroots    = -1;
6295fce210SBarry Smith   b->nleaves   = -1;
6329046d53SLisandro Dalcin   b->minleaf   = PETSC_MAX_INT;
6429046d53SLisandro Dalcin   b->maxleaf   = PETSC_MIN_INT;
6595fce210SBarry Smith   b->nranks    = -1;
6695fce210SBarry Smith   b->rankorder = PETSC_TRUE;
6795fce210SBarry Smith   b->ingroup   = MPI_GROUP_NULL;
6895fce210SBarry Smith   b->outgroup  = MPI_GROUP_NULL;
6995fce210SBarry Smith   b->graphset  = PETSC_FALSE;
7020c24465SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
7120c24465SJunchao Zhang   b->use_gpu_aware_mpi    = use_gpu_aware_mpi;
7220c24465SJunchao Zhang   b->use_stream_aware_mpi = PETSC_FALSE;
7371438e86SJunchao Zhang   b->unknown_input_stream= PETSC_FALSE;
7427f636e8SJunchao Zhang   #if defined(PETSC_HAVE_KOKKOS) /* Prefer kokkos over cuda*/
7520c24465SJunchao Zhang     b->backend = PETSCSF_BACKEND_KOKKOS;
7627f636e8SJunchao Zhang   #elif defined(PETSC_HAVE_CUDA)
7727f636e8SJunchao Zhang     b->backend = PETSCSF_BACKEND_CUDA;
7859af0bd3SScott Kruger   #elif defined(PETSC_HAVE_HIP)
7959af0bd3SScott Kruger     b->backend = PETSCSF_BACKEND_HIP;
8020c24465SJunchao Zhang   #endif
8171438e86SJunchao Zhang 
8271438e86SJunchao Zhang   #if defined(PETSC_HAVE_NVSHMEM)
8371438e86SJunchao Zhang     b->use_nvshmem     = PETSC_FALSE; /* Default is not to try NVSHMEM */
8471438e86SJunchao Zhang     b->use_nvshmem_get = PETSC_FALSE; /* Default is to use nvshmem_put based protocol */
859566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_nvshmem",&b->use_nvshmem,NULL));
869566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_nvshmem_get",&b->use_nvshmem_get,NULL));
8771438e86SJunchao Zhang   #endif
8820c24465SJunchao Zhang #endif
8960c22052SBarry Smith   b->vscat.from_n = -1;
9060c22052SBarry Smith   b->vscat.to_n   = -1;
9160c22052SBarry Smith   b->vscat.unit   = MPIU_SCALAR;
9295fce210SBarry Smith  *sf = b;
9395fce210SBarry Smith   PetscFunctionReturn(0);
9495fce210SBarry Smith }
9595fce210SBarry Smith 
9629046d53SLisandro Dalcin /*@
9795fce210SBarry Smith    PetscSFReset - Reset a star forest so that different sizes or neighbors can be used
9895fce210SBarry Smith 
9995fce210SBarry Smith    Collective
10095fce210SBarry Smith 
1014165533cSJose E. Roman    Input Parameter:
10295fce210SBarry Smith .  sf - star forest
10395fce210SBarry Smith 
10495fce210SBarry Smith    Level: advanced
10595fce210SBarry Smith 
106db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFSetGraph()`, `PetscSFDestroy()`
10795fce210SBarry Smith @*/
10895fce210SBarry Smith PetscErrorCode PetscSFReset(PetscSF sf)
10995fce210SBarry Smith {
11095fce210SBarry Smith   PetscFunctionBegin;
11195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
112*dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf,Reset);
11329046d53SLisandro Dalcin   sf->nroots   = -1;
11429046d53SLisandro Dalcin   sf->nleaves  = -1;
11529046d53SLisandro Dalcin   sf->minleaf  = PETSC_MAX_INT;
11629046d53SLisandro Dalcin   sf->maxleaf  = PETSC_MIN_INT;
11795fce210SBarry Smith   sf->mine     = NULL;
11895fce210SBarry Smith   sf->remote   = NULL;
11929046d53SLisandro Dalcin   sf->graphset = PETSC_FALSE;
1209566063dSJacob Faibussowitsch   PetscCall(PetscFree(sf->mine_alloc));
1219566063dSJacob Faibussowitsch   PetscCall(PetscFree(sf->remote_alloc));
12221c688dcSJed Brown   sf->nranks = -1;
1239566063dSJacob Faibussowitsch   PetscCall(PetscFree4(sf->ranks,sf->roffset,sf->rmine,sf->rremote));
12429046d53SLisandro Dalcin   sf->degreeknown = PETSC_FALSE;
1259566063dSJacob Faibussowitsch   PetscCall(PetscFree(sf->degree));
1269566063dSJacob Faibussowitsch   if (sf->ingroup  != MPI_GROUP_NULL) PetscCallMPI(MPI_Group_free(&sf->ingroup));
1279566063dSJacob Faibussowitsch   if (sf->outgroup != MPI_GROUP_NULL) PetscCallMPI(MPI_Group_free(&sf->outgroup));
128013b3241SStefano Zampini   if (sf->multi) sf->multi->multi = NULL;
1299566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf->multi));
1309566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&sf->map));
13171438e86SJunchao Zhang 
13271438e86SJunchao Zhang  #if defined(PETSC_HAVE_DEVICE)
1339566063dSJacob Faibussowitsch   for (PetscInt i=0; i<2; i++) PetscCall(PetscSFFree(sf,PETSC_MEMTYPE_DEVICE,sf->rmine_d[i]));
13471438e86SJunchao Zhang  #endif
13571438e86SJunchao Zhang 
13695fce210SBarry Smith   sf->setupcalled = PETSC_FALSE;
13795fce210SBarry Smith   PetscFunctionReturn(0);
13895fce210SBarry Smith }
13995fce210SBarry Smith 
14095fce210SBarry Smith /*@C
14129046d53SLisandro Dalcin    PetscSFSetType - Set the PetscSF communication implementation
14295fce210SBarry Smith 
14395fce210SBarry Smith    Collective on PetscSF
14495fce210SBarry Smith 
14595fce210SBarry Smith    Input Parameters:
14695fce210SBarry Smith +  sf - the PetscSF context
14795fce210SBarry Smith -  type - a known method
14895fce210SBarry Smith 
14995fce210SBarry Smith    Options Database Key:
15095fce210SBarry Smith .  -sf_type <type> - Sets the method; use -help for a list
15170616304SStefano Zampini    of available methods (for instance, window, basic, neighbor)
15295fce210SBarry Smith 
15395fce210SBarry Smith    Notes:
15495fce210SBarry Smith    See "include/petscsf.h" for available methods (for instance)
15595fce210SBarry Smith +    PETSCSFWINDOW - MPI-2/3 one-sided
15695fce210SBarry Smith -    PETSCSFBASIC - basic implementation using MPI-1 two-sided
15795fce210SBarry Smith 
15895fce210SBarry Smith   Level: intermediate
15995fce210SBarry Smith 
160db781477SPatrick Sanan .seealso: `PetscSFType`, `PetscSFCreate()`
16195fce210SBarry Smith @*/
16295fce210SBarry Smith PetscErrorCode PetscSFSetType(PetscSF sf,PetscSFType type)
16395fce210SBarry Smith {
16495fce210SBarry Smith   PetscBool      match;
1655f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(PetscSF);
16695fce210SBarry Smith 
16795fce210SBarry Smith   PetscFunctionBegin;
16895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
16995fce210SBarry Smith   PetscValidCharPointer(type,2);
17095fce210SBarry Smith 
1719566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)sf,type,&match));
17295fce210SBarry Smith   if (match) PetscFunctionReturn(0);
17395fce210SBarry Smith 
1749566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(PetscSFList,type,&r));
17528b400f6SJacob Faibussowitsch   PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PetscSF type %s",type);
17629046d53SLisandro Dalcin   /* Destroy the previous PetscSF implementation context */
177*dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf,Destroy);
1789566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(sf->ops,sizeof(*sf->ops)));
1799566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)sf,type));
1809566063dSJacob Faibussowitsch   PetscCall((*r)(sf));
18195fce210SBarry Smith   PetscFunctionReturn(0);
18295fce210SBarry Smith }
18395fce210SBarry Smith 
18429046d53SLisandro Dalcin /*@C
18529046d53SLisandro Dalcin   PetscSFGetType - Get the PetscSF communication implementation
18629046d53SLisandro Dalcin 
18729046d53SLisandro Dalcin   Not Collective
18829046d53SLisandro Dalcin 
18929046d53SLisandro Dalcin   Input Parameter:
19029046d53SLisandro Dalcin . sf  - the PetscSF context
19129046d53SLisandro Dalcin 
19229046d53SLisandro Dalcin   Output Parameter:
19329046d53SLisandro Dalcin . type - the PetscSF type name
19429046d53SLisandro Dalcin 
19529046d53SLisandro Dalcin   Level: intermediate
19629046d53SLisandro Dalcin 
197db781477SPatrick Sanan .seealso: `PetscSFSetType()`, `PetscSFCreate()`
19829046d53SLisandro Dalcin @*/
19929046d53SLisandro Dalcin PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type)
20029046d53SLisandro Dalcin {
20129046d53SLisandro Dalcin   PetscFunctionBegin;
20229046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID,1);
20329046d53SLisandro Dalcin   PetscValidPointer(type,2);
20429046d53SLisandro Dalcin   *type = ((PetscObject)sf)->type_name;
20529046d53SLisandro Dalcin   PetscFunctionReturn(0);
20629046d53SLisandro Dalcin }
20729046d53SLisandro Dalcin 
2081fb7b255SJunchao Zhang /*@C
20995fce210SBarry Smith    PetscSFDestroy - destroy star forest
21095fce210SBarry Smith 
21195fce210SBarry Smith    Collective
21295fce210SBarry Smith 
2134165533cSJose E. Roman    Input Parameter:
21495fce210SBarry Smith .  sf - address of star forest
21595fce210SBarry Smith 
21695fce210SBarry Smith    Level: intermediate
21795fce210SBarry Smith 
218db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFReset()`
21995fce210SBarry Smith @*/
22095fce210SBarry Smith PetscErrorCode PetscSFDestroy(PetscSF *sf)
22195fce210SBarry Smith {
22295fce210SBarry Smith   PetscFunctionBegin;
22395fce210SBarry Smith   if (!*sf) PetscFunctionReturn(0);
22495fce210SBarry Smith   PetscValidHeaderSpecific((*sf),PETSCSF_CLASSID,1);
22529046d53SLisandro Dalcin   if (--((PetscObject)(*sf))->refct > 0) {*sf = NULL; PetscFunctionReturn(0);}
2269566063dSJacob Faibussowitsch   PetscCall(PetscSFReset(*sf));
227*dbbe0bcdSBarry Smith   PetscTryTypeMethod((*sf),Destroy);
2289566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&(*sf)->vscat.lsf));
2299566063dSJacob Faibussowitsch   if ((*sf)->vscat.bs > 1) PetscCallMPI(MPI_Type_free(&(*sf)->vscat.unit));
2309566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(sf));
23195fce210SBarry Smith   PetscFunctionReturn(0);
23295fce210SBarry Smith }
23395fce210SBarry Smith 
234c4e6a40aSLawrence Mitchell static PetscErrorCode PetscSFCheckGraphValid_Private(PetscSF sf)
235c4e6a40aSLawrence Mitchell {
236c4e6a40aSLawrence Mitchell   PetscInt           i, nleaves;
237c4e6a40aSLawrence Mitchell   PetscMPIInt        size;
238c4e6a40aSLawrence Mitchell   const PetscInt    *ilocal;
239c4e6a40aSLawrence Mitchell   const PetscSFNode *iremote;
240c4e6a40aSLawrence Mitchell 
241c4e6a40aSLawrence Mitchell   PetscFunctionBegin;
24276bd3646SJed Brown   if (!sf->graphset || !PetscDefined(USE_DEBUG)) PetscFunctionReturn(0);
2439566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,&iremote));
2449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size));
245c4e6a40aSLawrence Mitchell   for (i = 0; i < nleaves; i++) {
246c4e6a40aSLawrence Mitchell     const PetscInt rank = iremote[i].rank;
247c4e6a40aSLawrence Mitchell     const PetscInt remote = iremote[i].index;
248c4e6a40aSLawrence Mitchell     const PetscInt leaf = ilocal ? ilocal[i] : i;
249c9cc58a2SBarry 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);
25008401ef6SPierre 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);
25108401ef6SPierre 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);
252c4e6a40aSLawrence Mitchell   }
253c4e6a40aSLawrence Mitchell   PetscFunctionReturn(0);
254c4e6a40aSLawrence Mitchell }
255c4e6a40aSLawrence Mitchell 
25695fce210SBarry Smith /*@
25795fce210SBarry Smith    PetscSFSetUp - set up communication structures
25895fce210SBarry Smith 
25995fce210SBarry Smith    Collective
26095fce210SBarry Smith 
2614165533cSJose E. Roman    Input Parameter:
26295fce210SBarry Smith .  sf - star forest communication object
26395fce210SBarry Smith 
26495fce210SBarry Smith    Level: beginner
26595fce210SBarry Smith 
266db781477SPatrick Sanan .seealso: `PetscSFSetFromOptions()`, `PetscSFSetType()`
26795fce210SBarry Smith @*/
26895fce210SBarry Smith PetscErrorCode PetscSFSetUp(PetscSF sf)
26995fce210SBarry Smith {
27095fce210SBarry Smith   PetscFunctionBegin;
27129046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
27229046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
27395fce210SBarry Smith   if (sf->setupcalled) PetscFunctionReturn(0);
2749566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_SetUp,sf,0,0,0));
2759566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckGraphValid_Private(sf));
2769566063dSJacob Faibussowitsch   if (!((PetscObject)sf)->type_name) PetscCall(PetscSFSetType(sf,PETSCSFBASIC)); /* Zero all sf->ops */
277*dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf,SetUp);
27820c24465SJunchao Zhang #if defined(PETSC_HAVE_CUDA)
27920c24465SJunchao Zhang   if (sf->backend == PETSCSF_BACKEND_CUDA) {
28071438e86SJunchao Zhang     sf->ops->Malloc = PetscSFMalloc_CUDA;
28171438e86SJunchao Zhang     sf->ops->Free   = PetscSFFree_CUDA;
28220c24465SJunchao Zhang   }
28320c24465SJunchao Zhang #endif
28459af0bd3SScott Kruger #if defined(PETSC_HAVE_HIP)
28559af0bd3SScott Kruger   if (sf->backend == PETSCSF_BACKEND_HIP) {
28659af0bd3SScott Kruger     sf->ops->Malloc = PetscSFMalloc_HIP;
28759af0bd3SScott Kruger     sf->ops->Free   = PetscSFFree_HIP;
28859af0bd3SScott Kruger   }
28959af0bd3SScott Kruger #endif
29020c24465SJunchao Zhang 
29159af0bd3SScott Kruger #
29220c24465SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
29320c24465SJunchao Zhang   if (sf->backend == PETSCSF_BACKEND_KOKKOS) {
29420c24465SJunchao Zhang     sf->ops->Malloc = PetscSFMalloc_Kokkos;
29520c24465SJunchao Zhang     sf->ops->Free   = PetscSFFree_Kokkos;
29620c24465SJunchao Zhang   }
29720c24465SJunchao Zhang #endif
2989566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_SetUp,sf,0,0,0));
29995fce210SBarry Smith   sf->setupcalled = PETSC_TRUE;
30095fce210SBarry Smith   PetscFunctionReturn(0);
30195fce210SBarry Smith }
30295fce210SBarry Smith 
3038af6ec1cSBarry Smith /*@
30495fce210SBarry Smith    PetscSFSetFromOptions - set PetscSF options using the options database
30595fce210SBarry Smith 
30695fce210SBarry Smith    Logically Collective
30795fce210SBarry Smith 
3084165533cSJose E. Roman    Input Parameter:
30995fce210SBarry Smith .  sf - star forest
31095fce210SBarry Smith 
31195fce210SBarry Smith    Options Database Keys:
31260263706SJed Brown +  -sf_type               - implementation type, see PetscSFSetType()
31351ccb202SJunchao Zhang .  -sf_rank_order         - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise
314b85e67b7SJunchao Zhang .  -sf_use_default_stream - Assume callers of SF computed the input root/leafdata with the default cuda stream. SF will also
315c2a741eeSJunchao Zhang                             use the default stream to process data. Therefore, no stream synchronization is needed between SF and its caller (default: true).
316c06a8e02SRichard Tran Mills                             If true, this option only works with -use_gpu_aware_mpi 1.
31720c24465SJunchao Zhang .  -sf_use_stream_aware_mpi  - Assume the underlying MPI is cuda-stream aware and SF won't sync streams for send/recv buffers passed to MPI (default: false).
318c06a8e02SRichard Tran Mills                                If true, this option only works with -use_gpu_aware_mpi 1.
31995fce210SBarry Smith 
32059af0bd3SScott Kruger -  -sf_backend cuda | hip | kokkos -Select the device backend SF uses. Currently SF has these backends: cuda, hip and Kokkos.
32159af0bd3SScott Kruger                               On CUDA (HIP) devices, one can choose cuda (hip) or kokkos with the default being kokkos. On other devices,
32220c24465SJunchao Zhang                               the only available is kokkos.
32320c24465SJunchao Zhang 
32495fce210SBarry Smith    Level: intermediate
32595fce210SBarry Smith @*/
32695fce210SBarry Smith PetscErrorCode PetscSFSetFromOptions(PetscSF sf)
32795fce210SBarry Smith {
32895fce210SBarry Smith   PetscSFType    deft;
32995fce210SBarry Smith   char           type[256];
33095fce210SBarry Smith   PetscBool      flg;
33195fce210SBarry Smith 
33295fce210SBarry Smith   PetscFunctionBegin;
33395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
334d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)sf);
33595fce210SBarry Smith   deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC;
3369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-sf_type","PetscSF implementation type","PetscSFSetType",PetscSFList,deft,type,sizeof(type),&flg));
3379566063dSJacob Faibussowitsch   PetscCall(PetscSFSetType(sf,flg ? type : deft));
3389566063dSJacob 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));
3397fd2d3dbSJunchao Zhang  #if defined(PETSC_HAVE_DEVICE)
34020c24465SJunchao Zhang   {
34120c24465SJunchao Zhang     char        backendstr[32] = {0};
34259af0bd3SScott Kruger     PetscBool   isCuda = PETSC_FALSE,isHip = PETSC_FALSE,isKokkos = PETSC_FALSE,set;
34320c24465SJunchao Zhang     /* Change the defaults set in PetscSFCreate() with command line options */
3449566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-sf_unknown_input_stream","SF root/leafdata is computed on arbitary streams unknown to SF","PetscSFSetFromOptions",sf->unknown_input_stream,&sf->unknown_input_stream,NULL));
3459566063dSJacob 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));
3469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsString("-sf_backend","Select the device backend SF uses","PetscSFSetFromOptions",NULL,backendstr,sizeof(backendstr),&set));
3479566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp("cuda",backendstr,&isCuda));
3489566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp("kokkos",backendstr,&isKokkos));
3499566063dSJacob Faibussowitsch     PetscCall(PetscStrcasecmp("hip",backendstr,&isHip));
35059af0bd3SScott Kruger   #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
35120c24465SJunchao Zhang     if (isCuda) sf->backend = PETSCSF_BACKEND_CUDA;
35220c24465SJunchao Zhang     else if (isKokkos) sf->backend = PETSCSF_BACKEND_KOKKOS;
35359af0bd3SScott Kruger     else if (isHip) sf->backend = PETSCSF_BACKEND_HIP;
35428b400f6SJacob 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);
35520c24465SJunchao Zhang   #elif defined(PETSC_HAVE_KOKKOS)
35608401ef6SPierre Jolivet     PetscCheck(!set || isKokkos,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"-sf_backend %s is not supported. You can only choose kokkos", backendstr);
35720c24465SJunchao Zhang    #endif
35820c24465SJunchao Zhang   }
359c2a741eeSJunchao Zhang  #endif
360*dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf,SetFromOptions,PetscOptionsObject);
361d0609cedSBarry Smith   PetscOptionsEnd();
36295fce210SBarry Smith   PetscFunctionReturn(0);
36395fce210SBarry Smith }
36495fce210SBarry Smith 
36529046d53SLisandro Dalcin /*@
36695fce210SBarry Smith    PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order
36795fce210SBarry Smith 
36895fce210SBarry Smith    Logically Collective
36995fce210SBarry Smith 
3704165533cSJose E. Roman    Input Parameters:
37195fce210SBarry Smith +  sf - star forest
37295fce210SBarry Smith -  flg - PETSC_TRUE to sort, PETSC_FALSE to skip sorting (lower setup cost, but non-deterministic)
37395fce210SBarry Smith 
37495fce210SBarry Smith    Level: advanced
37595fce210SBarry Smith 
376db781477SPatrick Sanan .seealso: `PetscSFGatherBegin()`, `PetscSFScatterBegin()`
37795fce210SBarry Smith @*/
37895fce210SBarry Smith PetscErrorCode PetscSFSetRankOrder(PetscSF sf,PetscBool flg)
37995fce210SBarry Smith {
38095fce210SBarry Smith   PetscFunctionBegin;
38195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
38295fce210SBarry Smith   PetscValidLogicalCollectiveBool(sf,flg,2);
38328b400f6SJacob Faibussowitsch   PetscCheck(!sf->multi,PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()");
38495fce210SBarry Smith   sf->rankorder = flg;
38595fce210SBarry Smith   PetscFunctionReturn(0);
38695fce210SBarry Smith }
38795fce210SBarry Smith 
3888dbb0df6SBarry Smith /*@C
38995fce210SBarry Smith    PetscSFSetGraph - Set a parallel star forest
39095fce210SBarry Smith 
39195fce210SBarry Smith    Collective
39295fce210SBarry Smith 
3934165533cSJose E. Roman    Input Parameters:
39495fce210SBarry Smith +  sf - star forest
39595fce210SBarry Smith .  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
39695fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
397c4e6a40aSLawrence Mitchell .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage (locations must be >= 0, enforced
398c4e6a40aSLawrence Mitchell during setup in debug mode)
39995fce210SBarry Smith .  localmode - copy mode for ilocal
400c4e6a40aSLawrence Mitchell .  iremote - remote locations of root vertices for each leaf on the current process (locations must be >= 0, enforced
401c4e6a40aSLawrence Mitchell during setup in debug mode)
40295fce210SBarry Smith -  remotemode - copy mode for iremote
40395fce210SBarry Smith 
40495fce210SBarry Smith    Level: intermediate
40595fce210SBarry Smith 
40695452b02SPatrick Sanan    Notes:
407db2b9530SVaclav Hapla    Leaf indices in ilocal must be unique, otherwise an error occurs.
40838ab3f8aSBarry Smith 
409db2b9530SVaclav Hapla    Input arrays ilocal and iremote follow the PetscCopyMode semantics.
410db2b9530SVaclav Hapla    In particular, if localmode/remotemode is PETSC_OWN_POINTER or PETSC_USE_POINTER,
411db2b9530SVaclav Hapla    PETSc might modify the respective array;
412db2b9530SVaclav Hapla    if PETSC_USE_POINTER, the user must delete the array after PetscSFDestroy().
413e97d3de3SVaclav Hapla    Only if PETSC_COPY_VALUES is used, the respective array is guaranteed to stay intact and a const array can be passed (but a cast to non-const is needed).
414db2b9530SVaclav Hapla 
415db2b9530SVaclav Hapla    Fortran Notes:
416db2b9530SVaclav Hapla    In Fortran you must use PETSC_COPY_VALUES for localmode and remotemode.
417c4e6a40aSLawrence Mitchell 
41872bf8598SVaclav Hapla    Developer Notes:
419db2b9530SVaclav Hapla    We sort leaves to check for duplicates and contiguousness and to find minleaf/maxleaf.
420db2b9530SVaclav Hapla    This also allows to compare leaf sets of two SFs easily.
42172bf8598SVaclav Hapla 
422db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFView()`, `PetscSFGetGraph()`
42395fce210SBarry Smith @*/
424db2b9530SVaclav Hapla PetscErrorCode PetscSFSetGraph(PetscSF sf,PetscInt nroots,PetscInt nleaves,PetscInt *ilocal,PetscCopyMode localmode,PetscSFNode *iremote,PetscCopyMode remotemode)
42595fce210SBarry Smith {
426db2b9530SVaclav Hapla   PetscBool       unique, contiguous;
42795fce210SBarry Smith 
42895fce210SBarry Smith   PetscFunctionBegin;
42995fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
43029046d53SLisandro Dalcin   if (nleaves > 0 && ilocal) PetscValidIntPointer(ilocal,4);
43129046d53SLisandro Dalcin   if (nleaves > 0) PetscValidPointer(iremote,6);
43208401ef6SPierre Jolivet   PetscCheck(nroots  >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nroots %" PetscInt_FMT ", cannot be negative",nroots);
43308401ef6SPierre Jolivet   PetscCheck(nleaves >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nleaves %" PetscInt_FMT ", cannot be negative",nleaves);
4348da24d32SBarry Smith   /* enums may be handled as unsigned by some compilers, NVHPC for example, the int cast
4358da24d32SBarry Smith    * below is to prevent NVHPC from warning about meaningless comparison of unsigned with zero */
4368da24d32SBarry Smith   PetscCheck((int)localmode  >= PETSC_COPY_VALUES && localmode  <= PETSC_USE_POINTER,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Wrong localmode %d",localmode);
4378da24d32SBarry Smith   PetscCheck((int)remotemode >= PETSC_COPY_VALUES && remotemode <= PETSC_USE_POINTER,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Wrong remotemode %d",remotemode);
43829046d53SLisandro Dalcin 
4392a67d2daSStefano Zampini   if (sf->nroots >= 0) { /* Reset only if graph already set */
4409566063dSJacob Faibussowitsch     PetscCall(PetscSFReset(sf));
4412a67d2daSStefano Zampini   }
4422a67d2daSStefano Zampini 
4439566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_SetGraph,sf,0,0,0));
44429046d53SLisandro Dalcin 
44595fce210SBarry Smith   sf->nroots  = nroots;
44695fce210SBarry Smith   sf->nleaves = nleaves;
44729046d53SLisandro Dalcin 
448db2b9530SVaclav Hapla   if (localmode == PETSC_COPY_VALUES && ilocal) {
449db2b9530SVaclav Hapla     PetscInt *tlocal = NULL;
450db2b9530SVaclav Hapla 
4519566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves,&tlocal));
4529566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tlocal,ilocal,nleaves));
453db2b9530SVaclav Hapla     ilocal = tlocal;
454db2b9530SVaclav Hapla   }
455db2b9530SVaclav Hapla   if (remotemode == PETSC_COPY_VALUES) {
456db2b9530SVaclav Hapla     PetscSFNode *tremote = NULL;
457db2b9530SVaclav Hapla 
4589566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves,&tremote));
4599566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(tremote,iremote,nleaves));
460db2b9530SVaclav Hapla     iremote = tremote;
461db2b9530SVaclav Hapla   }
462db2b9530SVaclav Hapla 
46329046d53SLisandro Dalcin   if (nleaves && ilocal) {
464db2b9530SVaclav Hapla     PetscSFNode   work;
465db2b9530SVaclav Hapla 
4669566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithDataArray(nleaves, ilocal, iremote, sizeof(PetscSFNode), &work));
4679566063dSJacob Faibussowitsch     PetscCall(PetscSortedCheckDupsInt(nleaves, ilocal, &unique));
468db2b9530SVaclav Hapla     unique = PetscNot(unique);
469db2b9530SVaclav 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");
470db2b9530SVaclav Hapla     sf->minleaf = ilocal[0];
471db2b9530SVaclav Hapla     sf->maxleaf = ilocal[nleaves-1];
472db2b9530SVaclav Hapla     contiguous = (PetscBool) (unique && ilocal[0] == 0 && ilocal[nleaves-1] == nleaves-1);
47329046d53SLisandro Dalcin   } else {
47429046d53SLisandro Dalcin     sf->minleaf = 0;
47529046d53SLisandro Dalcin     sf->maxleaf = nleaves - 1;
476db2b9530SVaclav Hapla     unique      = PETSC_TRUE;
477db2b9530SVaclav Hapla     contiguous  = PETSC_TRUE;
47829046d53SLisandro Dalcin   }
47929046d53SLisandro Dalcin 
480db2b9530SVaclav Hapla   if (contiguous) {
481db2b9530SVaclav Hapla     if (localmode == PETSC_USE_POINTER) {
482db2b9530SVaclav Hapla       ilocal = NULL;
483db2b9530SVaclav Hapla     } else {
4849566063dSJacob Faibussowitsch       PetscCall(PetscFree(ilocal));
485db2b9530SVaclav Hapla     }
486db2b9530SVaclav Hapla   }
487db2b9530SVaclav Hapla   sf->mine            = ilocal;
488db2b9530SVaclav Hapla   if (localmode == PETSC_USE_POINTER) {
48929046d53SLisandro Dalcin     sf->mine_alloc    = NULL;
490db2b9530SVaclav Hapla   } else {
491db2b9530SVaclav Hapla     sf->mine_alloc    = ilocal;
49295fce210SBarry Smith   }
493db2b9530SVaclav Hapla   sf->remote          = iremote;
494db2b9530SVaclav Hapla   if (remotemode == PETSC_USE_POINTER) {
49529046d53SLisandro Dalcin     sf->remote_alloc  = NULL;
496db2b9530SVaclav Hapla   } else {
497db2b9530SVaclav Hapla     sf->remote_alloc  = iremote;
49895fce210SBarry Smith   }
4999566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_SetGraph,sf,0,0,0));
50029046d53SLisandro Dalcin   sf->graphset = PETSC_TRUE;
50195fce210SBarry Smith   PetscFunctionReturn(0);
50295fce210SBarry Smith }
50395fce210SBarry Smith 
50429046d53SLisandro Dalcin /*@
505dd5b3ca6SJunchao Zhang   PetscSFSetGraphWithPattern - Sets the graph of an SF with a specific pattern
506dd5b3ca6SJunchao Zhang 
507dd5b3ca6SJunchao Zhang   Collective
508dd5b3ca6SJunchao Zhang 
509dd5b3ca6SJunchao Zhang   Input Parameters:
510dd5b3ca6SJunchao Zhang + sf      - The PetscSF
511dd5b3ca6SJunchao Zhang . map     - Layout of roots over all processes (insignificant when pattern is PETSCSF_PATTERN_ALLTOALL)
512dd5b3ca6SJunchao Zhang - pattern - One of PETSCSF_PATTERN_ALLGATHER, PETSCSF_PATTERN_GATHER, PETSCSF_PATTERN_ALLTOALL
513dd5b3ca6SJunchao Zhang 
514dd5b3ca6SJunchao Zhang   Notes:
515dd5b3ca6SJunchao Zhang   It is easier to explain PetscSFPattern using vectors. Suppose we have an MPI vector x and its layout is map.
516dd5b3ca6SJunchao Zhang   n and N are local and global sizes of x respectively.
517dd5b3ca6SJunchao Zhang 
518dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLGATHER, the routine creates a graph that if one does Bcast on it, it will copy x to
519dd5b3ca6SJunchao Zhang   sequential vectors y on all ranks.
520dd5b3ca6SJunchao Zhang 
521dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_GATHER, the routine creates a graph that if one does Bcast on it, it will copy x to a
522dd5b3ca6SJunchao Zhang   sequential vector y on rank 0.
523dd5b3ca6SJunchao Zhang 
524dd5b3ca6SJunchao Zhang   In above cases, entries of x are roots and entries of y are leaves.
525dd5b3ca6SJunchao Zhang 
526dd5b3ca6SJunchao Zhang   With PETSCSF_PATTERN_ALLTOALL, map is insignificant. Suppose NP is size of sf's communicator. The routine
527dd5b3ca6SJunchao Zhang   creates a graph that every rank has NP leaves and NP roots. On rank i, its leaf j is connected to root i
528dd5b3ca6SJunchao Zhang   of rank j. Here 0 <=i,j<NP. It is a kind of MPI_Alltoall with sendcount/recvcount being 1. Note that it does
529dd5b3ca6SJunchao Zhang   not mean one can not send multiple items. One just needs to create a new MPI datatype for the mulptiple data
530dd5b3ca6SJunchao Zhang   items with MPI_Type_contiguous() and use that as the <unit> argument in SF routines.
531dd5b3ca6SJunchao Zhang 
532dd5b3ca6SJunchao Zhang   In this case, roots and leaves are symmetric.
533dd5b3ca6SJunchao Zhang 
534dd5b3ca6SJunchao Zhang   Level: intermediate
535dd5b3ca6SJunchao Zhang  @*/
536dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFSetGraphWithPattern(PetscSF sf,PetscLayout map,PetscSFPattern pattern)
537dd5b3ca6SJunchao Zhang {
538dd5b3ca6SJunchao Zhang   MPI_Comm       comm;
539dd5b3ca6SJunchao Zhang   PetscInt       n,N,res[2];
540dd5b3ca6SJunchao Zhang   PetscMPIInt    rank,size;
541dd5b3ca6SJunchao Zhang   PetscSFType    type;
542dd5b3ca6SJunchao Zhang 
543dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
5442abc8c78SJacob Faibussowitsch   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
5452abc8c78SJacob Faibussowitsch   if (pattern != PETSCSF_PATTERN_ALLTOALL) PetscValidPointer(map,2);
5469566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
5472c71b3e2SJacob Faibussowitsch   PetscCheck(pattern >= PETSCSF_PATTERN_ALLGATHER && pattern <= PETSCSF_PATTERN_ALLTOALL,comm,PETSC_ERR_ARG_OUTOFRANGE,"Unsupported PetscSFPattern %d",pattern);
5489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
5499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
550dd5b3ca6SJunchao Zhang 
551dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLTOALL) {
552dd5b3ca6SJunchao Zhang     type = PETSCSFALLTOALL;
5539566063dSJacob Faibussowitsch     PetscCall(PetscLayoutCreate(comm,&sf->map));
5549566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetLocalSize(sf->map,size));
5559566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetSize(sf->map,((PetscInt)size)*size));
5569566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetUp(sf->map));
557dd5b3ca6SJunchao Zhang   } else {
5589566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetLocalSize(map,&n));
5599566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetSize(map,&N));
560dd5b3ca6SJunchao Zhang     res[0] = n;
561dd5b3ca6SJunchao Zhang     res[1] = -n;
562dd5b3ca6SJunchao Zhang     /* Check if n are same over all ranks so that we can optimize it */
5631c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(MPI_IN_PLACE,res,2,MPIU_INT,MPI_MAX,comm));
564dd5b3ca6SJunchao Zhang     if (res[0] == -res[1]) { /* same n */
565dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHER  : PETSCSFGATHER;
566dd5b3ca6SJunchao Zhang     } else {
567dd5b3ca6SJunchao Zhang       type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHERV : PETSCSFGATHERV;
568dd5b3ca6SJunchao Zhang     }
5699566063dSJacob Faibussowitsch     PetscCall(PetscLayoutReference(map,&sf->map));
570dd5b3ca6SJunchao Zhang   }
5719566063dSJacob Faibussowitsch   PetscCall(PetscSFSetType(sf,type));
572dd5b3ca6SJunchao Zhang 
573dd5b3ca6SJunchao Zhang   sf->pattern = pattern;
574dd5b3ca6SJunchao Zhang   sf->mine    = NULL; /* Contiguous */
575dd5b3ca6SJunchao Zhang 
576dd5b3ca6SJunchao Zhang   /* Set nleaves, nroots here in case user calls PetscSFGetGraph, which is legal to call even before PetscSFSetUp is called.
577dd5b3ca6SJunchao Zhang      Also set other easy stuff.
578dd5b3ca6SJunchao Zhang    */
579dd5b3ca6SJunchao Zhang   if (pattern == PETSCSF_PATTERN_ALLGATHER) {
580dd5b3ca6SJunchao Zhang     sf->nleaves      = N;
581dd5b3ca6SJunchao Zhang     sf->nroots       = n;
582dd5b3ca6SJunchao Zhang     sf->nranks       = size;
583dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
584dd5b3ca6SJunchao Zhang     sf->maxleaf      = N - 1;
585dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_GATHER) {
586dd5b3ca6SJunchao Zhang     sf->nleaves      = rank ? 0 : N;
587dd5b3ca6SJunchao Zhang     sf->nroots       = n;
588dd5b3ca6SJunchao Zhang     sf->nranks       = rank ? 0 : size;
589dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
590dd5b3ca6SJunchao Zhang     sf->maxleaf      = rank ? -1 : N - 1;
591dd5b3ca6SJunchao Zhang   } else if (pattern == PETSCSF_PATTERN_ALLTOALL) {
592dd5b3ca6SJunchao Zhang     sf->nleaves      = size;
593dd5b3ca6SJunchao Zhang     sf->nroots       = size;
594dd5b3ca6SJunchao Zhang     sf->nranks       = size;
595dd5b3ca6SJunchao Zhang     sf->minleaf      = 0;
596dd5b3ca6SJunchao Zhang     sf->maxleaf      = size - 1;
597dd5b3ca6SJunchao Zhang   }
598dd5b3ca6SJunchao Zhang   sf->ndranks  = 0; /* We do not need to separate out distinguished ranks for patterned graphs to improve communication performance */
599dd5b3ca6SJunchao Zhang   sf->graphset = PETSC_TRUE;
600dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
601dd5b3ca6SJunchao Zhang }
602dd5b3ca6SJunchao Zhang 
603dd5b3ca6SJunchao Zhang /*@
60495fce210SBarry Smith    PetscSFCreateInverseSF - given a PetscSF in which all vertices have degree 1, creates the inverse map
60595fce210SBarry Smith 
60695fce210SBarry Smith    Collective
60795fce210SBarry Smith 
6084165533cSJose E. Roman    Input Parameter:
60995fce210SBarry Smith .  sf - star forest to invert
61095fce210SBarry Smith 
6114165533cSJose E. Roman    Output Parameter:
61295fce210SBarry Smith .  isf - inverse of sf
6134165533cSJose E. Roman 
61495fce210SBarry Smith    Level: advanced
61595fce210SBarry Smith 
61695fce210SBarry Smith    Notes:
61795fce210SBarry Smith    All roots must have degree 1.
61895fce210SBarry Smith 
61995fce210SBarry Smith    The local space may be a permutation, but cannot be sparse.
62095fce210SBarry Smith 
621db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`
62295fce210SBarry Smith @*/
62395fce210SBarry Smith PetscErrorCode PetscSFCreateInverseSF(PetscSF sf,PetscSF *isf)
62495fce210SBarry Smith {
62595fce210SBarry Smith   PetscMPIInt    rank;
62695fce210SBarry Smith   PetscInt       i,nroots,nleaves,maxlocal,count,*newilocal;
62795fce210SBarry Smith   const PetscInt *ilocal;
62895fce210SBarry Smith   PetscSFNode    *roots,*leaves;
62995fce210SBarry Smith 
63095fce210SBarry Smith   PetscFunctionBegin;
63129046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
63229046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
63329046d53SLisandro Dalcin   PetscValidPointer(isf,2);
63429046d53SLisandro Dalcin 
6359566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL));
63629046d53SLisandro Dalcin   maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
63729046d53SLisandro Dalcin 
6389566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank));
6399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nroots,&roots,maxlocal,&leaves));
640ae9aee6dSMatthew G. Knepley   for (i=0; i<maxlocal; i++) {
64195fce210SBarry Smith     leaves[i].rank  = rank;
64295fce210SBarry Smith     leaves[i].index = i;
64395fce210SBarry Smith   }
64495fce210SBarry Smith   for (i=0; i <nroots; i++) {
64595fce210SBarry Smith     roots[i].rank  = -1;
64695fce210SBarry Smith     roots[i].index = -1;
64795fce210SBarry Smith   }
6489566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf,MPIU_2INT,leaves,roots,MPI_REPLACE));
6499566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sf,MPIU_2INT,leaves,roots,MPI_REPLACE));
65095fce210SBarry Smith 
65195fce210SBarry Smith   /* Check whether our leaves are sparse */
65295fce210SBarry Smith   for (i=0,count=0; i<nroots; i++) if (roots[i].rank >= 0) count++;
65395fce210SBarry Smith   if (count == nroots) newilocal = NULL;
65495fce210SBarry Smith   else {                        /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */
6559566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count,&newilocal));
65695fce210SBarry Smith     for (i=0,count=0; i<nroots; i++) {
65795fce210SBarry Smith       if (roots[i].rank >= 0) {
65895fce210SBarry Smith         newilocal[count]   = i;
65995fce210SBarry Smith         roots[count].rank  = roots[i].rank;
66095fce210SBarry Smith         roots[count].index = roots[i].index;
66195fce210SBarry Smith         count++;
66295fce210SBarry Smith       }
66395fce210SBarry Smith     }
66495fce210SBarry Smith   }
66595fce210SBarry Smith 
6669566063dSJacob Faibussowitsch   PetscCall(PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,isf));
6679566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*isf,maxlocal,count,newilocal,PETSC_OWN_POINTER,roots,PETSC_COPY_VALUES));
6689566063dSJacob Faibussowitsch   PetscCall(PetscFree2(roots,leaves));
66995fce210SBarry Smith   PetscFunctionReturn(0);
67095fce210SBarry Smith }
67195fce210SBarry Smith 
67295fce210SBarry Smith /*@
67395fce210SBarry Smith    PetscSFDuplicate - duplicate a PetscSF, optionally preserving rank connectivity and graph
67495fce210SBarry Smith 
67595fce210SBarry Smith    Collective
67695fce210SBarry Smith 
6774165533cSJose E. Roman    Input Parameters:
67895fce210SBarry Smith +  sf - communication object to duplicate
67995fce210SBarry Smith -  opt - PETSCSF_DUPLICATE_CONFONLY, PETSCSF_DUPLICATE_RANKS, or PETSCSF_DUPLICATE_GRAPH (see PetscSFDuplicateOption)
68095fce210SBarry Smith 
6814165533cSJose E. Roman    Output Parameter:
68295fce210SBarry Smith .  newsf - new communication object
68395fce210SBarry Smith 
68495fce210SBarry Smith    Level: beginner
68595fce210SBarry Smith 
686db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFSetType()`, `PetscSFSetGraph()`
68795fce210SBarry Smith @*/
68895fce210SBarry Smith PetscErrorCode PetscSFDuplicate(PetscSF sf,PetscSFDuplicateOption opt,PetscSF *newsf)
68995fce210SBarry Smith {
69029046d53SLisandro Dalcin   PetscSFType    type;
69197929ea7SJunchao Zhang   MPI_Datatype   dtype=MPIU_SCALAR;
69295fce210SBarry Smith 
69395fce210SBarry Smith   PetscFunctionBegin;
69429046d53SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
69529046d53SLisandro Dalcin   PetscValidLogicalCollectiveEnum(sf,opt,2);
69629046d53SLisandro Dalcin   PetscValidPointer(newsf,3);
6979566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sf),newsf));
6989566063dSJacob Faibussowitsch   PetscCall(PetscSFGetType(sf,&type));
6999566063dSJacob Faibussowitsch   if (type) PetscCall(PetscSFSetType(*newsf,type));
70054a24607SJunchao Zhang   (*newsf)->allow_multi_leaves = sf->allow_multi_leaves; /* Dup this flag ealier since PetscSFSetGraph() below checks on this flag */
70195fce210SBarry Smith   if (opt == PETSCSF_DUPLICATE_GRAPH) {
702dd5b3ca6SJunchao Zhang     PetscSFCheckGraphSet(sf,1);
703dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
70495fce210SBarry Smith       PetscInt          nroots,nleaves;
70595fce210SBarry Smith       const PetscInt    *ilocal;
70695fce210SBarry Smith       const PetscSFNode *iremote;
7079566063dSJacob Faibussowitsch       PetscCall(PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote));
7089566063dSJacob Faibussowitsch       PetscCall(PetscSFSetGraph(*newsf,nroots,nleaves,(PetscInt*)ilocal,PETSC_COPY_VALUES,(PetscSFNode*)iremote,PETSC_COPY_VALUES));
709dd5b3ca6SJunchao Zhang     } else {
7109566063dSJacob Faibussowitsch       PetscCall(PetscSFSetGraphWithPattern(*newsf,sf->map,sf->pattern));
711dd5b3ca6SJunchao Zhang     }
71295fce210SBarry Smith   }
71397929ea7SJunchao Zhang   /* Since oldtype is committed, so is newtype, according to MPI */
7149566063dSJacob Faibussowitsch   if (sf->vscat.bs > 1) PetscCallMPI(MPI_Type_dup(sf->vscat.unit,&dtype));
71597929ea7SJunchao Zhang   (*newsf)->vscat.bs     = sf->vscat.bs;
71697929ea7SJunchao Zhang   (*newsf)->vscat.unit   = dtype;
71797929ea7SJunchao Zhang   (*newsf)->vscat.to_n   = sf->vscat.to_n;
71897929ea7SJunchao Zhang   (*newsf)->vscat.from_n = sf->vscat.from_n;
71997929ea7SJunchao Zhang   /* Do not copy lsf. Build it on demand since it is rarely used */
72097929ea7SJunchao Zhang 
72120c24465SJunchao Zhang #if defined(PETSC_HAVE_DEVICE)
72220c24465SJunchao Zhang   (*newsf)->backend              = sf->backend;
72371438e86SJunchao Zhang   (*newsf)->unknown_input_stream= sf->unknown_input_stream;
72420c24465SJunchao Zhang   (*newsf)->use_gpu_aware_mpi    = sf->use_gpu_aware_mpi;
72520c24465SJunchao Zhang   (*newsf)->use_stream_aware_mpi = sf->use_stream_aware_mpi;
72620c24465SJunchao Zhang #endif
727*dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf,Duplicate,opt,*newsf);
72820c24465SJunchao Zhang   /* Don't do PetscSFSetUp() since the new sf's graph might have not been set. */
72995fce210SBarry Smith   PetscFunctionReturn(0);
73095fce210SBarry Smith }
73195fce210SBarry Smith 
73295fce210SBarry Smith /*@C
73395fce210SBarry Smith    PetscSFGetGraph - Get the graph specifying a parallel star forest
73495fce210SBarry Smith 
73595fce210SBarry Smith    Not Collective
73695fce210SBarry Smith 
7374165533cSJose E. Roman    Input Parameter:
73895fce210SBarry Smith .  sf - star forest
73995fce210SBarry Smith 
7404165533cSJose E. Roman    Output Parameters:
74195fce210SBarry Smith +  nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves)
74295fce210SBarry Smith .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
743bc6585dcSJunchao Zhang .  ilocal - locations of leaves in leafdata buffers (if returned value is NULL, it means leaves are in contiguous storage)
74495fce210SBarry Smith -  iremote - remote locations of root vertices for each leaf on the current process
74595fce210SBarry Smith 
746373e0d91SLisandro Dalcin    Notes:
747373e0d91SLisandro Dalcin      We are not currently requiring that the graph is set, thus returning nroots=-1 if it has not been set yet
748373e0d91SLisandro Dalcin 
749db2b9530SVaclav Hapla      The returned ilocal/iremote might contain values in different order than the input ones in PetscSFSetGraph(),
750db2b9530SVaclav Hapla      see its manpage for details.
751db2b9530SVaclav Hapla 
7528dbb0df6SBarry Smith    Fortran Notes:
7538dbb0df6SBarry Smith      The returned iremote array is a copy and must be deallocated after use. Consequently, if you
7548dbb0df6SBarry Smith      want to update the graph, you must call PetscSFSetGraph() after modifying the iremote array.
7558dbb0df6SBarry Smith 
7568dbb0df6SBarry Smith      To check for a NULL ilocal use
7578dbb0df6SBarry Smith $      if (loc(ilocal) == loc(PETSC_NULL_INTEGER)) then
758ca797d7aSLawrence Mitchell 
75995fce210SBarry Smith    Level: intermediate
76095fce210SBarry Smith 
761db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`
76295fce210SBarry Smith @*/
76395fce210SBarry Smith PetscErrorCode PetscSFGetGraph(PetscSF sf,PetscInt *nroots,PetscInt *nleaves,const PetscInt **ilocal,const PetscSFNode **iremote)
76495fce210SBarry Smith {
76595fce210SBarry Smith   PetscFunctionBegin;
76695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
767b8dee149SJunchao Zhang   if (sf->ops->GetGraph) {
7689566063dSJacob Faibussowitsch     PetscCall((sf->ops->GetGraph)(sf,nroots,nleaves,ilocal,iremote));
769b8dee149SJunchao Zhang   } else {
77095fce210SBarry Smith     if (nroots) *nroots = sf->nroots;
77195fce210SBarry Smith     if (nleaves) *nleaves = sf->nleaves;
77295fce210SBarry Smith     if (ilocal) *ilocal = sf->mine;
77395fce210SBarry Smith     if (iremote) *iremote = sf->remote;
774b8dee149SJunchao Zhang   }
77595fce210SBarry Smith   PetscFunctionReturn(0);
77695fce210SBarry Smith }
77795fce210SBarry Smith 
77829046d53SLisandro Dalcin /*@
77995fce210SBarry Smith    PetscSFGetLeafRange - Get the active leaf ranges
78095fce210SBarry Smith 
78195fce210SBarry Smith    Not Collective
78295fce210SBarry Smith 
7834165533cSJose E. Roman    Input Parameter:
78495fce210SBarry Smith .  sf - star forest
78595fce210SBarry Smith 
7864165533cSJose E. Roman    Output Parameters:
787dd5b3ca6SJunchao Zhang +  minleaf - minimum active leaf on this process. Return 0 if there are no leaves.
788dd5b3ca6SJunchao Zhang -  maxleaf - maximum active leaf on this process. Return -1 if there are no leaves.
78995fce210SBarry Smith 
79095fce210SBarry Smith    Level: developer
79195fce210SBarry Smith 
792db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
79395fce210SBarry Smith @*/
79495fce210SBarry Smith PetscErrorCode PetscSFGetLeafRange(PetscSF sf,PetscInt *minleaf,PetscInt *maxleaf)
79595fce210SBarry Smith {
79695fce210SBarry Smith   PetscFunctionBegin;
79795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
79829046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
79995fce210SBarry Smith   if (minleaf) *minleaf = sf->minleaf;
80095fce210SBarry Smith   if (maxleaf) *maxleaf = sf->maxleaf;
80195fce210SBarry Smith   PetscFunctionReturn(0);
80295fce210SBarry Smith }
80395fce210SBarry Smith 
80495fce210SBarry Smith /*@C
805fe2efc57SMark    PetscSFViewFromOptions - View from Options
806fe2efc57SMark 
807fe2efc57SMark    Collective on PetscSF
808fe2efc57SMark 
809fe2efc57SMark    Input Parameters:
810fe2efc57SMark +  A - the star forest
811736c3998SJose E. Roman .  obj - Optional object
812736c3998SJose E. Roman -  name - command line option
813fe2efc57SMark 
814fe2efc57SMark    Level: intermediate
815db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFView`, `PetscObjectViewFromOptions()`, `PetscSFCreate()`
816fe2efc57SMark @*/
817fe2efc57SMark PetscErrorCode  PetscSFViewFromOptions(PetscSF A,PetscObject obj,const char name[])
818fe2efc57SMark {
819fe2efc57SMark   PetscFunctionBegin;
820fe2efc57SMark   PetscValidHeaderSpecific(A,PETSCSF_CLASSID,1);
8219566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A,obj,name));
822fe2efc57SMark   PetscFunctionReturn(0);
823fe2efc57SMark }
824fe2efc57SMark 
825fe2efc57SMark /*@C
82695fce210SBarry Smith    PetscSFView - view a star forest
82795fce210SBarry Smith 
82895fce210SBarry Smith    Collective
82995fce210SBarry Smith 
8304165533cSJose E. Roman    Input Parameters:
83195fce210SBarry Smith +  sf - star forest
83295fce210SBarry Smith -  viewer - viewer to display graph, for example PETSC_VIEWER_STDOUT_WORLD
83395fce210SBarry Smith 
83495fce210SBarry Smith    Level: beginner
83595fce210SBarry Smith 
836db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `PetscSFSetGraph()`
83795fce210SBarry Smith @*/
83895fce210SBarry Smith PetscErrorCode PetscSFView(PetscSF sf,PetscViewer viewer)
83995fce210SBarry Smith {
84095fce210SBarry Smith   PetscBool         iascii;
84195fce210SBarry Smith   PetscViewerFormat format;
84295fce210SBarry Smith 
84395fce210SBarry Smith   PetscFunctionBegin;
84495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
8459566063dSJacob Faibussowitsch   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf),&viewer));
84695fce210SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
84795fce210SBarry Smith   PetscCheckSameComm(sf,1,viewer,2);
8489566063dSJacob Faibussowitsch   if (sf->graphset) PetscCall(PetscSFSetUp(sf));
8499566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
85053dd6d7dSJunchao Zhang   if (iascii && viewer->format != PETSC_VIEWER_ASCII_MATLAB) {
85195fce210SBarry Smith     PetscMPIInt rank;
85281bfa7aaSJed Brown     PetscInt    ii,i,j;
85395fce210SBarry Smith 
8549566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)sf,viewer));
8559566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
856dd5b3ca6SJunchao Zhang     if (sf->pattern == PETSCSF_PATTERN_GENERAL) {
85780153354SVaclav Hapla       if (!sf->graphset) {
8589566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"PetscSFSetGraph() has not been called yet\n"));
8599566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
86080153354SVaclav Hapla         PetscFunctionReturn(0);
86180153354SVaclav Hapla       }
8629566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank));
8639566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
8649566063dSJacob 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));
86595fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
8669566063dSJacob Faibussowitsch         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));
86795fce210SBarry Smith       }
8689566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
8699566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetFormat(viewer,&format));
87095fce210SBarry Smith       if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
87181bfa7aaSJed Brown         PetscMPIInt *tmpranks,*perm;
8729566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(sf->nranks,&tmpranks,sf->nranks,&perm));
8739566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(tmpranks,sf->ranks,sf->nranks));
87481bfa7aaSJed Brown         for (i=0; i<sf->nranks; i++) perm[i] = i;
8759566063dSJacob Faibussowitsch         PetscCall(PetscSortMPIIntWithArray(sf->nranks,tmpranks,perm));
8769566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Roots referenced by my leaves, by rank\n",rank));
87781bfa7aaSJed Brown         for (ii=0; ii<sf->nranks; ii++) {
87881bfa7aaSJed Brown           i = perm[ii];
8799566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d: %" PetscInt_FMT " edges\n",rank,sf->ranks[i],sf->roffset[i+1]-sf->roffset[i]));
88095fce210SBarry Smith           for (j=sf->roffset[i]; j<sf->roffset[i+1]; j++) {
8819566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"[%d]    %" PetscInt_FMT " <- %" PetscInt_FMT "\n",rank,sf->rmine[j],sf->rremote[j]));
88295fce210SBarry Smith           }
88395fce210SBarry Smith         }
8849566063dSJacob Faibussowitsch         PetscCall(PetscFree2(tmpranks,perm));
88595fce210SBarry Smith       }
8869566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
8879566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
888dd5b3ca6SJunchao Zhang     }
8899566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
89095fce210SBarry Smith   }
891*dbbe0bcdSBarry Smith   PetscTryTypeMethod(sf,View,viewer);
89295fce210SBarry Smith   PetscFunctionReturn(0);
89395fce210SBarry Smith }
89495fce210SBarry Smith 
89595fce210SBarry Smith /*@C
896dec1416fSJunchao Zhang    PetscSFGetRootRanks - Get root ranks and number of vertices referenced by leaves on this process
89795fce210SBarry Smith 
89895fce210SBarry Smith    Not Collective
89995fce210SBarry Smith 
9004165533cSJose E. Roman    Input Parameter:
90195fce210SBarry Smith .  sf - star forest
90295fce210SBarry Smith 
9034165533cSJose E. Roman    Output Parameters:
90495fce210SBarry Smith +  nranks - number of ranks referenced by local part
90595fce210SBarry Smith .  ranks - array of ranks
90695fce210SBarry Smith .  roffset - offset in rmine/rremote for each rank (length nranks+1)
90795fce210SBarry Smith .  rmine - concatenated array holding local indices referencing each remote rank
90895fce210SBarry Smith -  rremote - concatenated array holding remote indices referenced for each remote rank
90995fce210SBarry Smith 
91095fce210SBarry Smith    Level: developer
91195fce210SBarry Smith 
912db781477SPatrick Sanan .seealso: `PetscSFGetLeafRanks()`
91395fce210SBarry Smith @*/
914dec1416fSJunchao Zhang PetscErrorCode PetscSFGetRootRanks(PetscSF sf,PetscInt *nranks,const PetscMPIInt **ranks,const PetscInt **roffset,const PetscInt **rmine,const PetscInt **rremote)
91595fce210SBarry Smith {
91695fce210SBarry Smith   PetscFunctionBegin;
91795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
91828b400f6SJacob Faibussowitsch   PetscCheck(sf->setupcalled,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
919dec1416fSJunchao Zhang   if (sf->ops->GetRootRanks) {
9209566063dSJacob Faibussowitsch     PetscCall((sf->ops->GetRootRanks)(sf,nranks,ranks,roffset,rmine,rremote));
921dec1416fSJunchao Zhang   } else {
922dec1416fSJunchao Zhang     /* The generic implementation */
92395fce210SBarry Smith     if (nranks)  *nranks  = sf->nranks;
92495fce210SBarry Smith     if (ranks)   *ranks   = sf->ranks;
92595fce210SBarry Smith     if (roffset) *roffset = sf->roffset;
92695fce210SBarry Smith     if (rmine)   *rmine   = sf->rmine;
92795fce210SBarry Smith     if (rremote) *rremote = sf->rremote;
928dec1416fSJunchao Zhang   }
92995fce210SBarry Smith   PetscFunctionReturn(0);
93095fce210SBarry Smith }
93195fce210SBarry Smith 
9328750ddebSJunchao Zhang /*@C
9338750ddebSJunchao Zhang    PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process
9348750ddebSJunchao Zhang 
9358750ddebSJunchao Zhang    Not Collective
9368750ddebSJunchao Zhang 
9374165533cSJose E. Roman    Input Parameter:
9388750ddebSJunchao Zhang .  sf - star forest
9398750ddebSJunchao Zhang 
9404165533cSJose E. Roman    Output Parameters:
9418750ddebSJunchao Zhang +  niranks - number of leaf ranks referencing roots on this process
9428750ddebSJunchao Zhang .  iranks - array of ranks
9438750ddebSJunchao Zhang .  ioffset - offset in irootloc for each rank (length niranks+1)
9448750ddebSJunchao Zhang -  irootloc - concatenated array holding local indices of roots referenced by each leaf rank
9458750ddebSJunchao Zhang 
9468750ddebSJunchao Zhang    Level: developer
9478750ddebSJunchao Zhang 
948db781477SPatrick Sanan .seealso: `PetscSFGetRootRanks()`
9498750ddebSJunchao Zhang @*/
9508750ddebSJunchao Zhang PetscErrorCode PetscSFGetLeafRanks(PetscSF sf,PetscInt *niranks,const PetscMPIInt **iranks,const PetscInt **ioffset,const PetscInt **irootloc)
9518750ddebSJunchao Zhang {
9528750ddebSJunchao Zhang   PetscFunctionBegin;
9538750ddebSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
95428b400f6SJacob Faibussowitsch   PetscCheck(sf->setupcalled,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUp() before obtaining ranks");
9558750ddebSJunchao Zhang   if (sf->ops->GetLeafRanks) {
9569566063dSJacob Faibussowitsch     PetscCall((sf->ops->GetLeafRanks)(sf,niranks,iranks,ioffset,irootloc));
9578750ddebSJunchao Zhang   } else {
9588750ddebSJunchao Zhang     PetscSFType type;
9599566063dSJacob Faibussowitsch     PetscCall(PetscSFGetType(sf,&type));
96098921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type);
9618750ddebSJunchao Zhang   }
9628750ddebSJunchao Zhang   PetscFunctionReturn(0);
9638750ddebSJunchao Zhang }
9648750ddebSJunchao Zhang 
965b5a8e515SJed Brown static PetscBool InList(PetscMPIInt needle,PetscMPIInt n,const PetscMPIInt *list) {
966b5a8e515SJed Brown   PetscInt i;
967b5a8e515SJed Brown   for (i=0; i<n; i++) {
968b5a8e515SJed Brown     if (needle == list[i]) return PETSC_TRUE;
969b5a8e515SJed Brown   }
970b5a8e515SJed Brown   return PETSC_FALSE;
971b5a8e515SJed Brown }
972b5a8e515SJed Brown 
97395fce210SBarry Smith /*@C
97421c688dcSJed Brown    PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by PetscSF implementations.
97521c688dcSJed Brown 
97621c688dcSJed Brown    Collective
97721c688dcSJed Brown 
9784165533cSJose E. Roman    Input Parameters:
979b5a8e515SJed Brown +  sf - PetscSF to set up; PetscSFSetGraph() must have been called
980b5a8e515SJed Brown -  dgroup - MPI_Group of ranks to be distinguished (e.g., for self or shared memory exchange)
98121c688dcSJed Brown 
98221c688dcSJed Brown    Level: developer
98321c688dcSJed Brown 
984db781477SPatrick Sanan .seealso: `PetscSFGetRootRanks()`
98521c688dcSJed Brown @*/
986b5a8e515SJed Brown PetscErrorCode PetscSFSetUpRanks(PetscSF sf,MPI_Group dgroup)
98721c688dcSJed Brown {
98821c688dcSJed Brown   PetscTable         table;
98921c688dcSJed Brown   PetscTablePosition pos;
990b5a8e515SJed Brown   PetscMPIInt        size,groupsize,*groupranks;
991247e8311SStefano Zampini   PetscInt           *rcount,*ranks;
992247e8311SStefano Zampini   PetscInt           i, irank = -1,orank = -1;
99321c688dcSJed Brown 
99421c688dcSJed Brown   PetscFunctionBegin;
99521c688dcSJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
99629046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
9979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf),&size));
9989566063dSJacob Faibussowitsch   PetscCall(PetscTableCreate(10,size,&table));
99921c688dcSJed Brown   for (i=0; i<sf->nleaves; i++) {
100021c688dcSJed Brown     /* Log 1-based rank */
10019566063dSJacob Faibussowitsch     PetscCall(PetscTableAdd(table,sf->remote[i].rank+1,1,ADD_VALUES));
100221c688dcSJed Brown   }
10039566063dSJacob Faibussowitsch   PetscCall(PetscTableGetCount(table,&sf->nranks));
10049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(sf->nranks,&sf->ranks,sf->nranks+1,&sf->roffset,sf->nleaves,&sf->rmine,sf->nleaves,&sf->rremote));
10059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(sf->nranks,&rcount,sf->nranks,&ranks));
10069566063dSJacob Faibussowitsch   PetscCall(PetscTableGetHeadPosition(table,&pos));
100721c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
10089566063dSJacob Faibussowitsch     PetscCall(PetscTableGetNext(table,&pos,&ranks[i],&rcount[i]));
100921c688dcSJed Brown     ranks[i]--;             /* Convert back to 0-based */
101021c688dcSJed Brown   }
10119566063dSJacob Faibussowitsch   PetscCall(PetscTableDestroy(&table));
1012b5a8e515SJed Brown 
1013b5a8e515SJed Brown   /* We expect that dgroup is reliably "small" while nranks could be large */
1014b5a8e515SJed Brown   {
10157fb8a5e4SKarl Rupp     MPI_Group group = MPI_GROUP_NULL;
1016b5a8e515SJed Brown     PetscMPIInt *dgroupranks;
10179566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group));
10189566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_size(dgroup,&groupsize));
10199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(groupsize,&dgroupranks));
10209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(groupsize,&groupranks));
1021b5a8e515SJed Brown     for (i=0; i<groupsize; i++) dgroupranks[i] = i;
10229566063dSJacob Faibussowitsch     if (groupsize) PetscCallMPI(MPI_Group_translate_ranks(dgroup,groupsize,dgroupranks,group,groupranks));
10239566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_free(&group));
10249566063dSJacob Faibussowitsch     PetscCall(PetscFree(dgroupranks));
1025b5a8e515SJed Brown   }
1026b5a8e515SJed Brown 
1027b5a8e515SJed Brown   /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */
1028b5a8e515SJed Brown   for (sf->ndranks=0,i=sf->nranks; sf->ndranks<i;) {
1029b5a8e515SJed Brown     for (i--; sf->ndranks<i; i--) { /* Scan i backward looking for distinguished rank */
1030b5a8e515SJed Brown       if (InList(ranks[i],groupsize,groupranks)) break;
1031b5a8e515SJed Brown     }
1032b5a8e515SJed Brown     for (; sf->ndranks<=i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */
1033b5a8e515SJed Brown       if (!InList(ranks[sf->ndranks],groupsize,groupranks)) break;
1034b5a8e515SJed Brown     }
1035b5a8e515SJed Brown     if (sf->ndranks < i) {                         /* Swap ranks[sf->ndranks] with ranks[i] */
1036b5a8e515SJed Brown       PetscInt tmprank,tmpcount;
1037247e8311SStefano Zampini 
1038b5a8e515SJed Brown       tmprank             = ranks[i];
1039b5a8e515SJed Brown       tmpcount            = rcount[i];
1040b5a8e515SJed Brown       ranks[i]            = ranks[sf->ndranks];
1041b5a8e515SJed Brown       rcount[i]           = rcount[sf->ndranks];
1042b5a8e515SJed Brown       ranks[sf->ndranks]  = tmprank;
1043b5a8e515SJed Brown       rcount[sf->ndranks] = tmpcount;
1044b5a8e515SJed Brown       sf->ndranks++;
1045b5a8e515SJed Brown     }
1046b5a8e515SJed Brown   }
10479566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupranks));
10489566063dSJacob Faibussowitsch   PetscCall(PetscSortIntWithArray(sf->ndranks,ranks,rcount));
10499566063dSJacob Faibussowitsch   PetscCall(PetscSortIntWithArray(sf->nranks-sf->ndranks,ranks+sf->ndranks,rcount+sf->ndranks));
105021c688dcSJed Brown   sf->roffset[0] = 0;
105121c688dcSJed Brown   for (i=0; i<sf->nranks; i++) {
10529566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(ranks[i],sf->ranks+i));
105321c688dcSJed Brown     sf->roffset[i+1] = sf->roffset[i] + rcount[i];
105421c688dcSJed Brown     rcount[i]        = 0;
105521c688dcSJed Brown   }
1056247e8311SStefano Zampini   for (i=0, irank = -1, orank = -1; i<sf->nleaves; i++) {
1057247e8311SStefano Zampini     /* short circuit */
1058247e8311SStefano Zampini     if (orank != sf->remote[i].rank) {
105921c688dcSJed Brown       /* Search for index of iremote[i].rank in sf->ranks */
10609566063dSJacob Faibussowitsch       PetscCall(PetscFindMPIInt(sf->remote[i].rank,sf->ndranks,sf->ranks,&irank));
1061b5a8e515SJed Brown       if (irank < 0) {
10629566063dSJacob Faibussowitsch         PetscCall(PetscFindMPIInt(sf->remote[i].rank,sf->nranks-sf->ndranks,sf->ranks+sf->ndranks,&irank));
1063b5a8e515SJed Brown         if (irank >= 0) irank += sf->ndranks;
106421c688dcSJed Brown       }
1065247e8311SStefano Zampini       orank = sf->remote[i].rank;
1066247e8311SStefano Zampini     }
106708401ef6SPierre Jolivet     PetscCheck(irank >= 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Could not find rank %" PetscInt_FMT " in array",sf->remote[i].rank);
106821c688dcSJed Brown     sf->rmine[sf->roffset[irank] + rcount[irank]]   = sf->mine ? sf->mine[i] : i;
106921c688dcSJed Brown     sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index;
107021c688dcSJed Brown     rcount[irank]++;
107121c688dcSJed Brown   }
10729566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rcount,ranks));
107321c688dcSJed Brown   PetscFunctionReturn(0);
107421c688dcSJed Brown }
107521c688dcSJed Brown 
107621c688dcSJed Brown /*@C
107795fce210SBarry Smith    PetscSFGetGroups - gets incoming and outgoing process groups
107895fce210SBarry Smith 
107995fce210SBarry Smith    Collective
108095fce210SBarry Smith 
10814165533cSJose E. Roman    Input Parameter:
108295fce210SBarry Smith .  sf - star forest
108395fce210SBarry Smith 
10844165533cSJose E. Roman    Output Parameters:
108595fce210SBarry Smith +  incoming - group of origin processes for incoming edges (leaves that reference my roots)
108695fce210SBarry Smith -  outgoing - group of destination processes for outgoing edges (roots that I reference)
108795fce210SBarry Smith 
108895fce210SBarry Smith    Level: developer
108995fce210SBarry Smith 
1090db781477SPatrick Sanan .seealso: `PetscSFGetWindow()`, `PetscSFRestoreWindow()`
109195fce210SBarry Smith @*/
109295fce210SBarry Smith PetscErrorCode PetscSFGetGroups(PetscSF sf,MPI_Group *incoming,MPI_Group *outgoing)
109395fce210SBarry Smith {
10947fb8a5e4SKarl Rupp   MPI_Group      group = MPI_GROUP_NULL;
109595fce210SBarry Smith 
109695fce210SBarry Smith   PetscFunctionBegin;
109708401ef6SPierre Jolivet   PetscCheck(sf->nranks >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSFSetUpRanks() before obtaining groups");
109895fce210SBarry Smith   if (sf->ingroup == MPI_GROUP_NULL) {
109995fce210SBarry Smith     PetscInt       i;
110095fce210SBarry Smith     const PetscInt *indegree;
110195fce210SBarry Smith     PetscMPIInt    rank,*outranks,*inranks;
110295fce210SBarry Smith     PetscSFNode    *remote;
110395fce210SBarry Smith     PetscSF        bgcount;
110495fce210SBarry Smith 
110595fce210SBarry Smith     /* Compute the number of incoming ranks */
11069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(sf->nranks,&remote));
110795fce210SBarry Smith     for (i=0; i<sf->nranks; i++) {
110895fce210SBarry Smith       remote[i].rank  = sf->ranks[i];
110995fce210SBarry Smith       remote[i].index = 0;
111095fce210SBarry Smith     }
11119566063dSJacob Faibussowitsch     PetscCall(PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,&bgcount));
11129566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(bgcount,1,sf->nranks,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER));
11139566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(bgcount,&indegree));
11149566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(bgcount,&indegree));
111595fce210SBarry Smith     /* Enumerate the incoming ranks */
11169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(indegree[0],&inranks,sf->nranks,&outranks));
11179566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank));
111895fce210SBarry Smith     for (i=0; i<sf->nranks; i++) outranks[i] = rank;
11199566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherBegin(bgcount,MPI_INT,outranks,inranks));
11209566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherEnd(bgcount,MPI_INT,outranks,inranks));
11219566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group));
11229566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_incl(group,indegree[0],inranks,&sf->ingroup));
11239566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_free(&group));
11249566063dSJacob Faibussowitsch     PetscCall(PetscFree2(inranks,outranks));
11259566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&bgcount));
112695fce210SBarry Smith   }
112795fce210SBarry Smith   *incoming = sf->ingroup;
112895fce210SBarry Smith 
112995fce210SBarry Smith   if (sf->outgroup == MPI_GROUP_NULL) {
11309566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)sf),&group));
11319566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_incl(group,sf->nranks,sf->ranks,&sf->outgroup));
11329566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Group_free(&group));
113395fce210SBarry Smith   }
113495fce210SBarry Smith   *outgoing = sf->outgroup;
113595fce210SBarry Smith   PetscFunctionReturn(0);
113695fce210SBarry Smith }
113795fce210SBarry Smith 
113829046d53SLisandro Dalcin /*@
11393b8d980fSPierre Jolivet    PetscSFGetMultiSF - gets the inner SF implementing gathers and scatters
114095fce210SBarry Smith 
114195fce210SBarry Smith    Collective
114295fce210SBarry Smith 
11434165533cSJose E. Roman    Input Parameter:
114495fce210SBarry Smith .  sf - star forest that may contain roots with 0 or with more than 1 vertex
114595fce210SBarry Smith 
11464165533cSJose E. Roman    Output Parameter:
114795fce210SBarry Smith .  multi - star forest with split roots, such that each root has degree exactly 1
114895fce210SBarry Smith 
114995fce210SBarry Smith    Level: developer
115095fce210SBarry Smith 
115195fce210SBarry Smith    Notes:
115295fce210SBarry Smith 
115395fce210SBarry Smith    In most cases, users should use PetscSFGatherBegin() and PetscSFScatterBegin() instead of manipulating multi
115495fce210SBarry Smith    directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming
115595fce210SBarry Smith    edge, it is a candidate for future optimization that might involve its removal.
115695fce210SBarry Smith 
1157db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`, `PetscSFGatherBegin()`, `PetscSFScatterBegin()`, `PetscSFComputeMultiRootOriginalNumbering()`
115895fce210SBarry Smith @*/
115995fce210SBarry Smith PetscErrorCode PetscSFGetMultiSF(PetscSF sf,PetscSF *multi)
116095fce210SBarry Smith {
116195fce210SBarry Smith   PetscFunctionBegin;
116295fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
116395fce210SBarry Smith   PetscValidPointer(multi,2);
116495fce210SBarry Smith   if (sf->nroots < 0) {         /* Graph has not been set yet; why do we need this? */
11659566063dSJacob Faibussowitsch     PetscCall(PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi));
116695fce210SBarry Smith     *multi = sf->multi;
1167013b3241SStefano Zampini     sf->multi->multi = sf->multi;
116895fce210SBarry Smith     PetscFunctionReturn(0);
116995fce210SBarry Smith   }
117095fce210SBarry Smith   if (!sf->multi) {
117195fce210SBarry Smith     const PetscInt *indegree;
11729837ea96SMatthew G. Knepley     PetscInt       i,*inoffset,*outones,*outoffset,maxlocal;
117395fce210SBarry Smith     PetscSFNode    *remote;
117429046d53SLisandro Dalcin     maxlocal = sf->maxleaf+1; /* TODO: We should use PetscSFGetLeafRange() */
11759566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf,&indegree));
11769566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf,&indegree));
11779566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(sf->nroots+1,&inoffset,maxlocal,&outones,maxlocal,&outoffset));
117895fce210SBarry Smith     inoffset[0] = 0;
117995fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i+1] = inoffset[i] + indegree[i];
11809837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) outones[i] = 1;
11819566063dSJacob Faibussowitsch     PetscCall(PetscSFFetchAndOpBegin(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM));
11829566063dSJacob Faibussowitsch     PetscCall(PetscSFFetchAndOpEnd(sf,MPIU_INT,inoffset,outones,outoffset,MPI_SUM));
118395fce210SBarry Smith     for (i=0; i<sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */
118476bd3646SJed Brown     if (PetscDefined(USE_DEBUG)) { /* Check that the expected number of increments occurred */
118595fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
1186c9cc58a2SBarry Smith         PetscCheck(inoffset[i] + indegree[i] == inoffset[i+1],PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect result after PetscSFFetchAndOp");
118795fce210SBarry Smith       }
118876bd3646SJed Brown     }
11899566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(sf->nleaves,&remote));
119095fce210SBarry Smith     for (i=0; i<sf->nleaves; i++) {
119195fce210SBarry Smith       remote[i].rank  = sf->remote[i].rank;
119238e7336fSToby Isaac       remote[i].index = outoffset[sf->mine ? sf->mine[i] : i];
119395fce210SBarry Smith     }
11949566063dSJacob Faibussowitsch     PetscCall(PetscSFDuplicate(sf,PETSCSF_DUPLICATE_RANKS,&sf->multi));
1195013b3241SStefano Zampini     sf->multi->multi = sf->multi;
11969566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER));
119795fce210SBarry Smith     if (sf->rankorder) {        /* Sort the ranks */
119895fce210SBarry Smith       PetscMPIInt rank;
119995fce210SBarry Smith       PetscInt    *inranks,*newoffset,*outranks,*newoutoffset,*tmpoffset,maxdegree;
120095fce210SBarry Smith       PetscSFNode *newremote;
12019566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank));
120295fce210SBarry Smith       for (i=0,maxdegree=0; i<sf->nroots; i++) maxdegree = PetscMax(maxdegree,indegree[i]);
12039566063dSJacob Faibussowitsch       PetscCall(PetscMalloc5(sf->multi->nroots,&inranks,sf->multi->nroots,&newoffset,maxlocal,&outranks,maxlocal,&newoutoffset,maxdegree,&tmpoffset));
12049837ea96SMatthew G. Knepley       for (i=0; i<maxlocal; i++) outranks[i] = rank;
12059566063dSJacob Faibussowitsch       PetscCall(PetscSFReduceBegin(sf->multi,MPIU_INT,outranks,inranks,MPI_REPLACE));
12069566063dSJacob Faibussowitsch       PetscCall(PetscSFReduceEnd(sf->multi,MPIU_INT,outranks,inranks,MPI_REPLACE));
120795fce210SBarry Smith       /* Sort the incoming ranks at each vertex, build the inverse map */
120895fce210SBarry Smith       for (i=0; i<sf->nroots; i++) {
120995fce210SBarry Smith         PetscInt j;
121095fce210SBarry Smith         for (j=0; j<indegree[i]; j++) tmpoffset[j] = j;
12119566063dSJacob Faibussowitsch         PetscCall(PetscSortIntWithArray(indegree[i],inranks+inoffset[i],tmpoffset));
121295fce210SBarry Smith         for (j=0; j<indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j;
121395fce210SBarry Smith       }
12149566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(sf->multi,MPIU_INT,newoffset,newoutoffset,MPI_REPLACE));
12159566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(sf->multi,MPIU_INT,newoffset,newoutoffset,MPI_REPLACE));
12169566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(sf->nleaves,&newremote));
121795fce210SBarry Smith       for (i=0; i<sf->nleaves; i++) {
121895fce210SBarry Smith         newremote[i].rank  = sf->remote[i].rank;
121901365b40SToby Isaac         newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i];
122095fce210SBarry Smith       }
12219566063dSJacob Faibussowitsch       PetscCall(PetscSFSetGraph(sf->multi,inoffset[sf->nroots],sf->nleaves,sf->mine,PETSC_COPY_VALUES,newremote,PETSC_OWN_POINTER));
12229566063dSJacob Faibussowitsch       PetscCall(PetscFree5(inranks,newoffset,outranks,newoutoffset,tmpoffset));
122395fce210SBarry Smith     }
12249566063dSJacob Faibussowitsch     PetscCall(PetscFree3(inoffset,outones,outoffset));
122595fce210SBarry Smith   }
122695fce210SBarry Smith   *multi = sf->multi;
122795fce210SBarry Smith   PetscFunctionReturn(0);
122895fce210SBarry Smith }
122995fce210SBarry Smith 
123095fce210SBarry Smith /*@C
123172502a1fSJunchao Zhang    PetscSFCreateEmbeddedRootSF - removes edges from all but the selected roots, does not remap indices
123295fce210SBarry Smith 
123395fce210SBarry Smith    Collective
123495fce210SBarry Smith 
12354165533cSJose E. Roman    Input Parameters:
123695fce210SBarry Smith +  sf - original star forest
1237ba2a7774SJunchao Zhang .  nselected  - number of selected roots on this process
1238ba2a7774SJunchao Zhang -  selected   - indices of the selected roots on this process
123995fce210SBarry Smith 
12404165533cSJose E. Roman    Output Parameter:
1241cd620004SJunchao Zhang .  esf - new star forest
124295fce210SBarry Smith 
124395fce210SBarry Smith    Level: advanced
124495fce210SBarry Smith 
124595fce210SBarry Smith    Note:
124695fce210SBarry Smith    To use the new PetscSF, it may be necessary to know the indices of the leaves that are still participating. This can
124795fce210SBarry Smith    be done by calling PetscSFGetGraph().
124895fce210SBarry Smith 
1249db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`, `PetscSFGetGraph()`
125095fce210SBarry Smith @*/
125172502a1fSJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedRootSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *esf)
125295fce210SBarry Smith {
1253cd620004SJunchao Zhang   PetscInt          i,j,n,nroots,nleaves,esf_nleaves,*new_ilocal,minleaf,maxleaf,maxlocal;
1254cd620004SJunchao Zhang   const PetscInt    *ilocal;
1255cd620004SJunchao Zhang   signed char       *rootdata,*leafdata,*leafmem;
1256ba2a7774SJunchao Zhang   const PetscSFNode *iremote;
1257f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1258f659e5c7SJunchao Zhang   MPI_Comm          comm;
125995fce210SBarry Smith 
126095fce210SBarry Smith   PetscFunctionBegin;
126195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
126229046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1263dadcf809SJacob Faibussowitsch   if (nselected) PetscValidIntPointer(selected,3);
1264cd620004SJunchao Zhang   PetscValidPointer(esf,4);
12650511a646SMatthew G. Knepley 
12669566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
12679566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_EmbedSF,sf,0,0,0));
12689566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf,&comm));
12699566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote));
1270cd620004SJunchao Zhang 
127176bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {  /* Error out if selected[] has dups or  out of range indices */
1272cd620004SJunchao Zhang     PetscBool dups;
12739566063dSJacob Faibussowitsch     PetscCall(PetscCheckDupsInt(nselected,selected,&dups));
127428b400f6SJacob Faibussowitsch     PetscCheck(!dups,comm,PETSC_ERR_ARG_WRONG,"selected[] has dups");
1275cd620004SJunchao Zhang     for (i=0; i<nselected; i++)
1276c9cc58a2SBarry Smith       PetscCheck(selected[i] >= 0 && selected[i] < nroots,comm,PETSC_ERR_ARG_OUTOFRANGE,"selected root indice %" PetscInt_FMT " is out of [0,%" PetscInt_FMT ")",selected[i],nroots);
1277cd620004SJunchao Zhang   }
1278f659e5c7SJunchao Zhang 
1279*dbbe0bcdSBarry Smith   if (sf->ops->CreateEmbeddedRootSF) PetscUseTypeMethod(sf,CreateEmbeddedRootSF ,nselected,selected,esf);
1280*dbbe0bcdSBarry Smith   else {
1281cd620004SJunchao Zhang     /* A generic version of creating embedded sf */
12829566063dSJacob Faibussowitsch     PetscCall(PetscSFGetLeafRange(sf,&minleaf,&maxleaf));
1283cd620004SJunchao Zhang     maxlocal = maxleaf - minleaf + 1;
12849566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(nroots,&rootdata,maxlocal,&leafmem));
1285cd620004SJunchao Zhang     leafdata = leafmem - minleaf;
1286cd620004SJunchao Zhang     /* Tag selected roots and bcast to leaves */
1287cd620004SJunchao Zhang     for (i=0; i<nselected; i++) rootdata[selected[i]] = 1;
12889566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf,MPI_SIGNED_CHAR,rootdata,leafdata,MPI_REPLACE));
12899566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf,MPI_SIGNED_CHAR,rootdata,leafdata,MPI_REPLACE));
1290ba2a7774SJunchao Zhang 
1291cd620004SJunchao Zhang     /* Build esf with leaves that are still connected */
1292cd620004SJunchao Zhang     esf_nleaves = 0;
1293cd620004SJunchao Zhang     for (i=0; i<nleaves; i++) {
1294cd620004SJunchao Zhang       j = ilocal ? ilocal[i] : i;
1295cd620004SJunchao Zhang       /* esf_nleaves += leafdata[j] should work in theory, but failed with SFWindow bugs
1296cd620004SJunchao Zhang          with PetscSFBcast. See https://gitlab.com/petsc/petsc/issues/555
1297cd620004SJunchao Zhang       */
1298cd620004SJunchao Zhang       esf_nleaves += (leafdata[j] ? 1 : 0);
1299cd620004SJunchao Zhang     }
13009566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(esf_nleaves,&new_ilocal));
13019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(esf_nleaves,&new_iremote));
1302cd620004SJunchao Zhang     for (i=n=0; i<nleaves; i++) {
1303cd620004SJunchao Zhang       j = ilocal ? ilocal[i] : i;
1304cd620004SJunchao Zhang       if (leafdata[j]) {
1305cd620004SJunchao Zhang         new_ilocal[n]        = j;
1306cd620004SJunchao Zhang         new_iremote[n].rank  = iremote[i].rank;
1307cd620004SJunchao Zhang         new_iremote[n].index = iremote[i].index;
1308fc1ede2bSMatthew G. Knepley         ++n;
130995fce210SBarry Smith       }
131095fce210SBarry Smith     }
13119566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(comm,esf));
13129566063dSJacob Faibussowitsch     PetscCall(PetscSFSetFromOptions(*esf));
13139566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*esf,nroots,esf_nleaves,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER));
13149566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata,leafmem));
1315f659e5c7SJunchao Zhang   }
13169566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_EmbedSF,sf,0,0,0));
131795fce210SBarry Smith   PetscFunctionReturn(0);
131895fce210SBarry Smith }
131995fce210SBarry Smith 
13202f5fb4c2SMatthew G. Knepley /*@C
13212f5fb4c2SMatthew G. Knepley   PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves, does not remap indices
13222f5fb4c2SMatthew G. Knepley 
13232f5fb4c2SMatthew G. Knepley   Collective
13242f5fb4c2SMatthew G. Knepley 
13254165533cSJose E. Roman   Input Parameters:
13262f5fb4c2SMatthew G. Knepley + sf - original star forest
1327f659e5c7SJunchao Zhang . nselected  - number of selected leaves on this process
1328f659e5c7SJunchao Zhang - selected   - indices of the selected leaves on this process
13292f5fb4c2SMatthew G. Knepley 
13304165533cSJose E. Roman   Output Parameter:
13312f5fb4c2SMatthew G. Knepley .  newsf - new star forest
13322f5fb4c2SMatthew G. Knepley 
13332f5fb4c2SMatthew G. Knepley   Level: advanced
13342f5fb4c2SMatthew G. Knepley 
1335db781477SPatrick Sanan .seealso: `PetscSFCreateEmbeddedRootSF()`, `PetscSFSetGraph()`, `PetscSFGetGraph()`
13362f5fb4c2SMatthew G. Knepley @*/
1337f659e5c7SJunchao Zhang PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf,PetscInt nselected,const PetscInt *selected,PetscSF *newsf)
13382f5fb4c2SMatthew G. Knepley {
1339f659e5c7SJunchao Zhang   const PetscSFNode *iremote;
1340f659e5c7SJunchao Zhang   PetscSFNode       *new_iremote;
1341f659e5c7SJunchao Zhang   const PetscInt    *ilocal;
1342f659e5c7SJunchao Zhang   PetscInt          i,nroots,*leaves,*new_ilocal;
1343f659e5c7SJunchao Zhang   MPI_Comm          comm;
13442f5fb4c2SMatthew G. Knepley 
13452f5fb4c2SMatthew G. Knepley   PetscFunctionBegin;
13462f5fb4c2SMatthew G. Knepley   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
134729046d53SLisandro Dalcin   PetscSFCheckGraphSet(sf,1);
1348dadcf809SJacob Faibussowitsch   if (nselected) PetscValidIntPointer(selected,3);
13492f5fb4c2SMatthew G. Knepley   PetscValidPointer(newsf,4);
13502f5fb4c2SMatthew G. Knepley 
1351f659e5c7SJunchao Zhang   /* Uniq selected[] and put results in leaves[] */
13529566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf,&comm));
13539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nselected,&leaves));
13549566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(leaves,selected,nselected));
13559566063dSJacob Faibussowitsch   PetscCall(PetscSortedRemoveDupsInt(&nselected,leaves));
135608401ef6SPierre 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);
1357f659e5c7SJunchao Zhang 
1358f659e5c7SJunchao Zhang   /* Optimize the routine only when sf is setup and hence we can reuse sf's communication pattern */
1359*dbbe0bcdSBarry Smith   if (sf->setupcalled && sf->ops->CreateEmbeddedLeafSF) PetscUseTypeMethod(sf,CreateEmbeddedLeafSF ,nselected,leaves,newsf);
1360*dbbe0bcdSBarry Smith   else {
13619566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sf,&nroots,NULL,&ilocal,&iremote));
13629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nselected,&new_ilocal));
13639566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nselected,&new_iremote));
1364f659e5c7SJunchao Zhang     for (i=0; i<nselected; ++i) {
1365f659e5c7SJunchao Zhang       const PetscInt l     = leaves[i];
1366f659e5c7SJunchao Zhang       new_ilocal[i]        = ilocal ? ilocal[l] : l;
1367f659e5c7SJunchao Zhang       new_iremote[i].rank  = iremote[l].rank;
1368f659e5c7SJunchao Zhang       new_iremote[i].index = iremote[l].index;
13692f5fb4c2SMatthew G. Knepley     }
13709566063dSJacob Faibussowitsch     PetscCall(PetscSFDuplicate(sf,PETSCSF_DUPLICATE_CONFONLY,newsf));
13719566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*newsf,nroots,nselected,new_ilocal,PETSC_OWN_POINTER,new_iremote,PETSC_OWN_POINTER));
1372f659e5c7SJunchao Zhang   }
13739566063dSJacob Faibussowitsch   PetscCall(PetscFree(leaves));
13742f5fb4c2SMatthew G. Knepley   PetscFunctionReturn(0);
13752f5fb4c2SMatthew G. Knepley }
13762f5fb4c2SMatthew G. Knepley 
137795fce210SBarry Smith /*@C
1378ad227feaSJunchao Zhang    PetscSFBcastBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to PetscSFBcastEnd()
13793482bfa8SJunchao Zhang 
13803482bfa8SJunchao Zhang    Collective on PetscSF
13813482bfa8SJunchao Zhang 
13824165533cSJose E. Roman    Input Parameters:
13833482bfa8SJunchao Zhang +  sf - star forest on which to communicate
13843482bfa8SJunchao Zhang .  unit - data type associated with each node
13853482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
13863482bfa8SJunchao Zhang -  op - operation to use for reduction
13873482bfa8SJunchao Zhang 
13884165533cSJose E. Roman    Output Parameter:
13893482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
13903482bfa8SJunchao Zhang 
13913482bfa8SJunchao Zhang    Level: intermediate
13923482bfa8SJunchao Zhang 
1393d0295fc0SJunchao Zhang    Notes:
1394d0295fc0SJunchao Zhang     When petsc is configured with device support, it will use its own mechanism to figure out whether the given data pointers
1395d0295fc0SJunchao Zhang     are host pointers or device pointers, which may incur a noticable cost. If you already knew the info, you should
1396ad227feaSJunchao Zhang     use PetscSFBcastWithMemTypeBegin() instead.
1397db781477SPatrick Sanan .seealso: `PetscSFBcastEnd()`, `PetscSFBcastWithMemTypeBegin()`
13983482bfa8SJunchao Zhang @*/
1399ad227feaSJunchao Zhang PetscErrorCode PetscSFBcastBegin(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
14003482bfa8SJunchao Zhang {
1401eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
14023482bfa8SJunchao Zhang 
14033482bfa8SJunchao Zhang   PetscFunctionBegin;
14043482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
14059566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
14069566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0));
14079566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(rootdata,&rootmtype));
14089566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafdata,&leafmtype));
1409*dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf,BcastBegin ,unit,rootmtype,rootdata,leafmtype,leafdata,op);
14109566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0));
14113482bfa8SJunchao Zhang   PetscFunctionReturn(0);
14123482bfa8SJunchao Zhang }
14133482bfa8SJunchao Zhang 
14143482bfa8SJunchao Zhang /*@C
1415ad227feaSJunchao Zhang    PetscSFBcastWithMemTypeBegin - begin pointwise broadcast with root value being reduced to leaf value with explicit memory types, to be concluded with call to PetscSFBcastEnd()
1416d0295fc0SJunchao Zhang 
1417d0295fc0SJunchao Zhang    Collective on PetscSF
1418d0295fc0SJunchao Zhang 
14194165533cSJose E. Roman    Input Parameters:
1420d0295fc0SJunchao Zhang +  sf - star forest on which to communicate
1421d0295fc0SJunchao Zhang .  unit - data type associated with each node
1422d0295fc0SJunchao Zhang .  rootmtype - memory type of rootdata
1423d0295fc0SJunchao Zhang .  rootdata - buffer to broadcast
1424d0295fc0SJunchao Zhang .  leafmtype - memory type of leafdata
1425d0295fc0SJunchao Zhang -  op - operation to use for reduction
1426d0295fc0SJunchao Zhang 
14274165533cSJose E. Roman    Output Parameter:
1428d0295fc0SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
1429d0295fc0SJunchao Zhang 
1430d0295fc0SJunchao Zhang    Level: intermediate
1431d0295fc0SJunchao Zhang 
1432db781477SPatrick Sanan .seealso: `PetscSFBcastEnd()`, `PetscSFBcastBegin()`
1433d0295fc0SJunchao Zhang @*/
1434ad227feaSJunchao Zhang PetscErrorCode PetscSFBcastWithMemTypeBegin(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,const void *rootdata,PetscMemType leafmtype,void *leafdata,MPI_Op op)
1435d0295fc0SJunchao Zhang {
1436d0295fc0SJunchao Zhang   PetscFunctionBegin;
1437d0295fc0SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
14389566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
14399566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0));
1440*dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf,BcastBegin ,unit,rootmtype,rootdata,leafmtype,leafdata,op);
14419566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0));
1442d0295fc0SJunchao Zhang   PetscFunctionReturn(0);
1443d0295fc0SJunchao Zhang }
1444d0295fc0SJunchao Zhang 
1445d0295fc0SJunchao Zhang /*@C
1446ad227feaSJunchao Zhang    PetscSFBcastEnd - end a broadcast & reduce operation started with PetscSFBcastBegin()
14473482bfa8SJunchao Zhang 
14483482bfa8SJunchao Zhang    Collective
14493482bfa8SJunchao Zhang 
14504165533cSJose E. Roman    Input Parameters:
14513482bfa8SJunchao Zhang +  sf - star forest
14523482bfa8SJunchao Zhang .  unit - data type
14533482bfa8SJunchao Zhang .  rootdata - buffer to broadcast
14543482bfa8SJunchao Zhang -  op - operation to use for reduction
14553482bfa8SJunchao Zhang 
14564165533cSJose E. Roman    Output Parameter:
14573482bfa8SJunchao Zhang .  leafdata - buffer to be reduced with values from each leaf's respective root
14583482bfa8SJunchao Zhang 
14593482bfa8SJunchao Zhang    Level: intermediate
14603482bfa8SJunchao Zhang 
1461db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`, `PetscSFReduceEnd()`
14623482bfa8SJunchao Zhang @*/
1463ad227feaSJunchao Zhang PetscErrorCode PetscSFBcastEnd(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
14643482bfa8SJunchao Zhang {
14653482bfa8SJunchao Zhang   PetscFunctionBegin;
14663482bfa8SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
14679566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_BcastEnd,sf,0,0,0));
1468*dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf,BcastEnd ,unit,rootdata,leafdata,op);
14699566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_BcastEnd,sf,0,0,0));
14703482bfa8SJunchao Zhang   PetscFunctionReturn(0);
14713482bfa8SJunchao Zhang }
14723482bfa8SJunchao Zhang 
14733482bfa8SJunchao Zhang /*@C
147495fce210SBarry Smith    PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to PetscSFReduceEnd()
147595fce210SBarry Smith 
147695fce210SBarry Smith    Collective
147795fce210SBarry Smith 
14784165533cSJose E. Roman    Input Parameters:
147995fce210SBarry Smith +  sf - star forest
148095fce210SBarry Smith .  unit - data type
148195fce210SBarry Smith .  leafdata - values to reduce
148295fce210SBarry Smith -  op - reduction operation
148395fce210SBarry Smith 
14844165533cSJose E. Roman    Output Parameter:
148595fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
148695fce210SBarry Smith 
148795fce210SBarry Smith    Level: intermediate
148895fce210SBarry Smith 
1489d0295fc0SJunchao Zhang    Notes:
1490d0295fc0SJunchao Zhang     When petsc is configured with device support, it will use its own mechanism to figure out whether the given data pointers
1491d0295fc0SJunchao Zhang     are host pointers or device pointers, which may incur a noticable cost. If you already knew the info, you should
1492d0295fc0SJunchao Zhang     use PetscSFReduceWithMemTypeBegin() instead.
1493d0295fc0SJunchao Zhang 
1494db781477SPatrick Sanan .seealso: `PetscSFBcastBegin()`, `PetscSFReduceWithMemTypeBegin()`
149595fce210SBarry Smith @*/
149695fce210SBarry Smith PetscErrorCode PetscSFReduceBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
149795fce210SBarry Smith {
1498eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype;
149995fce210SBarry Smith 
150095fce210SBarry Smith   PetscFunctionBegin;
150195fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
15029566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
15039566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0));
15049566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(rootdata,&rootmtype));
15059566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafdata,&leafmtype));
15069566063dSJacob Faibussowitsch   PetscCall((sf->ops->ReduceBegin)(sf,unit,leafmtype,leafdata,rootmtype,rootdata,op));
15079566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_ReduceBegin,sf,0,0,0));
150895fce210SBarry Smith   PetscFunctionReturn(0);
150995fce210SBarry Smith }
151095fce210SBarry Smith 
151195fce210SBarry Smith /*@C
1512d0295fc0SJunchao Zhang    PetscSFReduceWithMemTypeBegin - begin reduction of leafdata into rootdata with explicit memory types, to be completed with call to PetscSFReduceEnd()
1513d0295fc0SJunchao Zhang 
1514d0295fc0SJunchao Zhang    Collective
1515d0295fc0SJunchao Zhang 
15164165533cSJose E. Roman    Input Parameters:
1517d0295fc0SJunchao Zhang +  sf - star forest
1518d0295fc0SJunchao Zhang .  unit - data type
1519d0295fc0SJunchao Zhang .  leafmtype - memory type of leafdata
1520d0295fc0SJunchao Zhang .  leafdata - values to reduce
1521d0295fc0SJunchao Zhang .  rootmtype - memory type of rootdata
1522d0295fc0SJunchao Zhang -  op - reduction operation
1523d0295fc0SJunchao Zhang 
15244165533cSJose E. Roman    Output Parameter:
1525d0295fc0SJunchao Zhang .  rootdata - result of reduction of values from all leaves of each root
1526d0295fc0SJunchao Zhang 
1527d0295fc0SJunchao Zhang    Level: intermediate
1528d0295fc0SJunchao Zhang 
1529db781477SPatrick Sanan .seealso: `PetscSFBcastBegin()`, `PetscSFReduceBegin()`
1530d0295fc0SJunchao Zhang @*/
1531d0295fc0SJunchao Zhang PetscErrorCode PetscSFReduceWithMemTypeBegin(PetscSF sf,MPI_Datatype unit,PetscMemType leafmtype,const void *leafdata,PetscMemType rootmtype,void *rootdata,MPI_Op op)
1532d0295fc0SJunchao Zhang {
1533d0295fc0SJunchao Zhang   PetscFunctionBegin;
1534d0295fc0SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
15359566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
15369566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_ReduceBegin,sf,0,0,0));
15379566063dSJacob Faibussowitsch   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));
1539d0295fc0SJunchao Zhang   PetscFunctionReturn(0);
1540d0295fc0SJunchao Zhang }
1541d0295fc0SJunchao Zhang 
1542d0295fc0SJunchao Zhang /*@C
154395fce210SBarry Smith    PetscSFReduceEnd - end a reduction operation started with PetscSFReduceBegin()
154495fce210SBarry Smith 
154595fce210SBarry Smith    Collective
154695fce210SBarry Smith 
15474165533cSJose E. Roman    Input Parameters:
154895fce210SBarry Smith +  sf - star forest
154995fce210SBarry Smith .  unit - data type
155095fce210SBarry Smith .  leafdata - values to reduce
155195fce210SBarry Smith -  op - reduction operation
155295fce210SBarry Smith 
15534165533cSJose E. Roman    Output Parameter:
155495fce210SBarry Smith .  rootdata - result of reduction of values from all leaves of each root
155595fce210SBarry Smith 
155695fce210SBarry Smith    Level: intermediate
155795fce210SBarry Smith 
1558db781477SPatrick Sanan .seealso: `PetscSFSetGraph()`, `PetscSFBcastEnd()`
155995fce210SBarry Smith @*/
156095fce210SBarry Smith PetscErrorCode PetscSFReduceEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
156195fce210SBarry Smith {
156295fce210SBarry Smith   PetscFunctionBegin;
156395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
15649566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_ReduceEnd,sf,0,0,0));
1565*dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf,ReduceEnd ,unit,leafdata,rootdata,op);
15669566063dSJacob Faibussowitsch   if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_ReduceEnd,sf,0,0,0));
156795fce210SBarry Smith   PetscFunctionReturn(0);
156895fce210SBarry Smith }
156995fce210SBarry Smith 
157095fce210SBarry Smith /*@C
1571a1729e3fSJunchao Zhang    PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
1572a1729e3fSJunchao Zhang 
1573a1729e3fSJunchao Zhang    Collective
1574a1729e3fSJunchao Zhang 
15754165533cSJose E. Roman    Input Parameters:
1576a1729e3fSJunchao Zhang +  sf - star forest
1577a1729e3fSJunchao Zhang .  unit - data type
1578a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1579a1729e3fSJunchao Zhang -  op - operation to use for reduction
1580a1729e3fSJunchao Zhang 
15814165533cSJose E. Roman    Output Parameters:
1582a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1583a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1584a1729e3fSJunchao Zhang 
1585a1729e3fSJunchao Zhang    Level: advanced
1586a1729e3fSJunchao Zhang 
1587a1729e3fSJunchao Zhang    Note:
1588a1729e3fSJunchao Zhang    The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process
1589a1729e3fSJunchao Zhang    might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is
1590a1729e3fSJunchao Zhang    not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as
1591a1729e3fSJunchao Zhang    integers.
1592a1729e3fSJunchao Zhang 
1593db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeBegin()`, `PetscSFReduceBegin()`, `PetscSFSetGraph()`
1594a1729e3fSJunchao Zhang @*/
1595a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1596a1729e3fSJunchao Zhang {
1597eb02082bSJunchao Zhang   PetscMemType   rootmtype,leafmtype,leafupdatemtype;
1598a1729e3fSJunchao Zhang 
1599a1729e3fSJunchao Zhang   PetscFunctionBegin;
1600a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
16019566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
16029566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0));
16039566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(rootdata,&rootmtype));
16049566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafdata,&leafmtype));
16059566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafupdate,&leafupdatemtype));
160608401ef6SPierre Jolivet   PetscCheck(leafmtype == leafupdatemtype,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for leafdata and leafupdate in different memory types");
1607*dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf,FetchAndOpBegin ,unit,rootmtype,rootdata,leafmtype,leafdata,leafupdate,op);
16089566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0));
1609a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1610a1729e3fSJunchao Zhang }
1611a1729e3fSJunchao Zhang 
1612a1729e3fSJunchao Zhang /*@C
1613d3b3e55cSJunchao Zhang    PetscSFFetchAndOpWithMemTypeBegin - begin operation with explicit memory types that fetches values from root and updates atomically by applying operation using my leaf value, to be completed with PetscSFFetchAndOpEnd()
1614d3b3e55cSJunchao Zhang 
1615d3b3e55cSJunchao Zhang    Collective
1616d3b3e55cSJunchao Zhang 
1617d3b3e55cSJunchao Zhang    Input Parameters:
1618d3b3e55cSJunchao Zhang +  sf - star forest
1619d3b3e55cSJunchao Zhang .  unit - data type
1620d3b3e55cSJunchao Zhang .  rootmtype - memory type of rootdata
1621d3b3e55cSJunchao Zhang .  leafmtype - memory type of leafdata
1622d3b3e55cSJunchao Zhang .  leafdata - leaf values to use in reduction
1623d3b3e55cSJunchao Zhang .  leafupdatemtype - memory type of leafupdate
1624d3b3e55cSJunchao Zhang -  op - operation to use for reduction
1625d3b3e55cSJunchao Zhang 
1626d3b3e55cSJunchao Zhang    Output Parameters:
1627d3b3e55cSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1628d3b3e55cSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1629d3b3e55cSJunchao Zhang 
1630d3b3e55cSJunchao Zhang    Level: advanced
1631d3b3e55cSJunchao Zhang 
1632d3b3e55cSJunchao Zhang    Note: See PetscSFFetchAndOpBegin() for more details.
1633d3b3e55cSJunchao Zhang 
1634c2e3fba1SPatrick Sanan .seealso: `PetscSFFetchAndOpBegin()`, `PetscSFComputeDegreeBegin()`, `PetscSFReduceBegin()`, `PetscSFSetGraph()`
1635d3b3e55cSJunchao Zhang @*/
1636d3b3e55cSJunchao Zhang PetscErrorCode PetscSFFetchAndOpWithMemTypeBegin(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,void *rootdata,PetscMemType leafmtype,const void *leafdata,PetscMemType leafupdatemtype,void *leafupdate,MPI_Op op)
1637d3b3e55cSJunchao Zhang {
1638d3b3e55cSJunchao Zhang   PetscFunctionBegin;
1639d3b3e55cSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
16409566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
16419566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_FetchAndOpBegin,sf,0,0,0));
164208401ef6SPierre Jolivet   PetscCheck(leafmtype == leafupdatemtype,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for leafdata and leafupdate in different memory types");
1643*dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf,FetchAndOpBegin ,unit,rootmtype,rootdata,leafmtype,leafdata,leafupdate,op);
16449566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_FetchAndOpBegin,sf,0,0,0));
1645d3b3e55cSJunchao Zhang   PetscFunctionReturn(0);
1646d3b3e55cSJunchao Zhang }
1647d3b3e55cSJunchao Zhang 
1648d3b3e55cSJunchao Zhang /*@C
1649a1729e3fSJunchao Zhang    PetscSFFetchAndOpEnd - end operation started in matching call to PetscSFFetchAndOpBegin() to fetch values from roots and update atomically by applying operation using my leaf value
1650a1729e3fSJunchao Zhang 
1651a1729e3fSJunchao Zhang    Collective
1652a1729e3fSJunchao Zhang 
16534165533cSJose E. Roman    Input Parameters:
1654a1729e3fSJunchao Zhang +  sf - star forest
1655a1729e3fSJunchao Zhang .  unit - data type
1656a1729e3fSJunchao Zhang .  leafdata - leaf values to use in reduction
1657a1729e3fSJunchao Zhang -  op - operation to use for reduction
1658a1729e3fSJunchao Zhang 
16594165533cSJose E. Roman    Output Parameters:
1660a1729e3fSJunchao Zhang +  rootdata - root values to be updated, input state is seen by first process to perform an update
1661a1729e3fSJunchao Zhang -  leafupdate - state at each leaf's respective root immediately prior to my atomic update
1662a1729e3fSJunchao Zhang 
1663a1729e3fSJunchao Zhang    Level: advanced
1664a1729e3fSJunchao Zhang 
1665db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeEnd()`, `PetscSFReduceEnd()`, `PetscSFSetGraph()`
1666a1729e3fSJunchao Zhang @*/
1667a1729e3fSJunchao Zhang PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
1668a1729e3fSJunchao Zhang {
1669a1729e3fSJunchao Zhang   PetscFunctionBegin;
1670a1729e3fSJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
16719566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_FetchAndOpEnd,sf,0,0,0));
1672*dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf,FetchAndOpEnd ,unit,rootdata,leafdata,leafupdate,op);
16739566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_FetchAndOpEnd,sf,0,0,0));
1674a1729e3fSJunchao Zhang   PetscFunctionReturn(0);
1675a1729e3fSJunchao Zhang }
1676a1729e3fSJunchao Zhang 
1677a1729e3fSJunchao Zhang /*@C
167895fce210SBarry Smith    PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with PetscSFComputeDegreeEnd()
167995fce210SBarry Smith 
168095fce210SBarry Smith    Collective
168195fce210SBarry Smith 
16824165533cSJose E. Roman    Input Parameter:
168395fce210SBarry Smith .  sf - star forest
168495fce210SBarry Smith 
16854165533cSJose E. Roman    Output Parameter:
168695fce210SBarry Smith .  degree - degree of each root vertex
168795fce210SBarry Smith 
168895fce210SBarry Smith    Level: advanced
168995fce210SBarry Smith 
1690ffe67aa5SVáclav Hapla    Notes:
1691ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1692ffe67aa5SVáclav Hapla 
1693db781477SPatrick Sanan .seealso: `PetscSFGatherBegin()`
169495fce210SBarry Smith @*/
169595fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf,const PetscInt **degree)
169695fce210SBarry Smith {
169795fce210SBarry Smith   PetscFunctionBegin;
169895fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
169995fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
170095fce210SBarry Smith   PetscValidPointer(degree,2);
1701803bd9e8SMatthew G. Knepley   if (!sf->degreeknown) {
17025b0d146aSStefano Zampini     PetscInt i, nroots = sf->nroots, maxlocal;
170328b400f6SJacob Faibussowitsch     PetscCheck(!sf->degree,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested.");
17045b0d146aSStefano Zampini     maxlocal = sf->maxleaf-sf->minleaf+1;
17059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nroots,&sf->degree));
17069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PetscMax(maxlocal,1),&sf->degreetmp)); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */
170729046d53SLisandro Dalcin     for (i=0; i<nroots; i++) sf->degree[i] = 0;
17089837ea96SMatthew G. Knepley     for (i=0; i<maxlocal; i++) sf->degreetmp[i] = 1;
17099566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf,MPIU_INT,sf->degreetmp-sf->minleaf,sf->degree,MPI_SUM));
171095fce210SBarry Smith   }
171195fce210SBarry Smith   *degree = NULL;
171295fce210SBarry Smith   PetscFunctionReturn(0);
171395fce210SBarry Smith }
171495fce210SBarry Smith 
171595fce210SBarry Smith /*@C
171695fce210SBarry Smith    PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with PetscSFComputeDegreeBegin()
171795fce210SBarry Smith 
171895fce210SBarry Smith    Collective
171995fce210SBarry Smith 
17204165533cSJose E. Roman    Input Parameter:
172195fce210SBarry Smith .  sf - star forest
172295fce210SBarry Smith 
17234165533cSJose E. Roman    Output Parameter:
172495fce210SBarry Smith .  degree - degree of each root vertex
172595fce210SBarry Smith 
172695fce210SBarry Smith    Level: developer
172795fce210SBarry Smith 
1728ffe67aa5SVáclav Hapla    Notes:
1729ffe67aa5SVáclav Hapla    The returned array is owned by PetscSF and automatically freed by PetscSFDestroy(). Hence no need to call PetscFree() on it.
1730ffe67aa5SVáclav Hapla 
173195fce210SBarry Smith .seealso:
173295fce210SBarry Smith @*/
173395fce210SBarry Smith PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf,const PetscInt **degree)
173495fce210SBarry Smith {
173595fce210SBarry Smith   PetscFunctionBegin;
173695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
173795fce210SBarry Smith   PetscSFCheckGraphSet(sf,1);
173829046d53SLisandro Dalcin   PetscValidPointer(degree,2);
173995fce210SBarry Smith   if (!sf->degreeknown) {
174028b400f6SJacob Faibussowitsch     PetscCheck(sf->degreetmp,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()");
17419566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf,MPIU_INT,sf->degreetmp-sf->minleaf,sf->degree,MPI_SUM));
17429566063dSJacob Faibussowitsch     PetscCall(PetscFree(sf->degreetmp));
174395fce210SBarry Smith     sf->degreeknown = PETSC_TRUE;
174495fce210SBarry Smith   }
174595fce210SBarry Smith   *degree = sf->degree;
174695fce210SBarry Smith   PetscFunctionReturn(0);
174795fce210SBarry Smith }
174895fce210SBarry Smith 
1749673100f5SVaclav Hapla /*@C
175066dfcd1aSVaclav Hapla    PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-SF returned by PetscSFGetMultiSF()).
175166dfcd1aSVaclav Hapla    Each multi-root is assigned index of the corresponding original root.
1752673100f5SVaclav Hapla 
1753673100f5SVaclav Hapla    Collective
1754673100f5SVaclav Hapla 
17554165533cSJose E. Roman    Input Parameters:
1756673100f5SVaclav Hapla +  sf - star forest
1757673100f5SVaclav Hapla -  degree - degree of each root vertex, computed with PetscSFComputeDegreeBegin()/PetscSFComputeDegreeEnd()
1758673100f5SVaclav Hapla 
17594165533cSJose E. Roman    Output Parameters:
176066dfcd1aSVaclav Hapla +  nMultiRoots - (optional) number of multi-roots (roots of multi-SF)
176166dfcd1aSVaclav Hapla -  multiRootsOrigNumbering - original indices of multi-roots; length of this array is nMultiRoots
1762673100f5SVaclav Hapla 
1763673100f5SVaclav Hapla    Level: developer
1764673100f5SVaclav Hapla 
1765ffe67aa5SVáclav Hapla    Notes:
1766ffe67aa5SVáclav Hapla    The returned array multiRootsOrigNumbering is newly allocated and should be destroyed with PetscFree() when no longer needed.
1767ffe67aa5SVáclav Hapla 
1768db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeBegin()`, `PetscSFComputeDegreeEnd()`, `PetscSFGetMultiSF()`
1769673100f5SVaclav Hapla @*/
177066dfcd1aSVaclav Hapla PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[])
1771673100f5SVaclav Hapla {
1772673100f5SVaclav Hapla   PetscSF             msf;
1773673100f5SVaclav Hapla   PetscInt            i, j, k, nroots, nmroots;
1774673100f5SVaclav Hapla 
1775673100f5SVaclav Hapla   PetscFunctionBegin;
1776673100f5SVaclav Hapla   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
17779566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
177829328920SVaclav Hapla   if (nroots) PetscValidIntPointer(degree,2);
177966dfcd1aSVaclav Hapla   if (nMultiRoots) PetscValidIntPointer(nMultiRoots,3);
178066dfcd1aSVaclav Hapla   PetscValidPointer(multiRootsOrigNumbering,4);
17819566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf,&msf));
17829566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL));
17839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nmroots, multiRootsOrigNumbering));
1784673100f5SVaclav Hapla   for (i=0,j=0,k=0; i<nroots; i++) {
1785673100f5SVaclav Hapla     if (!degree[i]) continue;
1786673100f5SVaclav Hapla     for (j=0; j<degree[i]; j++,k++) {
178766dfcd1aSVaclav Hapla       (*multiRootsOrigNumbering)[k] = i;
1788673100f5SVaclav Hapla     }
1789673100f5SVaclav Hapla   }
179008401ef6SPierre Jolivet   PetscCheck(k == nmroots,PETSC_COMM_SELF,PETSC_ERR_PLIB,"sanity check fail");
179166dfcd1aSVaclav Hapla   if (nMultiRoots) *nMultiRoots = nmroots;
1792673100f5SVaclav Hapla   PetscFunctionReturn(0);
1793673100f5SVaclav Hapla }
1794673100f5SVaclav Hapla 
179595fce210SBarry Smith /*@C
179695fce210SBarry Smith    PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with PetscSFGatherEnd()
179795fce210SBarry Smith 
179895fce210SBarry Smith    Collective
179995fce210SBarry Smith 
18004165533cSJose E. Roman    Input Parameters:
180195fce210SBarry Smith +  sf - star forest
180295fce210SBarry Smith .  unit - data type
180395fce210SBarry Smith -  leafdata - leaf data to gather to roots
180495fce210SBarry Smith 
18054165533cSJose E. Roman    Output Parameter:
180695fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
180795fce210SBarry Smith 
180895fce210SBarry Smith    Level: intermediate
180995fce210SBarry Smith 
1810db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeBegin()`, `PetscSFScatterBegin()`
181195fce210SBarry Smith @*/
181295fce210SBarry Smith PetscErrorCode PetscSFGatherBegin(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
181395fce210SBarry Smith {
1814a5526d50SJunchao Zhang   PetscSF        multi = NULL;
181595fce210SBarry Smith 
181695fce210SBarry Smith   PetscFunctionBegin;
181795fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
18189566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
18199566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf,&multi));
18209566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(multi,unit,leafdata,multirootdata,MPI_REPLACE));
182195fce210SBarry Smith   PetscFunctionReturn(0);
182295fce210SBarry Smith }
182395fce210SBarry Smith 
182495fce210SBarry Smith /*@C
182595fce210SBarry Smith    PetscSFGatherEnd - ends pointwise gather operation that was started with PetscSFGatherBegin()
182695fce210SBarry Smith 
182795fce210SBarry Smith    Collective
182895fce210SBarry Smith 
18294165533cSJose E. Roman    Input Parameters:
183095fce210SBarry Smith +  sf - star forest
183195fce210SBarry Smith .  unit - data type
183295fce210SBarry Smith -  leafdata - leaf data to gather to roots
183395fce210SBarry Smith 
18344165533cSJose E. Roman    Output Parameter:
183595fce210SBarry Smith .  multirootdata - root buffer to gather into, amount of space per root is equal to its degree
183695fce210SBarry Smith 
183795fce210SBarry Smith    Level: intermediate
183895fce210SBarry Smith 
1839db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeEnd()`, `PetscSFScatterEnd()`
184095fce210SBarry Smith @*/
184195fce210SBarry Smith PetscErrorCode PetscSFGatherEnd(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *multirootdata)
184295fce210SBarry Smith {
1843a5526d50SJunchao Zhang   PetscSF        multi = NULL;
184495fce210SBarry Smith 
184595fce210SBarry Smith   PetscFunctionBegin;
184695fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
18479566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf,&multi));
18489566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(multi,unit,leafdata,multirootdata,MPI_REPLACE));
184995fce210SBarry Smith   PetscFunctionReturn(0);
185095fce210SBarry Smith }
185195fce210SBarry Smith 
185295fce210SBarry Smith /*@C
185395fce210SBarry Smith    PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with PetscSFScatterEnd()
185495fce210SBarry Smith 
185595fce210SBarry Smith    Collective
185695fce210SBarry Smith 
18574165533cSJose E. Roman    Input Parameters:
185895fce210SBarry Smith +  sf - star forest
185995fce210SBarry Smith .  unit - data type
186095fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
186195fce210SBarry Smith 
18624165533cSJose E. Roman    Output Parameter:
186395fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
186495fce210SBarry Smith 
186595fce210SBarry Smith    Level: intermediate
186695fce210SBarry Smith 
1867db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeBegin()`, `PetscSFScatterBegin()`
186895fce210SBarry Smith @*/
186995fce210SBarry Smith PetscErrorCode PetscSFScatterBegin(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
187095fce210SBarry Smith {
1871a5526d50SJunchao Zhang   PetscSF        multi = NULL;
187295fce210SBarry Smith 
187395fce210SBarry Smith   PetscFunctionBegin;
187495fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
18759566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
18769566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf,&multi));
18779566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(multi,unit,multirootdata,leafdata,MPI_REPLACE));
187895fce210SBarry Smith   PetscFunctionReturn(0);
187995fce210SBarry Smith }
188095fce210SBarry Smith 
188195fce210SBarry Smith /*@C
188295fce210SBarry Smith    PetscSFScatterEnd - ends pointwise scatter operation that was started with PetscSFScatterBegin()
188395fce210SBarry Smith 
188495fce210SBarry Smith    Collective
188595fce210SBarry Smith 
18864165533cSJose E. Roman    Input Parameters:
188795fce210SBarry Smith +  sf - star forest
188895fce210SBarry Smith .  unit - data type
188995fce210SBarry Smith -  multirootdata - root buffer to send to each leaf, one unit of data per leaf
189095fce210SBarry Smith 
18914165533cSJose E. Roman    Output Parameter:
189295fce210SBarry Smith .  leafdata - leaf data to be update with personal data from each respective root
189395fce210SBarry Smith 
189495fce210SBarry Smith    Level: intermediate
189595fce210SBarry Smith 
1896db781477SPatrick Sanan .seealso: `PetscSFComputeDegreeEnd()`, `PetscSFScatterEnd()`
189795fce210SBarry Smith @*/
189895fce210SBarry Smith PetscErrorCode PetscSFScatterEnd(PetscSF sf,MPI_Datatype unit,const void *multirootdata,void *leafdata)
189995fce210SBarry Smith {
1900a5526d50SJunchao Zhang   PetscSF        multi = NULL;
190195fce210SBarry Smith 
190295fce210SBarry Smith   PetscFunctionBegin;
190395fce210SBarry Smith   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
19049566063dSJacob Faibussowitsch   PetscCall(PetscSFGetMultiSF(sf,&multi));
19059566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(multi,unit,multirootdata,leafdata,MPI_REPLACE));
190695fce210SBarry Smith   PetscFunctionReturn(0);
190795fce210SBarry Smith }
1908a7b3aa13SAta Mesgarnejad 
1909a072220fSLawrence Mitchell static PetscErrorCode PetscSFCheckLeavesUnique_Private(PetscSF sf)
1910a072220fSLawrence Mitchell {
1911a072220fSLawrence Mitchell   PetscInt        i, n, nleaves;
1912a072220fSLawrence Mitchell   const PetscInt *ilocal = NULL;
1913a072220fSLawrence Mitchell   PetscHSetI      seen;
1914a072220fSLawrence Mitchell 
1915a072220fSLawrence Mitchell   PetscFunctionBegin;
1916b458e8f1SJose E. Roman   if (PetscDefined(USE_DEBUG)) {
19179566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sf,NULL,&nleaves,&ilocal,NULL));
19189566063dSJacob Faibussowitsch     PetscCall(PetscHSetICreate(&seen));
1919a072220fSLawrence Mitchell     for (i = 0; i < nleaves; i++) {
1920a072220fSLawrence Mitchell       const PetscInt leaf = ilocal ? ilocal[i] : i;
19219566063dSJacob Faibussowitsch       PetscCall(PetscHSetIAdd(seen,leaf));
1922a072220fSLawrence Mitchell     }
19239566063dSJacob Faibussowitsch     PetscCall(PetscHSetIGetSize(seen,&n));
192408401ef6SPierre Jolivet     PetscCheck(n == nleaves,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Provided leaves have repeated values: all leaves must be unique");
19259566063dSJacob Faibussowitsch     PetscCall(PetscHSetIDestroy(&seen));
1926b458e8f1SJose E. Roman   }
1927a072220fSLawrence Mitchell   PetscFunctionReturn(0);
1928a072220fSLawrence Mitchell }
192954729392SStefano Zampini 
1930a7b3aa13SAta Mesgarnejad /*@
193104c0ada0SJunchao Zhang   PetscSFCompose - Compose a new PetscSF by putting the second SF under the first one in a top (roots) down (leaves) view
1932a7b3aa13SAta Mesgarnejad 
1933a7b3aa13SAta Mesgarnejad   Input Parameters:
193454729392SStefano Zampini + sfA - The first PetscSF
193554729392SStefano Zampini - sfB - The second PetscSF
1936a7b3aa13SAta Mesgarnejad 
1937a7b3aa13SAta Mesgarnejad   Output Parameters:
193854729392SStefano Zampini . sfBA - The composite SF
1939a7b3aa13SAta Mesgarnejad 
1940a7b3aa13SAta Mesgarnejad   Level: developer
1941a7b3aa13SAta Mesgarnejad 
1942a072220fSLawrence Mitchell   Notes:
194354729392SStefano Zampini   Currently, the two SFs must be defined on congruent communicators and they must be true star
194454729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots.
194554729392SStefano Zampini 
194654729392SStefano Zampini   sfA's leaf space and sfB's root space might be partially overlapped. The composition builds
194754729392SStefano Zampini   a graph with sfA's roots and sfB's leaves only when there is a path between them. Unconnected
194854729392SStefano Zampini   nodes (roots or leaves) are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a
194954729392SStefano Zampini   Bcast on sfA, then a Bcast on sfB, on connected nodes.
1950a072220fSLawrence Mitchell 
1951db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFComposeInverse()`, `PetscSFGetGraph()`, `PetscSFSetGraph()`
1952a7b3aa13SAta Mesgarnejad @*/
1953a7b3aa13SAta Mesgarnejad PetscErrorCode PetscSFCompose(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
1954a7b3aa13SAta Mesgarnejad {
1955a7b3aa13SAta Mesgarnejad   const PetscSFNode *remotePointsA,*remotePointsB;
1956d41018fbSJunchao Zhang   PetscSFNode       *remotePointsBA=NULL,*reorderedRemotePointsA = NULL,*leafdataB;
195754729392SStefano Zampini   const PetscInt    *localPointsA,*localPointsB;
195854729392SStefano Zampini   PetscInt          *localPointsBA;
195954729392SStefano Zampini   PetscInt          i,numRootsA,numLeavesA,numRootsB,numLeavesB,minleaf,maxleaf,numLeavesBA;
196054729392SStefano Zampini   PetscBool         denseB;
1961a7b3aa13SAta Mesgarnejad 
1962a7b3aa13SAta Mesgarnejad   PetscFunctionBegin;
1963a7b3aa13SAta Mesgarnejad   PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1);
196429046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfA,1);
196529046d53SLisandro Dalcin   PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2);
196629046d53SLisandro Dalcin   PetscSFCheckGraphSet(sfB,2);
196754729392SStefano Zampini   PetscCheckSameComm(sfA,1,sfB,2);
196829046d53SLisandro Dalcin   PetscValidPointer(sfBA,3);
19699566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckLeavesUnique_Private(sfA));
19709566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckLeavesUnique_Private(sfB));
197154729392SStefano Zampini 
19729566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfA,&numRootsA,&numLeavesA,&localPointsA,&remotePointsA));
19739566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfB,&numRootsB,&numLeavesB,&localPointsB,&remotePointsB));
1974d41018fbSJunchao Zhang   if (localPointsA) {
19759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numRootsB,&reorderedRemotePointsA));
197654729392SStefano Zampini     for (i=0; i<numRootsB; i++) {
197754729392SStefano Zampini       reorderedRemotePointsA[i].rank = -1;
197854729392SStefano Zampini       reorderedRemotePointsA[i].index = -1;
197954729392SStefano Zampini     }
198054729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
198154729392SStefano Zampini       if (localPointsA[i] >= numRootsB) continue;
198254729392SStefano Zampini       reorderedRemotePointsA[localPointsA[i]] = remotePointsA[i];
198354729392SStefano Zampini     }
1984d41018fbSJunchao Zhang     remotePointsA = reorderedRemotePointsA;
1985d41018fbSJunchao Zhang   }
19869566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(sfB,&minleaf,&maxleaf));
19879566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxleaf-minleaf+1,&leafdataB));
19889566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf,MPI_REPLACE));
19899566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sfB,MPIU_2INT,remotePointsA,leafdataB-minleaf,MPI_REPLACE));
19909566063dSJacob Faibussowitsch   PetscCall(PetscFree(reorderedRemotePointsA));
1991d41018fbSJunchao Zhang 
199254729392SStefano Zampini   denseB = (PetscBool)!localPointsB;
199354729392SStefano Zampini   for (i=0,numLeavesBA=0; i<numLeavesB; i++) {
199454729392SStefano Zampini     if (leafdataB[localPointsB ? localPointsB[i]-minleaf : i].rank == -1) denseB = PETSC_FALSE;
199554729392SStefano Zampini     else numLeavesBA++;
199654729392SStefano Zampini   }
199754729392SStefano Zampini   if (denseB) {
1998d41018fbSJunchao Zhang     localPointsBA  = NULL;
1999d41018fbSJunchao Zhang     remotePointsBA = leafdataB;
2000d41018fbSJunchao Zhang   } else {
20019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numLeavesBA,&localPointsBA));
20029566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numLeavesBA,&remotePointsBA));
200354729392SStefano Zampini     for (i=0,numLeavesBA=0; i<numLeavesB; i++) {
200454729392SStefano Zampini       const PetscInt l = localPointsB ? localPointsB[i] : i;
200554729392SStefano Zampini 
200654729392SStefano Zampini       if (leafdataB[l-minleaf].rank == -1) continue;
200754729392SStefano Zampini       remotePointsBA[numLeavesBA] = leafdataB[l-minleaf];
200854729392SStefano Zampini       localPointsBA[numLeavesBA] = l;
200954729392SStefano Zampini       numLeavesBA++;
201054729392SStefano Zampini     }
20119566063dSJacob Faibussowitsch     PetscCall(PetscFree(leafdataB));
2012d41018fbSJunchao Zhang   }
20139566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA));
20149566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(*sfBA));
20159566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER));
2016a7b3aa13SAta Mesgarnejad   PetscFunctionReturn(0);
2017a7b3aa13SAta Mesgarnejad }
20181c6ba672SJunchao Zhang 
201904c0ada0SJunchao Zhang /*@
202054729392SStefano Zampini   PetscSFComposeInverse - Compose a new PetscSF by putting the inverse of the second SF under the first one
202104c0ada0SJunchao Zhang 
202204c0ada0SJunchao Zhang   Input Parameters:
202354729392SStefano Zampini + sfA - The first PetscSF
202454729392SStefano Zampini - sfB - The second PetscSF
202504c0ada0SJunchao Zhang 
202604c0ada0SJunchao Zhang   Output Parameters:
202754729392SStefano Zampini . sfBA - The composite SF.
202804c0ada0SJunchao Zhang 
202904c0ada0SJunchao Zhang   Level: developer
203004c0ada0SJunchao Zhang 
203154729392SStefano Zampini   Notes:
203254729392SStefano Zampini   Currently, the two SFs must be defined on congruent communicators and they must be true star
203354729392SStefano Zampini   forests, i.e. the same leaf is not connected with different roots. Even more, all roots of the
203454729392SStefano Zampini   second SF must have a degree of 1, i.e., no roots have more than one leaf connected.
203554729392SStefano Zampini 
203654729392SStefano Zampini   sfA's leaf space and sfB's leaf space might be partially overlapped. The composition builds
203754729392SStefano Zampini   a graph with sfA's roots and sfB's roots only when there is a path between them. Unconnected
203854729392SStefano Zampini   roots are not in sfBA. Doing a Bcast on the new SF is equivalent to doing a Bcast on sfA, then
203954729392SStefano Zampini   a Reduce on sfB, on connected roots.
204054729392SStefano Zampini 
2041db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFCompose()`, `PetscSFGetGraph()`, `PetscSFSetGraph()`, `PetscSFCreateInverseSF()`
204204c0ada0SJunchao Zhang @*/
204304c0ada0SJunchao Zhang PetscErrorCode PetscSFComposeInverse(PetscSF sfA,PetscSF sfB,PetscSF *sfBA)
204404c0ada0SJunchao Zhang {
204504c0ada0SJunchao Zhang   const PetscSFNode *remotePointsA,*remotePointsB;
204604c0ada0SJunchao Zhang   PetscSFNode       *remotePointsBA;
204704c0ada0SJunchao Zhang   const PetscInt    *localPointsA,*localPointsB;
204854729392SStefano Zampini   PetscSFNode       *reorderedRemotePointsA = NULL;
204954729392SStefano Zampini   PetscInt          i,numRootsA,numLeavesA,numLeavesBA,numRootsB,numLeavesB,minleaf,maxleaf,*localPointsBA;
20505b0d146aSStefano Zampini   MPI_Op            op;
20515b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
20525b0d146aSStefano Zampini   PetscBool         iswin;
20535b0d146aSStefano Zampini #endif
205404c0ada0SJunchao Zhang 
205504c0ada0SJunchao Zhang   PetscFunctionBegin;
205604c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfA,PETSCSF_CLASSID,1);
205704c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfA,1);
205804c0ada0SJunchao Zhang   PetscValidHeaderSpecific(sfB,PETSCSF_CLASSID,2);
205904c0ada0SJunchao Zhang   PetscSFCheckGraphSet(sfB,2);
206054729392SStefano Zampini   PetscCheckSameComm(sfA,1,sfB,2);
206104c0ada0SJunchao Zhang   PetscValidPointer(sfBA,3);
20629566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckLeavesUnique_Private(sfA));
20639566063dSJacob Faibussowitsch   PetscCall(PetscSFCheckLeavesUnique_Private(sfB));
206454729392SStefano Zampini 
20659566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA));
20669566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB));
20675b0d146aSStefano Zampini 
20685b0d146aSStefano Zampini   /* TODO: Check roots of sfB have degree of 1 */
20695b0d146aSStefano Zampini   /* Once we implement it, we can replace the MPI_MAXLOC
207083df288dSJunchao Zhang      with MPI_REPLACE. In that case, MPI_MAXLOC and MPI_REPLACE have the same effect.
20715b0d146aSStefano Zampini      We use MPI_MAXLOC only to have a deterministic output from this routine if
20725b0d146aSStefano Zampini      the root condition is not meet.
20735b0d146aSStefano Zampini    */
20745b0d146aSStefano Zampini   op = MPI_MAXLOC;
20755b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES)
20765b0d146aSStefano Zampini   /* we accept a non-deterministic output (if any) with PETSCSFWINDOW, since MPI_MAXLOC cannot operate on MPIU_2INT with MPI_Accumulate */
20779566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)sfB,PETSCSFWINDOW,&iswin));
207883df288dSJunchao Zhang   if (iswin) op = MPI_REPLACE;
20795b0d146aSStefano Zampini #endif
20805b0d146aSStefano Zampini 
20819566063dSJacob Faibussowitsch   PetscCall(PetscSFGetLeafRange(sfB, &minleaf, &maxleaf));
20829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxleaf - minleaf + 1,&reorderedRemotePointsA));
208354729392SStefano Zampini   for (i=0; i<maxleaf - minleaf + 1; i++) {
208454729392SStefano Zampini     reorderedRemotePointsA[i].rank = -1;
208554729392SStefano Zampini     reorderedRemotePointsA[i].index = -1;
208654729392SStefano Zampini   }
208754729392SStefano Zampini   if (localPointsA) {
208854729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
208954729392SStefano Zampini       if (localPointsA[i] > maxleaf || localPointsA[i] < minleaf) continue;
209054729392SStefano Zampini       reorderedRemotePointsA[localPointsA[i] - minleaf] = remotePointsA[i];
209154729392SStefano Zampini     }
209254729392SStefano Zampini   } else {
209354729392SStefano Zampini     for (i=0; i<numLeavesA; i++) {
209454729392SStefano Zampini       if (i > maxleaf || i < minleaf) continue;
209554729392SStefano Zampini       reorderedRemotePointsA[i - minleaf] = remotePointsA[i];
209654729392SStefano Zampini     }
209754729392SStefano Zampini   }
209854729392SStefano Zampini 
20999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numRootsB,&localPointsBA));
21009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numRootsB,&remotePointsBA));
210154729392SStefano Zampini   for (i=0; i<numRootsB; i++) {
210254729392SStefano Zampini     remotePointsBA[i].rank = -1;
210354729392SStefano Zampini     remotePointsBA[i].index = -1;
210454729392SStefano Zampini   }
210554729392SStefano Zampini 
21069566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,op));
21079566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sfB,MPIU_2INT,reorderedRemotePointsA-minleaf,remotePointsBA,op));
21089566063dSJacob Faibussowitsch   PetscCall(PetscFree(reorderedRemotePointsA));
210954729392SStefano Zampini   for (i=0,numLeavesBA=0; i<numRootsB; i++) {
211054729392SStefano Zampini     if (remotePointsBA[i].rank == -1) continue;
211154729392SStefano Zampini     remotePointsBA[numLeavesBA].rank = remotePointsBA[i].rank;
211254729392SStefano Zampini     remotePointsBA[numLeavesBA].index = remotePointsBA[i].index;
211354729392SStefano Zampini     localPointsBA[numLeavesBA] = i;
211454729392SStefano Zampini     numLeavesBA++;
211554729392SStefano Zampini   }
21169566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfA),sfBA));
21179566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(*sfBA));
21189566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*sfBA,numRootsA,numLeavesBA,localPointsBA,PETSC_OWN_POINTER,remotePointsBA,PETSC_OWN_POINTER));
211904c0ada0SJunchao Zhang   PetscFunctionReturn(0);
212004c0ada0SJunchao Zhang }
212104c0ada0SJunchao Zhang 
21221c6ba672SJunchao Zhang /*
21231c6ba672SJunchao Zhang   PetscSFCreateLocalSF_Private - Creates a local PetscSF that only has intra-process edges of the global PetscSF
21241c6ba672SJunchao Zhang 
21251c6ba672SJunchao Zhang   Input Parameters:
21261c6ba672SJunchao Zhang . sf - The global PetscSF
21271c6ba672SJunchao Zhang 
21281c6ba672SJunchao Zhang   Output Parameters:
21291c6ba672SJunchao Zhang . out - The local PetscSF
21301c6ba672SJunchao Zhang  */
21311c6ba672SJunchao Zhang PetscErrorCode PetscSFCreateLocalSF_Private(PetscSF sf,PetscSF *out)
21321c6ba672SJunchao Zhang {
21331c6ba672SJunchao Zhang   MPI_Comm           comm;
21341c6ba672SJunchao Zhang   PetscMPIInt        myrank;
21351c6ba672SJunchao Zhang   const PetscInt     *ilocal;
21361c6ba672SJunchao Zhang   const PetscSFNode  *iremote;
21371c6ba672SJunchao Zhang   PetscInt           i,j,nroots,nleaves,lnleaves,*lilocal;
21381c6ba672SJunchao Zhang   PetscSFNode        *liremote;
21391c6ba672SJunchao Zhang   PetscSF            lsf;
21401c6ba672SJunchao Zhang 
21411c6ba672SJunchao Zhang   PetscFunctionBegin;
21421c6ba672SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
2143*dbbe0bcdSBarry Smith   if (sf->ops->CreateLocalSF) PetscUseTypeMethod(sf,CreateLocalSF ,out);
2144*dbbe0bcdSBarry Smith   else {
21451c6ba672SJunchao Zhang     /* Could use PetscSFCreateEmbeddedLeafSF, but since we know the comm is PETSC_COMM_SELF, we can make it fast */
21469566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)sf,&comm));
21479566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm,&myrank));
21481c6ba672SJunchao Zhang 
21491c6ba672SJunchao Zhang     /* Find out local edges and build a local SF */
21509566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote));
21511c6ba672SJunchao Zhang     for (i=lnleaves=0; i<nleaves; i++) {if (iremote[i].rank == (PetscInt)myrank) lnleaves++;}
21529566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lnleaves,&lilocal));
21539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lnleaves,&liremote));
21541c6ba672SJunchao Zhang 
21551c6ba672SJunchao Zhang     for (i=j=0; i<nleaves; i++) {
21561c6ba672SJunchao Zhang       if (iremote[i].rank == (PetscInt)myrank) {
21571c6ba672SJunchao Zhang         lilocal[j]        = ilocal? ilocal[i] : i; /* ilocal=NULL for contiguous storage */
21581c6ba672SJunchao Zhang         liremote[j].rank  = 0; /* rank in PETSC_COMM_SELF */
21591c6ba672SJunchao Zhang         liremote[j].index = iremote[i].index;
21601c6ba672SJunchao Zhang         j++;
21611c6ba672SJunchao Zhang       }
21621c6ba672SJunchao Zhang     }
21639566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PETSC_COMM_SELF,&lsf));
21649566063dSJacob Faibussowitsch     PetscCall(PetscSFSetFromOptions(lsf));
21659566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(lsf,nroots,lnleaves,lilocal,PETSC_OWN_POINTER,liremote,PETSC_OWN_POINTER));
21669566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(lsf));
21671c6ba672SJunchao Zhang     *out = lsf;
21681c6ba672SJunchao Zhang   }
21691c6ba672SJunchao Zhang   PetscFunctionReturn(0);
21701c6ba672SJunchao Zhang }
2171dd5b3ca6SJunchao Zhang 
2172dd5b3ca6SJunchao Zhang /* Similar to PetscSFBcast, but only Bcast to leaves on rank 0 */
2173dd5b3ca6SJunchao Zhang PetscErrorCode PetscSFBcastToZero_Private(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata)
2174dd5b3ca6SJunchao Zhang {
2175eb02082bSJunchao Zhang   PetscMemType       rootmtype,leafmtype;
2176dd5b3ca6SJunchao Zhang 
2177dd5b3ca6SJunchao Zhang   PetscFunctionBegin;
2178dd5b3ca6SJunchao Zhang   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
21799566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
21809566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCSF_BcastBegin,sf,0,0,0));
21819566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(rootdata,&rootmtype));
21829566063dSJacob Faibussowitsch   PetscCall(PetscGetMemType(leafdata,&leafmtype));
2183*dbbe0bcdSBarry Smith   PetscUseTypeMethod(sf,BcastToZero ,unit,rootmtype,rootdata,leafmtype,leafdata);
21849566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCSF_BcastBegin,sf,0,0,0));
2185dd5b3ca6SJunchao Zhang   PetscFunctionReturn(0);
2186dd5b3ca6SJunchao Zhang }
2187dd5b3ca6SJunchao Zhang 
2188157edd7aSVaclav Hapla /*@
2189157edd7aSVaclav Hapla   PetscSFConcatenate - concatenate multiple SFs into one
2190157edd7aSVaclav Hapla 
2191157edd7aSVaclav Hapla   Input Parameters:
2192157edd7aSVaclav Hapla + comm - the communicator
2193157edd7aSVaclav Hapla . nsfs - the number of input PetscSF
2194157edd7aSVaclav Hapla . sfs  - the array of input PetscSF
2195157edd7aSVaclav Hapla . shareRoots - the flag whether roots of input PetscSFs are taken as shared (PETSC_TRUE), or separate and concatenated (PETSC_FALSE)
2196157edd7aSVaclav Hapla - leafOffsets - the array of local leaf offsets, one for each input PetscSF, or NULL for contiguous storage
2197157edd7aSVaclav Hapla 
2198157edd7aSVaclav Hapla   Output Parameters:
2199157edd7aSVaclav Hapla . newsf - The resulting PetscSF
2200157edd7aSVaclav Hapla 
2201157edd7aSVaclav Hapla   Level: developer
2202157edd7aSVaclav Hapla 
2203157edd7aSVaclav Hapla   Notes:
2204157edd7aSVaclav Hapla   The communicator of all SFs in sfs must be comm.
2205157edd7aSVaclav Hapla 
2206157edd7aSVaclav Hapla   The offsets in leafOffsets are added to the original leaf indices.
2207157edd7aSVaclav Hapla 
2208157edd7aSVaclav Hapla   If all input SFs use contiguous leaf storage (ilocal = NULL), leafOffsets can be passed as NULL as well.
2209157edd7aSVaclav Hapla   In this case, NULL is also passed as ilocal to the resulting SF.
2210157edd7aSVaclav Hapla 
2211157edd7aSVaclav Hapla   If any input SF has non-null ilocal, leafOffsets is needed to distinguish leaves from different input SFs.
2212157edd7aSVaclav Hapla   In this case, user is responsible to provide correct offsets so that the resulting leaves are unique (otherwise an error occurs).
2213157edd7aSVaclav Hapla 
2214db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFCompose()`, `PetscSFGetGraph()`, `PetscSFSetGraph()`
2215157edd7aSVaclav Hapla @*/
2216157edd7aSVaclav Hapla PetscErrorCode PetscSFConcatenate(MPI_Comm comm, PetscInt nsfs, PetscSF sfs[], PetscBool shareRoots, PetscInt leafOffsets[], PetscSF *newsf)
2217157edd7aSVaclav Hapla {
2218157edd7aSVaclav Hapla   PetscInt            i, s, nLeaves, nRoots;
2219157edd7aSVaclav Hapla   PetscInt           *leafArrayOffsets;
2220157edd7aSVaclav Hapla   PetscInt           *ilocal_new;
2221157edd7aSVaclav Hapla   PetscSFNode        *iremote_new;
2222157edd7aSVaclav Hapla   PetscInt           *rootOffsets;
2223157edd7aSVaclav Hapla   PetscBool           all_ilocal_null = PETSC_FALSE;
2224157edd7aSVaclav Hapla   PetscMPIInt         rank;
2225157edd7aSVaclav Hapla 
2226157edd7aSVaclav Hapla   PetscFunctionBegin;
2227157edd7aSVaclav Hapla   {
2228157edd7aSVaclav Hapla     PetscSF dummy; /* just to have a PetscObject on comm for input validation */
2229157edd7aSVaclav Hapla 
22309566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(comm,&dummy));
2231157edd7aSVaclav Hapla     PetscValidLogicalCollectiveInt(dummy,nsfs,2);
2232157edd7aSVaclav Hapla     PetscValidPointer(sfs,3);
2233157edd7aSVaclav Hapla     for (i=0; i<nsfs; i++) {
2234157edd7aSVaclav Hapla       PetscValidHeaderSpecific(sfs[i],PETSCSF_CLASSID,3);
2235157edd7aSVaclav Hapla       PetscCheckSameComm(dummy,1,sfs[i],3);
2236157edd7aSVaclav Hapla     }
2237157edd7aSVaclav Hapla     PetscValidLogicalCollectiveBool(dummy,shareRoots,4);
2238157edd7aSVaclav Hapla     if (leafOffsets) PetscValidIntPointer(leafOffsets,5);
2239157edd7aSVaclav Hapla     PetscValidPointer(newsf,6);
22409566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&dummy));
2241157edd7aSVaclav Hapla   }
2242157edd7aSVaclav Hapla   if (!nsfs) {
22439566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(comm, newsf));
22449566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*newsf, 0, 0, NULL, PETSC_OWN_POINTER, NULL, PETSC_OWN_POINTER));
2245157edd7aSVaclav Hapla     PetscFunctionReturn(0);
2246157edd7aSVaclav Hapla   }
22479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2248157edd7aSVaclav Hapla 
22499566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nsfs+1, &rootOffsets));
2250157edd7aSVaclav Hapla   if (shareRoots) {
22519566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfs[0], &nRoots, NULL, NULL, NULL));
2252157edd7aSVaclav Hapla     if (PetscDefined(USE_DEBUG)) {
2253157edd7aSVaclav Hapla       for (s=1; s<nsfs; s++) {
2254157edd7aSVaclav Hapla         PetscInt nr;
2255157edd7aSVaclav Hapla 
22569566063dSJacob Faibussowitsch         PetscCall(PetscSFGetGraph(sfs[s], &nr, NULL, NULL, NULL));
2257157edd7aSVaclav Hapla         PetscCheck(nr == nRoots, comm, PETSC_ERR_ARG_SIZ, "shareRoots = PETSC_TRUE but sfs[%" PetscInt_FMT "] has a different number of roots (%" PetscInt_FMT ") than sfs[0] (%" PetscInt_FMT ")", s, nr, nRoots);
2258157edd7aSVaclav Hapla       }
2259157edd7aSVaclav Hapla     }
2260157edd7aSVaclav Hapla   } else {
2261157edd7aSVaclav Hapla     for (s=0; s<nsfs; s++) {
2262157edd7aSVaclav Hapla       PetscInt nr;
2263157edd7aSVaclav Hapla 
22649566063dSJacob Faibussowitsch       PetscCall(PetscSFGetGraph(sfs[s], &nr, NULL, NULL, NULL));
2265157edd7aSVaclav Hapla       rootOffsets[s+1] = rootOffsets[s] + nr;
2266157edd7aSVaclav Hapla     }
2267157edd7aSVaclav Hapla     nRoots = rootOffsets[nsfs];
2268157edd7aSVaclav Hapla   }
2269157edd7aSVaclav Hapla 
2270157edd7aSVaclav Hapla   /* Calculate leaf array offsets and automatic root offsets */
22719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsfs+1,&leafArrayOffsets));
2272157edd7aSVaclav Hapla   leafArrayOffsets[0] = 0;
2273157edd7aSVaclav Hapla   for (s=0; s<nsfs; s++) {
2274157edd7aSVaclav Hapla     PetscInt        nl;
2275157edd7aSVaclav Hapla 
22769566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfs[s], NULL, &nl, NULL, NULL));
2277157edd7aSVaclav Hapla     leafArrayOffsets[s+1] = leafArrayOffsets[s] + nl;
2278157edd7aSVaclav Hapla   }
2279157edd7aSVaclav Hapla   nLeaves = leafArrayOffsets[nsfs];
2280157edd7aSVaclav Hapla 
2281157edd7aSVaclav Hapla   if (!leafOffsets) {
2282157edd7aSVaclav Hapla     all_ilocal_null = PETSC_TRUE;
2283157edd7aSVaclav Hapla     for (s=0; s<nsfs; s++) {
2284157edd7aSVaclav Hapla       const PetscInt *ilocal;
2285157edd7aSVaclav Hapla 
22869566063dSJacob Faibussowitsch       PetscCall(PetscSFGetGraph(sfs[s], NULL, NULL, &ilocal, NULL));
2287157edd7aSVaclav Hapla       if (ilocal) {
2288157edd7aSVaclav Hapla         all_ilocal_null = PETSC_FALSE;
2289157edd7aSVaclav Hapla         break;
2290157edd7aSVaclav Hapla       }
2291157edd7aSVaclav Hapla     }
2292157edd7aSVaclav 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");
2293157edd7aSVaclav Hapla   }
2294157edd7aSVaclav Hapla 
2295157edd7aSVaclav Hapla   /* Renumber and concatenate local leaves */
2296157edd7aSVaclav Hapla   ilocal_new = NULL;
2297157edd7aSVaclav Hapla   if (!all_ilocal_null) {
22989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nLeaves, &ilocal_new));
2299157edd7aSVaclav Hapla     for (i = 0; i<nLeaves; i++) ilocal_new[i] = -1;
2300157edd7aSVaclav Hapla     for (s = 0; s<nsfs; s++) {
2301157edd7aSVaclav Hapla       const PetscInt   *ilocal;
2302157edd7aSVaclav Hapla       PetscInt         *ilocal_l = &ilocal_new[leafArrayOffsets[s]];
2303157edd7aSVaclav Hapla       PetscInt          i, nleaves_l;
2304157edd7aSVaclav Hapla 
23059566063dSJacob Faibussowitsch       PetscCall(PetscSFGetGraph(sfs[s], NULL, &nleaves_l, &ilocal, NULL));
2306157edd7aSVaclav Hapla       for (i=0; i<nleaves_l; i++) ilocal_l[i] = (ilocal ? ilocal[i] : i) + leafOffsets[s];
2307157edd7aSVaclav Hapla     }
2308157edd7aSVaclav Hapla   }
2309157edd7aSVaclav Hapla 
2310157edd7aSVaclav Hapla   /* Renumber and concatenate remote roots */
23119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nLeaves, &iremote_new));
2312157edd7aSVaclav Hapla   for (i = 0; i < nLeaves; i++) {
2313157edd7aSVaclav Hapla     iremote_new[i].rank   = -1;
2314157edd7aSVaclav Hapla     iremote_new[i].index  = -1;
2315157edd7aSVaclav Hapla   }
2316157edd7aSVaclav Hapla   for (s = 0; s<nsfs; s++) {
2317157edd7aSVaclav Hapla     PetscInt            i, nl, nr;
2318157edd7aSVaclav Hapla     PetscSF             tmp_sf;
2319157edd7aSVaclav Hapla     const PetscSFNode  *iremote;
2320157edd7aSVaclav Hapla     PetscSFNode        *tmp_rootdata;
2321157edd7aSVaclav Hapla     PetscSFNode        *tmp_leafdata = &iremote_new[leafArrayOffsets[s]];
2322157edd7aSVaclav Hapla 
23239566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfs[s], &nr, &nl, NULL, &iremote));
23249566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(comm, &tmp_sf));
2325157edd7aSVaclav Hapla     /* create helper SF with contiguous leaves */
23269566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(tmp_sf, nr, nl, NULL, PETSC_USE_POINTER, (PetscSFNode*) iremote, PETSC_COPY_VALUES));
23279566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(tmp_sf));
23289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nr, &tmp_rootdata));
2329157edd7aSVaclav Hapla     for (i = 0; i < nr; i++) {
2330157edd7aSVaclav Hapla       tmp_rootdata[i].index = i + rootOffsets[s];
2331157edd7aSVaclav Hapla       tmp_rootdata[i].rank  = (PetscInt) rank;
2332157edd7aSVaclav Hapla     }
23339566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(tmp_sf, MPIU_2INT, tmp_rootdata, tmp_leafdata, MPI_REPLACE));
23349566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(  tmp_sf, MPIU_2INT, tmp_rootdata, tmp_leafdata, MPI_REPLACE));
23359566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&tmp_sf));
23369566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp_rootdata));
2337157edd7aSVaclav Hapla   }
2338157edd7aSVaclav Hapla 
2339157edd7aSVaclav Hapla   /* Build the new SF */
23409566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, newsf));
23419566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*newsf, nRoots, nLeaves, ilocal_new, PETSC_OWN_POINTER, iremote_new, PETSC_OWN_POINTER));
23429566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(*newsf));
23439566063dSJacob Faibussowitsch   PetscCall(PetscFree(rootOffsets));
23449566063dSJacob Faibussowitsch   PetscCall(PetscFree(leafArrayOffsets));
2345157edd7aSVaclav Hapla   PetscFunctionReturn(0);
2346157edd7aSVaclav Hapla }
2347