1af0996ceSBarry Smith #include <petsc/private/sfimpl.h> /*I "petscsf.h" I*/ 2c4e6a40aSLawrence Mitchell #include <petsc/private/hashseti.h> 353dd6d7dSJunchao Zhang #include <petsc/private/viewerimpl.h> 4eec179cfSJacob Faibussowitsch #include <petsc/private/hashmapi.h> 595fce210SBarry Smith 67fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_CUDA) 77fd2d3dbSJunchao Zhang #include <cuda_runtime.h> 87fd2d3dbSJunchao Zhang #endif 97fd2d3dbSJunchao Zhang 107fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_HIP) 117fd2d3dbSJunchao Zhang #include <hip/hip_runtime.h> 127fd2d3dbSJunchao Zhang #endif 137fd2d3dbSJunchao Zhang 142abc8c78SJacob Faibussowitsch #if defined(PETSC_CLANG_STATIC_ANALYZER) 152abc8c78SJacob Faibussowitsch void PetscSFCheckGraphSet(PetscSF, int); 162abc8c78SJacob Faibussowitsch #else 1795fce210SBarry Smith #if defined(PETSC_USE_DEBUG) 187a46b595SBarry Smith #define PetscSFCheckGraphSet(sf, arg) PetscCheck((sf)->graphset, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFSetGraph() or PetscSFSetGraphWithPattern() on argument %d \"%s\" before %s()", (arg), #sf, PETSC_FUNCTION_NAME); 1995fce210SBarry Smith #else 209371c9d4SSatish Balay #define PetscSFCheckGraphSet(sf, arg) \ 219371c9d4SSatish Balay do { \ 229371c9d4SSatish Balay } while (0) 2395fce210SBarry Smith #endif 242abc8c78SJacob Faibussowitsch #endif 2595fce210SBarry Smith 264c8fdceaSLisandro Dalcin const char *const PetscSFDuplicateOptions[] = {"CONFONLY", "RANKS", "GRAPH", "PetscSFDuplicateOption", "PETSCSF_DUPLICATE_", NULL}; 271f40158dSVaclav Hapla const char *const PetscSFConcatenateRootModes[] = {"local", "shared", "global", "PetscSFConcatenateRootMode", "PETSCSF_CONCATENATE_ROOTMODE_", NULL}; 2895fce210SBarry Smith 298af6ec1cSBarry Smith /*@ 3095fce210SBarry Smith PetscSFCreate - create a star forest communication context 3195fce210SBarry Smith 32d083f849SBarry Smith Collective 3395fce210SBarry Smith 344165533cSJose E. Roman Input Parameter: 3595fce210SBarry Smith . comm - communicator on which the star forest will operate 3695fce210SBarry Smith 374165533cSJose E. Roman Output Parameter: 3895fce210SBarry Smith . sf - new star forest context 3995fce210SBarry Smith 4020662ed9SBarry Smith Options Database Key: 41cab54364SBarry Smith . -sf_type type - value of type may be 42cab54364SBarry Smith .vb 43cab54364SBarry Smith basic -Use MPI persistent Isend/Irecv for communication (Default) 44cab54364SBarry Smith window -Use MPI-3 one-sided window for communication 45cab54364SBarry Smith neighbor -Use MPI-3 neighborhood collectives for communication 46cab54364SBarry Smith .ve 47dd5b3ca6SJunchao Zhang 4895fce210SBarry Smith Level: intermediate 4995fce210SBarry Smith 50cab54364SBarry Smith Note: 51cab54364SBarry Smith When one knows the communication graph is one of the predefined graph, such as `MPI_Alltoall()`, `MPI_Allgatherv()`, 52cab54364SBarry Smith `MPI_Gatherv()`, one can create a `PetscSF` and then set its graph with `PetscSFSetGraphWithPattern()`. These special 5320662ed9SBarry Smith `SF`s are optimized and they have better performance than the general `SF`s. 54dd5b3ca6SJunchao Zhang 5520662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFSetType`, PetscSFSetGraph()`, `PetscSFSetGraphWithPattern()`, `PetscSFDestroy()` 5695fce210SBarry Smith @*/ 57d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreate(MPI_Comm comm, PetscSF *sf) 58d71ae5a4SJacob Faibussowitsch { 5995fce210SBarry Smith PetscSF b; 6095fce210SBarry Smith 6195fce210SBarry Smith PetscFunctionBegin; 6295fce210SBarry Smith PetscValidPointer(sf, 2); 639566063dSJacob Faibussowitsch PetscCall(PetscSFInitializePackage()); 6495fce210SBarry Smith 659566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(b, PETSCSF_CLASSID, "PetscSF", "Star Forest", "PetscSF", comm, PetscSFDestroy, PetscSFView)); 6695fce210SBarry Smith 6795fce210SBarry Smith b->nroots = -1; 6895fce210SBarry Smith b->nleaves = -1; 6929046d53SLisandro Dalcin b->minleaf = PETSC_MAX_INT; 7029046d53SLisandro Dalcin b->maxleaf = PETSC_MIN_INT; 7195fce210SBarry Smith b->nranks = -1; 7295fce210SBarry Smith b->rankorder = PETSC_TRUE; 7395fce210SBarry Smith b->ingroup = MPI_GROUP_NULL; 7495fce210SBarry Smith b->outgroup = MPI_GROUP_NULL; 7595fce210SBarry Smith b->graphset = PETSC_FALSE; 7620c24465SJunchao Zhang #if defined(PETSC_HAVE_DEVICE) 7720c24465SJunchao Zhang b->use_gpu_aware_mpi = use_gpu_aware_mpi; 7820c24465SJunchao Zhang b->use_stream_aware_mpi = PETSC_FALSE; 7971438e86SJunchao Zhang b->unknown_input_stream = PETSC_FALSE; 8027f636e8SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS) /* Prefer kokkos over cuda*/ 8120c24465SJunchao Zhang b->backend = PETSCSF_BACKEND_KOKKOS; 8227f636e8SJunchao Zhang #elif defined(PETSC_HAVE_CUDA) 8327f636e8SJunchao Zhang b->backend = PETSCSF_BACKEND_CUDA; 8459af0bd3SScott Kruger #elif defined(PETSC_HAVE_HIP) 8559af0bd3SScott Kruger b->backend = PETSCSF_BACKEND_HIP; 8620c24465SJunchao Zhang #endif 8771438e86SJunchao Zhang 8871438e86SJunchao Zhang #if defined(PETSC_HAVE_NVSHMEM) 8971438e86SJunchao Zhang b->use_nvshmem = PETSC_FALSE; /* Default is not to try NVSHMEM */ 9071438e86SJunchao Zhang b->use_nvshmem_get = PETSC_FALSE; /* Default is to use nvshmem_put based protocol */ 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_nvshmem", &b->use_nvshmem, NULL)); 929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_nvshmem_get", &b->use_nvshmem_get, NULL)); 9371438e86SJunchao Zhang #endif 9420c24465SJunchao Zhang #endif 9560c22052SBarry Smith b->vscat.from_n = -1; 9660c22052SBarry Smith b->vscat.to_n = -1; 9760c22052SBarry Smith b->vscat.unit = MPIU_SCALAR; 9895fce210SBarry Smith *sf = b; 993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10095fce210SBarry Smith } 10195fce210SBarry Smith 10229046d53SLisandro Dalcin /*@ 10395fce210SBarry Smith PetscSFReset - Reset a star forest so that different sizes or neighbors can be used 10495fce210SBarry Smith 10595fce210SBarry Smith Collective 10695fce210SBarry Smith 1074165533cSJose E. Roman Input Parameter: 10895fce210SBarry Smith . sf - star forest 10995fce210SBarry Smith 11095fce210SBarry Smith Level: advanced 11195fce210SBarry Smith 112cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`, `PetscSFSetGraph()`, `PetscSFDestroy()` 11395fce210SBarry Smith @*/ 114d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFReset(PetscSF sf) 115d71ae5a4SJacob Faibussowitsch { 11695fce210SBarry Smith PetscFunctionBegin; 11795fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 118dbbe0bcdSBarry Smith PetscTryTypeMethod(sf, Reset); 11929046d53SLisandro Dalcin sf->nroots = -1; 12029046d53SLisandro Dalcin sf->nleaves = -1; 12129046d53SLisandro Dalcin sf->minleaf = PETSC_MAX_INT; 12229046d53SLisandro Dalcin sf->maxleaf = PETSC_MIN_INT; 12395fce210SBarry Smith sf->mine = NULL; 12495fce210SBarry Smith sf->remote = NULL; 12529046d53SLisandro Dalcin sf->graphset = PETSC_FALSE; 1269566063dSJacob Faibussowitsch PetscCall(PetscFree(sf->mine_alloc)); 1279566063dSJacob Faibussowitsch PetscCall(PetscFree(sf->remote_alloc)); 12821c688dcSJed Brown sf->nranks = -1; 1299566063dSJacob Faibussowitsch PetscCall(PetscFree4(sf->ranks, sf->roffset, sf->rmine, sf->rremote)); 13029046d53SLisandro Dalcin sf->degreeknown = PETSC_FALSE; 1319566063dSJacob Faibussowitsch PetscCall(PetscFree(sf->degree)); 1329566063dSJacob Faibussowitsch if (sf->ingroup != MPI_GROUP_NULL) PetscCallMPI(MPI_Group_free(&sf->ingroup)); 1339566063dSJacob Faibussowitsch if (sf->outgroup != MPI_GROUP_NULL) PetscCallMPI(MPI_Group_free(&sf->outgroup)); 134013b3241SStefano Zampini if (sf->multi) sf->multi->multi = NULL; 1359566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf->multi)); 1369566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&sf->map)); 13771438e86SJunchao Zhang 13871438e86SJunchao Zhang #if defined(PETSC_HAVE_DEVICE) 1399566063dSJacob Faibussowitsch for (PetscInt i = 0; i < 2; i++) PetscCall(PetscSFFree(sf, PETSC_MEMTYPE_DEVICE, sf->rmine_d[i])); 14071438e86SJunchao Zhang #endif 14171438e86SJunchao Zhang 14295fce210SBarry Smith sf->setupcalled = PETSC_FALSE; 1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14495fce210SBarry Smith } 14595fce210SBarry Smith 14695fce210SBarry Smith /*@C 147cab54364SBarry Smith PetscSFSetType - Set the `PetscSF` communication implementation 14895fce210SBarry Smith 149c3339decSBarry Smith Collective 15095fce210SBarry Smith 15195fce210SBarry Smith Input Parameters: 152cab54364SBarry Smith + sf - the `PetscSF` context 15395fce210SBarry Smith - type - a known method 154cab54364SBarry Smith .vb 155cab54364SBarry Smith PETSCSFWINDOW - MPI-2/3 one-sided 156cab54364SBarry Smith PETSCSFBASIC - basic implementation using MPI-1 two-sided 157cab54364SBarry Smith .ve 15895fce210SBarry Smith 15995fce210SBarry Smith Options Database Key: 16020662ed9SBarry Smith . -sf_type <type> - Sets the method; for example `basic` or `window` use -help for a list of available methods 161cab54364SBarry Smith 162cab54364SBarry Smith Level: intermediate 16395fce210SBarry Smith 16495fce210SBarry Smith Notes: 16520662ed9SBarry Smith See `PetscSFType` for possible values 16695fce210SBarry Smith 16720662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFCreate()` 16895fce210SBarry Smith @*/ 169d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetType(PetscSF sf, PetscSFType type) 170d71ae5a4SJacob Faibussowitsch { 17195fce210SBarry Smith PetscBool match; 1725f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(PetscSF); 17395fce210SBarry Smith 17495fce210SBarry Smith PetscFunctionBegin; 17595fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 17695fce210SBarry Smith PetscValidCharPointer(type, 2); 17795fce210SBarry Smith 1789566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)sf, type, &match)); 1793ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 18095fce210SBarry Smith 1819566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(PetscSFList, type, &r)); 18228b400f6SJacob Faibussowitsch PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested PetscSF type %s", type); 18329046d53SLisandro Dalcin /* Destroy the previous PetscSF implementation context */ 184dbbe0bcdSBarry Smith PetscTryTypeMethod(sf, Destroy); 1859566063dSJacob Faibussowitsch PetscCall(PetscMemzero(sf->ops, sizeof(*sf->ops))); 1869566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)sf, type)); 1879566063dSJacob Faibussowitsch PetscCall((*r)(sf)); 1883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18995fce210SBarry Smith } 19095fce210SBarry Smith 19129046d53SLisandro Dalcin /*@C 192cab54364SBarry Smith PetscSFGetType - Get the `PetscSF` communication implementation 19329046d53SLisandro Dalcin 19429046d53SLisandro Dalcin Not Collective 19529046d53SLisandro Dalcin 19629046d53SLisandro Dalcin Input Parameter: 197cab54364SBarry Smith . sf - the `PetscSF` context 19829046d53SLisandro Dalcin 19929046d53SLisandro Dalcin Output Parameter: 200cab54364SBarry Smith . type - the `PetscSF` type name 20129046d53SLisandro Dalcin 20229046d53SLisandro Dalcin Level: intermediate 20329046d53SLisandro Dalcin 20420662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFSetType()`, `PetscSFCreate()` 20529046d53SLisandro Dalcin @*/ 206d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetType(PetscSF sf, PetscSFType *type) 207d71ae5a4SJacob Faibussowitsch { 20829046d53SLisandro Dalcin PetscFunctionBegin; 20929046d53SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 21029046d53SLisandro Dalcin PetscValidPointer(type, 2); 21129046d53SLisandro Dalcin *type = ((PetscObject)sf)->type_name; 2123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21329046d53SLisandro Dalcin } 21429046d53SLisandro Dalcin 2151fb7b255SJunchao Zhang /*@C 21620662ed9SBarry Smith PetscSFDestroy - destroy a star forest 21795fce210SBarry Smith 21895fce210SBarry Smith Collective 21995fce210SBarry Smith 2204165533cSJose E. Roman Input Parameter: 22195fce210SBarry Smith . sf - address of star forest 22295fce210SBarry Smith 22395fce210SBarry Smith Level: intermediate 22495fce210SBarry Smith 22520662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFCreate()`, `PetscSFReset()` 22695fce210SBarry Smith @*/ 227d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFDestroy(PetscSF *sf) 228d71ae5a4SJacob Faibussowitsch { 22995fce210SBarry Smith PetscFunctionBegin; 2303ba16761SJacob Faibussowitsch if (!*sf) PetscFunctionReturn(PETSC_SUCCESS); 23195fce210SBarry Smith PetscValidHeaderSpecific((*sf), PETSCSF_CLASSID, 1); 2329371c9d4SSatish Balay if (--((PetscObject)(*sf))->refct > 0) { 2339371c9d4SSatish Balay *sf = NULL; 2343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2359371c9d4SSatish Balay } 2369566063dSJacob Faibussowitsch PetscCall(PetscSFReset(*sf)); 237dbbe0bcdSBarry Smith PetscTryTypeMethod((*sf), Destroy); 2389566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&(*sf)->vscat.lsf)); 2399566063dSJacob Faibussowitsch if ((*sf)->vscat.bs > 1) PetscCallMPI(MPI_Type_free(&(*sf)->vscat.unit)); 2409566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(sf)); 2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24295fce210SBarry Smith } 24395fce210SBarry Smith 244d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFCheckGraphValid_Private(PetscSF sf) 245d71ae5a4SJacob Faibussowitsch { 246c4e6a40aSLawrence Mitchell PetscInt i, nleaves; 247c4e6a40aSLawrence Mitchell PetscMPIInt size; 248c4e6a40aSLawrence Mitchell const PetscInt *ilocal; 249c4e6a40aSLawrence Mitchell const PetscSFNode *iremote; 250c4e6a40aSLawrence Mitchell 251c4e6a40aSLawrence Mitchell PetscFunctionBegin; 2523ba16761SJacob Faibussowitsch if (!sf->graphset || !PetscDefined(USE_DEBUG)) PetscFunctionReturn(PETSC_SUCCESS); 2539566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, &ilocal, &iremote)); 2549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf), &size)); 255c4e6a40aSLawrence Mitchell for (i = 0; i < nleaves; i++) { 256c4e6a40aSLawrence Mitchell const PetscInt rank = iremote[i].rank; 257c4e6a40aSLawrence Mitchell const PetscInt remote = iremote[i].index; 258c4e6a40aSLawrence Mitchell const PetscInt leaf = ilocal ? ilocal[i] : i; 259c9cc58a2SBarry 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); 26008401ef6SPierre 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); 26108401ef6SPierre 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); 262c4e6a40aSLawrence Mitchell } 2633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 264c4e6a40aSLawrence Mitchell } 265c4e6a40aSLawrence Mitchell 26695fce210SBarry Smith /*@ 26720662ed9SBarry Smith PetscSFSetUp - set up communication structures for a `PetscSF`, after this is done it may be used to perform communication 26895fce210SBarry Smith 26995fce210SBarry Smith Collective 27095fce210SBarry Smith 2714165533cSJose E. Roman Input Parameter: 27295fce210SBarry Smith . sf - star forest communication object 27395fce210SBarry Smith 27495fce210SBarry Smith Level: beginner 27595fce210SBarry Smith 27620662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFSetFromOptions()`, `PetscSFSetType()` 27795fce210SBarry Smith @*/ 278d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetUp(PetscSF sf) 279d71ae5a4SJacob Faibussowitsch { 28095fce210SBarry Smith PetscFunctionBegin; 28129046d53SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 28229046d53SLisandro Dalcin PetscSFCheckGraphSet(sf, 1); 2833ba16761SJacob Faibussowitsch if (sf->setupcalled) PetscFunctionReturn(PETSC_SUCCESS); 2849566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_SetUp, sf, 0, 0, 0)); 2859566063dSJacob Faibussowitsch PetscCall(PetscSFCheckGraphValid_Private(sf)); 2869566063dSJacob Faibussowitsch if (!((PetscObject)sf)->type_name) PetscCall(PetscSFSetType(sf, PETSCSFBASIC)); /* Zero all sf->ops */ 287dbbe0bcdSBarry Smith PetscTryTypeMethod(sf, SetUp); 28820c24465SJunchao Zhang #if defined(PETSC_HAVE_CUDA) 28920c24465SJunchao Zhang if (sf->backend == PETSCSF_BACKEND_CUDA) { 29071438e86SJunchao Zhang sf->ops->Malloc = PetscSFMalloc_CUDA; 29171438e86SJunchao Zhang sf->ops->Free = PetscSFFree_CUDA; 29220c24465SJunchao Zhang } 29320c24465SJunchao Zhang #endif 29459af0bd3SScott Kruger #if defined(PETSC_HAVE_HIP) 29559af0bd3SScott Kruger if (sf->backend == PETSCSF_BACKEND_HIP) { 29659af0bd3SScott Kruger sf->ops->Malloc = PetscSFMalloc_HIP; 29759af0bd3SScott Kruger sf->ops->Free = PetscSFFree_HIP; 29859af0bd3SScott Kruger } 29959af0bd3SScott Kruger #endif 30020c24465SJunchao Zhang 30159af0bd3SScott Kruger # 30220c24465SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS) 30320c24465SJunchao Zhang if (sf->backend == PETSCSF_BACKEND_KOKKOS) { 30420c24465SJunchao Zhang sf->ops->Malloc = PetscSFMalloc_Kokkos; 30520c24465SJunchao Zhang sf->ops->Free = PetscSFFree_Kokkos; 30620c24465SJunchao Zhang } 30720c24465SJunchao Zhang #endif 3089566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_SetUp, sf, 0, 0, 0)); 30995fce210SBarry Smith sf->setupcalled = PETSC_TRUE; 3103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31195fce210SBarry Smith } 31295fce210SBarry Smith 3138af6ec1cSBarry Smith /*@ 314cab54364SBarry Smith PetscSFSetFromOptions - set `PetscSF` options using the options database 31595fce210SBarry Smith 31695fce210SBarry Smith Logically Collective 31795fce210SBarry Smith 3184165533cSJose E. Roman Input Parameter: 31995fce210SBarry Smith . sf - star forest 32095fce210SBarry Smith 32195fce210SBarry Smith Options Database Keys: 32220662ed9SBarry Smith + -sf_type - implementation type, see `PetscSFSetType()` 32351ccb202SJunchao Zhang . -sf_rank_order - sort composite points for gathers and scatters in rank order, gathers are non-deterministic otherwise 32420662ed9SBarry Smith . -sf_use_default_stream - Assume callers of `PetscSF` computed the input root/leafdata with the default CUDA stream. `PetscSF` will also 32520662ed9SBarry Smith use the default stream to process data. Therefore, no stream synchronization is needed between `PetscSF` and its caller (default: true). 32620662ed9SBarry Smith If true, this option only works with `-use_gpu_aware_mpi 1`. 32720662ed9SBarry Smith . -sf_use_stream_aware_mpi - Assume the underlying MPI is CUDA-stream aware and `PetscSF` won't sync streams for send/recv buffers passed to MPI (default: false). 32820662ed9SBarry Smith If true, this option only works with `-use_gpu_aware_mpi 1`. 32995fce210SBarry Smith 33020662ed9SBarry Smith - -sf_backend cuda | hip | kokkos -Select the device backend SF uses. Currently `PetscSF` has these backends: cuda, hip and Kokkos. 33159af0bd3SScott Kruger On CUDA (HIP) devices, one can choose cuda (hip) or kokkos with the default being kokkos. On other devices, 33220c24465SJunchao Zhang the only available is kokkos. 33320c24465SJunchao Zhang 33495fce210SBarry Smith Level: intermediate 335cab54364SBarry Smith 336cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`, `PetscSFSetType()` 33795fce210SBarry Smith @*/ 338d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetFromOptions(PetscSF sf) 339d71ae5a4SJacob Faibussowitsch { 34095fce210SBarry Smith PetscSFType deft; 34195fce210SBarry Smith char type[256]; 34295fce210SBarry Smith PetscBool flg; 34395fce210SBarry Smith 34495fce210SBarry Smith PetscFunctionBegin; 34595fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 346d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)sf); 34795fce210SBarry Smith deft = ((PetscObject)sf)->type_name ? ((PetscObject)sf)->type_name : PETSCSFBASIC; 3489566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-sf_type", "PetscSF implementation type", "PetscSFSetType", PetscSFList, deft, type, sizeof(type), &flg)); 3499566063dSJacob Faibussowitsch PetscCall(PetscSFSetType(sf, flg ? type : deft)); 3509566063dSJacob 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)); 3517fd2d3dbSJunchao Zhang #if defined(PETSC_HAVE_DEVICE) 35220c24465SJunchao Zhang { 35320c24465SJunchao Zhang char backendstr[32] = {0}; 35459af0bd3SScott Kruger PetscBool isCuda = PETSC_FALSE, isHip = PETSC_FALSE, isKokkos = PETSC_FALSE, set; 35520c24465SJunchao Zhang /* Change the defaults set in PetscSFCreate() with command line options */ 356d5b43468SJose E. Roman PetscCall(PetscOptionsBool("-sf_unknown_input_stream", "SF root/leafdata is computed on arbitrary streams unknown to SF", "PetscSFSetFromOptions", sf->unknown_input_stream, &sf->unknown_input_stream, NULL)); 3579566063dSJacob 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)); 3589566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-sf_backend", "Select the device backend SF uses", "PetscSFSetFromOptions", NULL, backendstr, sizeof(backendstr), &set)); 3599566063dSJacob Faibussowitsch PetscCall(PetscStrcasecmp("cuda", backendstr, &isCuda)); 3609566063dSJacob Faibussowitsch PetscCall(PetscStrcasecmp("kokkos", backendstr, &isKokkos)); 3619566063dSJacob Faibussowitsch PetscCall(PetscStrcasecmp("hip", backendstr, &isHip)); 36259af0bd3SScott Kruger #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) 36320c24465SJunchao Zhang if (isCuda) sf->backend = PETSCSF_BACKEND_CUDA; 36420c24465SJunchao Zhang else if (isKokkos) sf->backend = PETSCSF_BACKEND_KOKKOS; 36559af0bd3SScott Kruger else if (isHip) sf->backend = PETSCSF_BACKEND_HIP; 36628b400f6SJacob 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); 36720c24465SJunchao Zhang #elif defined(PETSC_HAVE_KOKKOS) 36808401ef6SPierre Jolivet PetscCheck(!set || isKokkos, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "-sf_backend %s is not supported. You can only choose kokkos", backendstr); 36920c24465SJunchao Zhang #endif 37020c24465SJunchao Zhang } 371c2a741eeSJunchao Zhang #endif 372dbbe0bcdSBarry Smith PetscTryTypeMethod(sf, SetFromOptions, PetscOptionsObject); 373d0609cedSBarry Smith PetscOptionsEnd(); 3743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 37595fce210SBarry Smith } 37695fce210SBarry Smith 37729046d53SLisandro Dalcin /*@ 37895fce210SBarry Smith PetscSFSetRankOrder - sort multi-points for gathers and scatters by rank order 37995fce210SBarry Smith 38095fce210SBarry Smith Logically Collective 38195fce210SBarry Smith 3824165533cSJose E. Roman Input Parameters: 38395fce210SBarry Smith + sf - star forest 384cab54364SBarry Smith - flg - `PETSC_TRUE` to sort, `PETSC_FALSE` to skip sorting (lower setup cost, but non-deterministic) 38595fce210SBarry Smith 38695fce210SBarry Smith Level: advanced 38795fce210SBarry Smith 38820662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFGatherBegin()`, `PetscSFScatterBegin()` 38995fce210SBarry Smith @*/ 390d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetRankOrder(PetscSF sf, PetscBool flg) 391d71ae5a4SJacob Faibussowitsch { 39295fce210SBarry Smith PetscFunctionBegin; 39395fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 39495fce210SBarry Smith PetscValidLogicalCollectiveBool(sf, flg, 2); 39528b400f6SJacob Faibussowitsch PetscCheck(!sf->multi, PetscObjectComm((PetscObject)sf), PETSC_ERR_ARG_WRONGSTATE, "Rank ordering must be set before first call to PetscSFGatherBegin() or PetscSFScatterBegin()"); 39695fce210SBarry Smith sf->rankorder = flg; 3973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39895fce210SBarry Smith } 39995fce210SBarry Smith 4008dbb0df6SBarry Smith /*@C 40195fce210SBarry Smith PetscSFSetGraph - Set a parallel star forest 40295fce210SBarry Smith 40395fce210SBarry Smith Collective 40495fce210SBarry Smith 4054165533cSJose E. Roman Input Parameters: 40695fce210SBarry Smith + sf - star forest 40795fce210SBarry Smith . nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) 40895fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process 40920662ed9SBarry Smith . ilocal - locations of leaves in leafdata buffers, pass `NULL` for contiguous storage (locations must be >= 0, enforced 410c4e6a40aSLawrence Mitchell during setup in debug mode) 41120662ed9SBarry Smith . localmode - copy mode for `ilocal` 412c4e6a40aSLawrence Mitchell . iremote - remote locations of root vertices for each leaf on the current process (locations must be >= 0, enforced 413c4e6a40aSLawrence Mitchell during setup in debug mode) 41420662ed9SBarry Smith - remotemode - copy mode for `iremote` 41595fce210SBarry Smith 41695fce210SBarry Smith Level: intermediate 41795fce210SBarry Smith 41895452b02SPatrick Sanan Notes: 41920662ed9SBarry Smith Leaf indices in `ilocal` must be unique, otherwise an error occurs. 42038ab3f8aSBarry Smith 42120662ed9SBarry Smith Input arrays `ilocal` and `iremote` follow the `PetscCopyMode` semantics. 42220662ed9SBarry Smith In particular, if `localmode` or `remotemode` is `PETSC_OWN_POINTER` or `PETSC_USE_POINTER`, 423db2b9530SVaclav Hapla PETSc might modify the respective array; 42420662ed9SBarry Smith if `PETSC_USE_POINTER`, the user must delete the array after `PetscSFDestroy()`. 425cab54364SBarry Smith Only if `PETSC_COPY_VALUES` is used, the respective array is guaranteed to stay intact and a const array can be passed (but a cast to non-const is needed). 426db2b9530SVaclav Hapla 427cab54364SBarry Smith Fortran Note: 42820662ed9SBarry Smith In Fortran you must use `PETSC_COPY_VALUES` for `localmode` and `remotemode`. 429c4e6a40aSLawrence Mitchell 430cab54364SBarry Smith Developer Note: 431db2b9530SVaclav Hapla We sort leaves to check for duplicates and contiguousness and to find minleaf/maxleaf. 43220662ed9SBarry Smith This also allows to compare leaf sets of two `PetscSF`s easily. 43372bf8598SVaclav Hapla 43420662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFGetGraph()` 43595fce210SBarry Smith @*/ 436d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraph(PetscSF sf, PetscInt nroots, PetscInt nleaves, PetscInt *ilocal, PetscCopyMode localmode, PetscSFNode *iremote, PetscCopyMode remotemode) 437d71ae5a4SJacob Faibussowitsch { 438db2b9530SVaclav Hapla PetscBool unique, contiguous; 43995fce210SBarry Smith 44095fce210SBarry Smith PetscFunctionBegin; 44195fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 44229046d53SLisandro Dalcin if (nleaves > 0 && ilocal) PetscValidIntPointer(ilocal, 4); 44329046d53SLisandro Dalcin if (nleaves > 0) PetscValidPointer(iremote, 6); 44408401ef6SPierre Jolivet PetscCheck(nroots >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nroots %" PetscInt_FMT ", cannot be negative", nroots); 44508401ef6SPierre Jolivet PetscCheck(nleaves >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nleaves %" PetscInt_FMT ", cannot be negative", nleaves); 4468da24d32SBarry Smith /* enums may be handled as unsigned by some compilers, NVHPC for example, the int cast 4478da24d32SBarry Smith * below is to prevent NVHPC from warning about meaningless comparison of unsigned with zero */ 4488da24d32SBarry Smith PetscCheck((int)localmode >= PETSC_COPY_VALUES && localmode <= PETSC_USE_POINTER, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wrong localmode %d", localmode); 4498da24d32SBarry Smith PetscCheck((int)remotemode >= PETSC_COPY_VALUES && remotemode <= PETSC_USE_POINTER, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wrong remotemode %d", remotemode); 45029046d53SLisandro Dalcin 4512a67d2daSStefano Zampini if (sf->nroots >= 0) { /* Reset only if graph already set */ 4529566063dSJacob Faibussowitsch PetscCall(PetscSFReset(sf)); 4532a67d2daSStefano Zampini } 4542a67d2daSStefano Zampini 4559566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_SetGraph, sf, 0, 0, 0)); 45629046d53SLisandro Dalcin 45795fce210SBarry Smith sf->nroots = nroots; 45895fce210SBarry Smith sf->nleaves = nleaves; 45929046d53SLisandro Dalcin 460db2b9530SVaclav Hapla if (localmode == PETSC_COPY_VALUES && ilocal) { 461db2b9530SVaclav Hapla PetscInt *tlocal = NULL; 462db2b9530SVaclav Hapla 4639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &tlocal)); 4649566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tlocal, ilocal, nleaves)); 465db2b9530SVaclav Hapla ilocal = tlocal; 466db2b9530SVaclav Hapla } 467db2b9530SVaclav Hapla if (remotemode == PETSC_COPY_VALUES) { 468db2b9530SVaclav Hapla PetscSFNode *tremote = NULL; 469db2b9530SVaclav Hapla 4709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &tremote)); 4719566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tremote, iremote, nleaves)); 472db2b9530SVaclav Hapla iremote = tremote; 473db2b9530SVaclav Hapla } 474db2b9530SVaclav Hapla 47529046d53SLisandro Dalcin if (nleaves && ilocal) { 476db2b9530SVaclav Hapla PetscSFNode work; 477db2b9530SVaclav Hapla 4789566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithDataArray(nleaves, ilocal, iremote, sizeof(PetscSFNode), &work)); 4799566063dSJacob Faibussowitsch PetscCall(PetscSortedCheckDupsInt(nleaves, ilocal, &unique)); 480db2b9530SVaclav Hapla unique = PetscNot(unique); 481db2b9530SVaclav 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"); 482db2b9530SVaclav Hapla sf->minleaf = ilocal[0]; 483db2b9530SVaclav Hapla sf->maxleaf = ilocal[nleaves - 1]; 484db2b9530SVaclav Hapla contiguous = (PetscBool)(unique && ilocal[0] == 0 && ilocal[nleaves - 1] == nleaves - 1); 48529046d53SLisandro Dalcin } else { 48629046d53SLisandro Dalcin sf->minleaf = 0; 48729046d53SLisandro Dalcin sf->maxleaf = nleaves - 1; 488db2b9530SVaclav Hapla unique = PETSC_TRUE; 489db2b9530SVaclav Hapla contiguous = PETSC_TRUE; 49029046d53SLisandro Dalcin } 49129046d53SLisandro Dalcin 492db2b9530SVaclav Hapla if (contiguous) { 493db2b9530SVaclav Hapla if (localmode == PETSC_USE_POINTER) { 494db2b9530SVaclav Hapla ilocal = NULL; 495db2b9530SVaclav Hapla } else { 4969566063dSJacob Faibussowitsch PetscCall(PetscFree(ilocal)); 497db2b9530SVaclav Hapla } 498db2b9530SVaclav Hapla } 499db2b9530SVaclav Hapla sf->mine = ilocal; 500db2b9530SVaclav Hapla if (localmode == PETSC_USE_POINTER) { 50129046d53SLisandro Dalcin sf->mine_alloc = NULL; 502db2b9530SVaclav Hapla } else { 503db2b9530SVaclav Hapla sf->mine_alloc = ilocal; 50495fce210SBarry Smith } 505db2b9530SVaclav Hapla sf->remote = iremote; 506db2b9530SVaclav Hapla if (remotemode == PETSC_USE_POINTER) { 50729046d53SLisandro Dalcin sf->remote_alloc = NULL; 508db2b9530SVaclav Hapla } else { 509db2b9530SVaclav Hapla sf->remote_alloc = iremote; 51095fce210SBarry Smith } 5119566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_SetGraph, sf, 0, 0, 0)); 51229046d53SLisandro Dalcin sf->graphset = PETSC_TRUE; 5133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51495fce210SBarry Smith } 51595fce210SBarry Smith 51629046d53SLisandro Dalcin /*@ 517cab54364SBarry Smith PetscSFSetGraphWithPattern - Sets the graph of a `PetscSF` with a specific pattern 518dd5b3ca6SJunchao Zhang 519dd5b3ca6SJunchao Zhang Collective 520dd5b3ca6SJunchao Zhang 521dd5b3ca6SJunchao Zhang Input Parameters: 522cab54364SBarry Smith + sf - The `PetscSF` 523cab54364SBarry Smith . map - Layout of roots over all processes (insignificant when pattern is `PETSCSF_PATTERN_ALLTOALL`) 524cab54364SBarry Smith - pattern - One of `PETSCSF_PATTERN_ALLGATHER`, `PETSCSF_PATTERN_GATHER`, `PETSCSF_PATTERN_ALLTOALL` 525cab54364SBarry Smith 526cab54364SBarry Smith Level: intermediate 527dd5b3ca6SJunchao Zhang 528dd5b3ca6SJunchao Zhang Notes: 52920662ed9SBarry Smith It is easier to explain `PetscSFPattern` using vectors. Suppose we have an MPI vector `x` and its `PetscLayout` is `map`. 53020662ed9SBarry Smith `n` and `N` are the local and global sizes of `x` respectively. 531dd5b3ca6SJunchao Zhang 53220662ed9SBarry Smith With `PETSCSF_PATTERN_ALLGATHER`, the routine creates a graph that if one does `PetscSFBcastBegin()`/`PetscSFBcastEnd()` on it, it will copy `x` to 53320662ed9SBarry Smith sequential vectors `y` on all MPI processes. 534dd5b3ca6SJunchao Zhang 53520662ed9SBarry Smith With `PETSCSF_PATTERN_GATHER`, the routine creates a graph that if one does `PetscSFBcastBegin()`/`PetscSFBcastEnd()` on it, it will copy `x` to a 53620662ed9SBarry Smith sequential vector `y` on rank 0. 537dd5b3ca6SJunchao Zhang 53820662ed9SBarry Smith In above cases, entries of `x` are roots and entries of `y` are leaves. 539dd5b3ca6SJunchao Zhang 54020662ed9SBarry Smith With `PETSCSF_PATTERN_ALLTOALL`, map is insignificant. Suppose NP is size of `sf`'s communicator. The routine 541dd5b3ca6SJunchao Zhang creates a graph that every rank has NP leaves and NP roots. On rank i, its leaf j is connected to root i 542cab54364SBarry Smith of rank j. Here 0 <=i,j<NP. It is a kind of `MPI_Alltoall()` with sendcount/recvcount being 1. Note that it does 543dd5b3ca6SJunchao Zhang not mean one can not send multiple items. One just needs to create a new MPI datatype for the mulptiple data 544cab54364SBarry Smith items with `MPI_Type_contiguous` and use that as the <unit> argument in SF routines. 545dd5b3ca6SJunchao Zhang 546dd5b3ca6SJunchao Zhang In this case, roots and leaves are symmetric. 547dd5b3ca6SJunchao Zhang 548cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFGetGraph()` 549dd5b3ca6SJunchao Zhang @*/ 550d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetGraphWithPattern(PetscSF sf, PetscLayout map, PetscSFPattern pattern) 551d71ae5a4SJacob Faibussowitsch { 552dd5b3ca6SJunchao Zhang MPI_Comm comm; 553dd5b3ca6SJunchao Zhang PetscInt n, N, res[2]; 554dd5b3ca6SJunchao Zhang PetscMPIInt rank, size; 555dd5b3ca6SJunchao Zhang PetscSFType type; 556dd5b3ca6SJunchao Zhang 557dd5b3ca6SJunchao Zhang PetscFunctionBegin; 5582abc8c78SJacob Faibussowitsch PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 5592abc8c78SJacob Faibussowitsch if (pattern != PETSCSF_PATTERN_ALLTOALL) PetscValidPointer(map, 2); 5609566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 5612c71b3e2SJacob Faibussowitsch PetscCheck(pattern >= PETSCSF_PATTERN_ALLGATHER && pattern <= PETSCSF_PATTERN_ALLTOALL, comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PetscSFPattern %d", pattern); 5629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 564dd5b3ca6SJunchao Zhang 565dd5b3ca6SJunchao Zhang if (pattern == PETSCSF_PATTERN_ALLTOALL) { 566dd5b3ca6SJunchao Zhang type = PETSCSFALLTOALL; 5679566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &sf->map)); 5689566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(sf->map, size)); 5699566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetSize(sf->map, ((PetscInt)size) * size)); 5709566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(sf->map)); 571dd5b3ca6SJunchao Zhang } else { 5729566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(map, &n)); 5739566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(map, &N)); 574dd5b3ca6SJunchao Zhang res[0] = n; 575dd5b3ca6SJunchao Zhang res[1] = -n; 576dd5b3ca6SJunchao Zhang /* Check if n are same over all ranks so that we can optimize it */ 5771c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE, res, 2, MPIU_INT, MPI_MAX, comm)); 578dd5b3ca6SJunchao Zhang if (res[0] == -res[1]) { /* same n */ 579dd5b3ca6SJunchao Zhang type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHER : PETSCSFGATHER; 580dd5b3ca6SJunchao Zhang } else { 581dd5b3ca6SJunchao Zhang type = (pattern == PETSCSF_PATTERN_ALLGATHER) ? PETSCSFALLGATHERV : PETSCSFGATHERV; 582dd5b3ca6SJunchao Zhang } 5839566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(map, &sf->map)); 584dd5b3ca6SJunchao Zhang } 5859566063dSJacob Faibussowitsch PetscCall(PetscSFSetType(sf, type)); 586dd5b3ca6SJunchao Zhang 587dd5b3ca6SJunchao Zhang sf->pattern = pattern; 588dd5b3ca6SJunchao Zhang sf->mine = NULL; /* Contiguous */ 589dd5b3ca6SJunchao Zhang 590dd5b3ca6SJunchao Zhang /* Set nleaves, nroots here in case user calls PetscSFGetGraph, which is legal to call even before PetscSFSetUp is called. 591dd5b3ca6SJunchao Zhang Also set other easy stuff. 592dd5b3ca6SJunchao Zhang */ 593dd5b3ca6SJunchao Zhang if (pattern == PETSCSF_PATTERN_ALLGATHER) { 594dd5b3ca6SJunchao Zhang sf->nleaves = N; 595dd5b3ca6SJunchao Zhang sf->nroots = n; 596dd5b3ca6SJunchao Zhang sf->nranks = size; 597dd5b3ca6SJunchao Zhang sf->minleaf = 0; 598dd5b3ca6SJunchao Zhang sf->maxleaf = N - 1; 599dd5b3ca6SJunchao Zhang } else if (pattern == PETSCSF_PATTERN_GATHER) { 600dd5b3ca6SJunchao Zhang sf->nleaves = rank ? 0 : N; 601dd5b3ca6SJunchao Zhang sf->nroots = n; 602dd5b3ca6SJunchao Zhang sf->nranks = rank ? 0 : size; 603dd5b3ca6SJunchao Zhang sf->minleaf = 0; 604dd5b3ca6SJunchao Zhang sf->maxleaf = rank ? -1 : N - 1; 605dd5b3ca6SJunchao Zhang } else if (pattern == PETSCSF_PATTERN_ALLTOALL) { 606dd5b3ca6SJunchao Zhang sf->nleaves = size; 607dd5b3ca6SJunchao Zhang sf->nroots = size; 608dd5b3ca6SJunchao Zhang sf->nranks = size; 609dd5b3ca6SJunchao Zhang sf->minleaf = 0; 610dd5b3ca6SJunchao Zhang sf->maxleaf = size - 1; 611dd5b3ca6SJunchao Zhang } 612dd5b3ca6SJunchao Zhang sf->ndranks = 0; /* We do not need to separate out distinguished ranks for patterned graphs to improve communication performance */ 613dd5b3ca6SJunchao Zhang sf->graphset = PETSC_TRUE; 6143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 615dd5b3ca6SJunchao Zhang } 616dd5b3ca6SJunchao Zhang 617dd5b3ca6SJunchao Zhang /*@ 618cab54364SBarry Smith PetscSFCreateInverseSF - given a `PetscSF` in which all vertices have degree 1, creates the inverse map 61995fce210SBarry Smith 62095fce210SBarry Smith Collective 62195fce210SBarry Smith 6224165533cSJose E. Roman Input Parameter: 62395fce210SBarry Smith . sf - star forest to invert 62495fce210SBarry Smith 6254165533cSJose E. Roman Output Parameter: 62620662ed9SBarry Smith . isf - inverse of `sf` 6274165533cSJose E. Roman 62895fce210SBarry Smith Level: advanced 62995fce210SBarry Smith 63095fce210SBarry Smith Notes: 63195fce210SBarry Smith All roots must have degree 1. 63295fce210SBarry Smith 63395fce210SBarry Smith The local space may be a permutation, but cannot be sparse. 63495fce210SBarry Smith 63520662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFSetGraph()` 63695fce210SBarry Smith @*/ 637d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateInverseSF(PetscSF sf, PetscSF *isf) 638d71ae5a4SJacob Faibussowitsch { 63995fce210SBarry Smith PetscMPIInt rank; 64095fce210SBarry Smith PetscInt i, nroots, nleaves, maxlocal, count, *newilocal; 64195fce210SBarry Smith const PetscInt *ilocal; 64295fce210SBarry Smith PetscSFNode *roots, *leaves; 64395fce210SBarry Smith 64495fce210SBarry Smith PetscFunctionBegin; 64529046d53SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 64629046d53SLisandro Dalcin PetscSFCheckGraphSet(sf, 1); 64729046d53SLisandro Dalcin PetscValidPointer(isf, 2); 64829046d53SLisandro Dalcin 6499566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, NULL)); 65029046d53SLisandro Dalcin maxlocal = sf->maxleaf + 1; /* TODO: We should use PetscSFGetLeafRange() */ 65129046d53SLisandro Dalcin 6529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank)); 6539566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nroots, &roots, maxlocal, &leaves)); 654ae9aee6dSMatthew G. Knepley for (i = 0; i < maxlocal; i++) { 65595fce210SBarry Smith leaves[i].rank = rank; 65695fce210SBarry Smith leaves[i].index = i; 65795fce210SBarry Smith } 65895fce210SBarry Smith for (i = 0; i < nroots; i++) { 65995fce210SBarry Smith roots[i].rank = -1; 66095fce210SBarry Smith roots[i].index = -1; 66195fce210SBarry Smith } 6629566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_2INT, leaves, roots, MPI_REPLACE)); 6639566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_2INT, leaves, roots, MPI_REPLACE)); 66495fce210SBarry Smith 66595fce210SBarry Smith /* Check whether our leaves are sparse */ 6669371c9d4SSatish Balay for (i = 0, count = 0; i < nroots; i++) 6679371c9d4SSatish Balay if (roots[i].rank >= 0) count++; 66895fce210SBarry Smith if (count == nroots) newilocal = NULL; 6699371c9d4SSatish Balay else { /* Index for sparse leaves and compact "roots" array (which is to become our leaves). */ PetscCall(PetscMalloc1(count, &newilocal)); 67095fce210SBarry Smith for (i = 0, count = 0; i < nroots; i++) { 67195fce210SBarry Smith if (roots[i].rank >= 0) { 67295fce210SBarry Smith newilocal[count] = i; 67395fce210SBarry Smith roots[count].rank = roots[i].rank; 67495fce210SBarry Smith roots[count].index = roots[i].index; 67595fce210SBarry Smith count++; 67695fce210SBarry Smith } 67795fce210SBarry Smith } 67895fce210SBarry Smith } 67995fce210SBarry Smith 6809566063dSJacob Faibussowitsch PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_CONFONLY, isf)); 6819566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*isf, maxlocal, count, newilocal, PETSC_OWN_POINTER, roots, PETSC_COPY_VALUES)); 6829566063dSJacob Faibussowitsch PetscCall(PetscFree2(roots, leaves)); 6833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68495fce210SBarry Smith } 68595fce210SBarry Smith 68695fce210SBarry Smith /*@ 687cab54364SBarry Smith PetscSFDuplicate - duplicate a `PetscSF`, optionally preserving rank connectivity and graph 68895fce210SBarry Smith 68995fce210SBarry Smith Collective 69095fce210SBarry Smith 6914165533cSJose E. Roman Input Parameters: 69295fce210SBarry Smith + sf - communication object to duplicate 693cab54364SBarry Smith - opt - `PETSCSF_DUPLICATE_CONFONLY`, `PETSCSF_DUPLICATE_RANKS`, or `PETSCSF_DUPLICATE_GRAPH` (see `PetscSFDuplicateOption`) 69495fce210SBarry Smith 6954165533cSJose E. Roman Output Parameter: 69695fce210SBarry Smith . newsf - new communication object 69795fce210SBarry Smith 69895fce210SBarry Smith Level: beginner 69995fce210SBarry Smith 70020662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFCreate()`, `PetscSFSetType()`, `PetscSFSetGraph()` 70195fce210SBarry Smith @*/ 702d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFDuplicate(PetscSF sf, PetscSFDuplicateOption opt, PetscSF *newsf) 703d71ae5a4SJacob Faibussowitsch { 70429046d53SLisandro Dalcin PetscSFType type; 70597929ea7SJunchao Zhang MPI_Datatype dtype = MPIU_SCALAR; 70695fce210SBarry Smith 70795fce210SBarry Smith PetscFunctionBegin; 70829046d53SLisandro Dalcin PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 70929046d53SLisandro Dalcin PetscValidLogicalCollectiveEnum(sf, opt, 2); 71029046d53SLisandro Dalcin PetscValidPointer(newsf, 3); 7119566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sf), newsf)); 7129566063dSJacob Faibussowitsch PetscCall(PetscSFGetType(sf, &type)); 7139566063dSJacob Faibussowitsch if (type) PetscCall(PetscSFSetType(*newsf, type)); 71435cb6cd3SPierre Jolivet (*newsf)->allow_multi_leaves = sf->allow_multi_leaves; /* Dup this flag earlier since PetscSFSetGraph() below checks on this flag */ 71595fce210SBarry Smith if (opt == PETSCSF_DUPLICATE_GRAPH) { 716dd5b3ca6SJunchao Zhang PetscSFCheckGraphSet(sf, 1); 717dd5b3ca6SJunchao Zhang if (sf->pattern == PETSCSF_PATTERN_GENERAL) { 71895fce210SBarry Smith PetscInt nroots, nleaves; 71995fce210SBarry Smith const PetscInt *ilocal; 72095fce210SBarry Smith const PetscSFNode *iremote; 7219566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, &iremote)); 7229566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*newsf, nroots, nleaves, (PetscInt *)ilocal, PETSC_COPY_VALUES, (PetscSFNode *)iremote, PETSC_COPY_VALUES)); 723dd5b3ca6SJunchao Zhang } else { 7249566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphWithPattern(*newsf, sf->map, sf->pattern)); 725dd5b3ca6SJunchao Zhang } 72695fce210SBarry Smith } 72797929ea7SJunchao Zhang /* Since oldtype is committed, so is newtype, according to MPI */ 7289566063dSJacob Faibussowitsch if (sf->vscat.bs > 1) PetscCallMPI(MPI_Type_dup(sf->vscat.unit, &dtype)); 72997929ea7SJunchao Zhang (*newsf)->vscat.bs = sf->vscat.bs; 73097929ea7SJunchao Zhang (*newsf)->vscat.unit = dtype; 73197929ea7SJunchao Zhang (*newsf)->vscat.to_n = sf->vscat.to_n; 73297929ea7SJunchao Zhang (*newsf)->vscat.from_n = sf->vscat.from_n; 73397929ea7SJunchao Zhang /* Do not copy lsf. Build it on demand since it is rarely used */ 73497929ea7SJunchao Zhang 73520c24465SJunchao Zhang #if defined(PETSC_HAVE_DEVICE) 73620c24465SJunchao Zhang (*newsf)->backend = sf->backend; 73771438e86SJunchao Zhang (*newsf)->unknown_input_stream = sf->unknown_input_stream; 73820c24465SJunchao Zhang (*newsf)->use_gpu_aware_mpi = sf->use_gpu_aware_mpi; 73920c24465SJunchao Zhang (*newsf)->use_stream_aware_mpi = sf->use_stream_aware_mpi; 74020c24465SJunchao Zhang #endif 741dbbe0bcdSBarry Smith PetscTryTypeMethod(sf, Duplicate, opt, *newsf); 74220c24465SJunchao Zhang /* Don't do PetscSFSetUp() since the new sf's graph might have not been set. */ 7433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 74495fce210SBarry Smith } 74595fce210SBarry Smith 74695fce210SBarry Smith /*@C 74795fce210SBarry Smith PetscSFGetGraph - Get the graph specifying a parallel star forest 74895fce210SBarry Smith 74995fce210SBarry Smith Not Collective 75095fce210SBarry Smith 7514165533cSJose E. Roman Input Parameter: 75295fce210SBarry Smith . sf - star forest 75395fce210SBarry Smith 7544165533cSJose E. Roman Output Parameters: 75595fce210SBarry Smith + nroots - number of root vertices on the current process (these are possible targets for other process to attach leaves) 75695fce210SBarry Smith . nleaves - number of leaf vertices on the current process, each of these references a root on any process 75720662ed9SBarry Smith . ilocal - locations of leaves in leafdata buffers (if returned value is `NULL`, it means leaves are in contiguous storage) 75895fce210SBarry Smith - iremote - remote locations of root vertices for each leaf on the current process 75995fce210SBarry Smith 760cab54364SBarry Smith Level: intermediate 761cab54364SBarry Smith 762373e0d91SLisandro Dalcin Notes: 76320662ed9SBarry Smith We are not currently requiring that the graph is set, thus returning `nroots` = -1 if it has not been set yet 764373e0d91SLisandro Dalcin 76520662ed9SBarry Smith The returned `ilocal` and `iremote` might contain values in different order than the input ones in `PetscSFSetGraph()` 766db2b9530SVaclav Hapla 7678dbb0df6SBarry Smith Fortran Notes: 76820662ed9SBarry Smith The returned `iremote` array is a copy and must be deallocated after use. Consequently, if you 76920662ed9SBarry Smith want to update the graph, you must call `PetscSFSetGraph()` after modifying the `iremote` array. 7708dbb0df6SBarry Smith 77120662ed9SBarry Smith To check for a `NULL` `ilocal` use 7728dbb0df6SBarry Smith $ if (loc(ilocal) == loc(PETSC_NULL_INTEGER)) then 773ca797d7aSLawrence Mitchell 77420662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()` 77595fce210SBarry Smith @*/ 776d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGraph(PetscSF sf, PetscInt *nroots, PetscInt *nleaves, const PetscInt **ilocal, const PetscSFNode **iremote) 777d71ae5a4SJacob Faibussowitsch { 77895fce210SBarry Smith PetscFunctionBegin; 77995fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 780b8dee149SJunchao Zhang if (sf->ops->GetGraph) { 7819566063dSJacob Faibussowitsch PetscCall((sf->ops->GetGraph)(sf, nroots, nleaves, ilocal, iremote)); 782b8dee149SJunchao Zhang } else { 78395fce210SBarry Smith if (nroots) *nroots = sf->nroots; 78495fce210SBarry Smith if (nleaves) *nleaves = sf->nleaves; 78595fce210SBarry Smith if (ilocal) *ilocal = sf->mine; 78695fce210SBarry Smith if (iremote) *iremote = sf->remote; 787b8dee149SJunchao Zhang } 7883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 78995fce210SBarry Smith } 79095fce210SBarry Smith 79129046d53SLisandro Dalcin /*@ 79295fce210SBarry Smith PetscSFGetLeafRange - Get the active leaf ranges 79395fce210SBarry Smith 79495fce210SBarry Smith Not Collective 79595fce210SBarry Smith 7964165533cSJose E. Roman Input Parameter: 79795fce210SBarry Smith . sf - star forest 79895fce210SBarry Smith 7994165533cSJose E. Roman Output Parameters: 80020662ed9SBarry Smith + minleaf - minimum active leaf on this process. Returns 0 if there are no leaves. 80120662ed9SBarry Smith - maxleaf - maximum active leaf on this process. Returns -1 if there are no leaves. 80295fce210SBarry Smith 80395fce210SBarry Smith Level: developer 80495fce210SBarry Smith 80520662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFType`, `PetscSFCreate()`, `PetscSFView()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 80695fce210SBarry Smith @*/ 807d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetLeafRange(PetscSF sf, PetscInt *minleaf, PetscInt *maxleaf) 808d71ae5a4SJacob Faibussowitsch { 80995fce210SBarry Smith PetscFunctionBegin; 81095fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 81129046d53SLisandro Dalcin PetscSFCheckGraphSet(sf, 1); 81295fce210SBarry Smith if (minleaf) *minleaf = sf->minleaf; 81395fce210SBarry Smith if (maxleaf) *maxleaf = sf->maxleaf; 8143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 81595fce210SBarry Smith } 81695fce210SBarry Smith 81795fce210SBarry Smith /*@C 818cab54364SBarry Smith PetscSFViewFromOptions - View a `PetscSF` based on arguments in the options database 819fe2efc57SMark 820*20f4b53cSBarry Smith Collective 821fe2efc57SMark 822fe2efc57SMark Input Parameters: 823fe2efc57SMark + A - the star forest 824cab54364SBarry Smith . obj - Optional object that provides the prefix for the option names 825736c3998SJose E. Roman - name - command line option 826fe2efc57SMark 827fe2efc57SMark Level: intermediate 828cab54364SBarry Smith 82920662ed9SBarry Smith Note: 83020662ed9SBarry Smith See `PetscObjectViewFromOptions()` for possible `PetscViewer` and `PetscViewerFormat` 83120662ed9SBarry Smith 832db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFView`, `PetscObjectViewFromOptions()`, `PetscSFCreate()` 833fe2efc57SMark @*/ 834d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFViewFromOptions(PetscSF A, PetscObject obj, const char name[]) 835d71ae5a4SJacob Faibussowitsch { 836fe2efc57SMark PetscFunctionBegin; 837fe2efc57SMark PetscValidHeaderSpecific(A, PETSCSF_CLASSID, 1); 8389566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 8393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 840fe2efc57SMark } 841fe2efc57SMark 842fe2efc57SMark /*@C 84395fce210SBarry Smith PetscSFView - view a star forest 84495fce210SBarry Smith 84595fce210SBarry Smith Collective 84695fce210SBarry Smith 8474165533cSJose E. Roman Input Parameters: 84895fce210SBarry Smith + sf - star forest 849cab54364SBarry Smith - viewer - viewer to display graph, for example `PETSC_VIEWER_STDOUT_WORLD` 85095fce210SBarry Smith 85195fce210SBarry Smith Level: beginner 85295fce210SBarry Smith 853cab54364SBarry Smith .seealso: `PetscSF`, `PetscViewer`, `PetscSFCreate()`, `PetscSFSetGraph()` 85495fce210SBarry Smith @*/ 855d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFView(PetscSF sf, PetscViewer viewer) 856d71ae5a4SJacob Faibussowitsch { 85795fce210SBarry Smith PetscBool iascii; 85895fce210SBarry Smith PetscViewerFormat format; 85995fce210SBarry Smith 86095fce210SBarry Smith PetscFunctionBegin; 86195fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 8629566063dSJacob Faibussowitsch if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sf), &viewer)); 86395fce210SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 86495fce210SBarry Smith PetscCheckSameComm(sf, 1, viewer, 2); 8659566063dSJacob Faibussowitsch if (sf->graphset) PetscCall(PetscSFSetUp(sf)); 8669566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 86753dd6d7dSJunchao Zhang if (iascii && viewer->format != PETSC_VIEWER_ASCII_MATLAB) { 86895fce210SBarry Smith PetscMPIInt rank; 86981bfa7aaSJed Brown PetscInt ii, i, j; 87095fce210SBarry Smith 8719566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)sf, viewer)); 8729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 873dd5b3ca6SJunchao Zhang if (sf->pattern == PETSCSF_PATTERN_GENERAL) { 87480153354SVaclav Hapla if (!sf->graphset) { 8759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "PetscSFSetGraph() has not been called yet\n")); 8769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 8773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 87880153354SVaclav Hapla } 8799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank)); 8809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 8819566063dSJacob 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)); 88248a46eb9SPierre Jolivet for (i = 0; i < sf->nleaves; i++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " <- (%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, sf->mine ? sf->mine[i] : i, sf->remote[i].rank, sf->remote[i].index)); 8839566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 8849566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 88595fce210SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 88681bfa7aaSJed Brown PetscMPIInt *tmpranks, *perm; 8879566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(sf->nranks, &tmpranks, sf->nranks, &perm)); 8889566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmpranks, sf->ranks, sf->nranks)); 88981bfa7aaSJed Brown for (i = 0; i < sf->nranks; i++) perm[i] = i; 8909566063dSJacob Faibussowitsch PetscCall(PetscSortMPIIntWithArray(sf->nranks, tmpranks, perm)); 8919566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Roots referenced by my leaves, by rank\n", rank)); 89281bfa7aaSJed Brown for (ii = 0; ii < sf->nranks; ii++) { 89381bfa7aaSJed Brown i = perm[ii]; 8949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %d: %" PetscInt_FMT " edges\n", rank, sf->ranks[i], sf->roffset[i + 1] - sf->roffset[i])); 89548a46eb9SPierre Jolivet for (j = sf->roffset[i]; j < sf->roffset[i + 1]; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " <- %" PetscInt_FMT "\n", rank, sf->rmine[j], sf->rremote[j])); 89695fce210SBarry Smith } 8979566063dSJacob Faibussowitsch PetscCall(PetscFree2(tmpranks, perm)); 89895fce210SBarry Smith } 8999566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 9009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 901dd5b3ca6SJunchao Zhang } 9029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 90395fce210SBarry Smith } 904dbbe0bcdSBarry Smith PetscTryTypeMethod(sf, View, viewer); 9053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 90695fce210SBarry Smith } 90795fce210SBarry Smith 90895fce210SBarry Smith /*@C 909dec1416fSJunchao Zhang PetscSFGetRootRanks - Get root ranks and number of vertices referenced by leaves on this process 91095fce210SBarry Smith 91195fce210SBarry Smith Not Collective 91295fce210SBarry Smith 9134165533cSJose E. Roman Input Parameter: 91495fce210SBarry Smith . sf - star forest 91595fce210SBarry Smith 9164165533cSJose E. Roman Output Parameters: 91795fce210SBarry Smith + nranks - number of ranks referenced by local part 91820662ed9SBarry Smith . ranks - [`nranks`] array of ranks 91920662ed9SBarry Smith . roffset - [`nranks`+1] offset in `rmine`/`rremote` for each rank 92020662ed9SBarry Smith . rmine - [`roffset`[`nranks`]] concatenated array holding local indices referencing each remote rank 92120662ed9SBarry Smith - rremote - [`roffset`[`nranks`]] concatenated array holding remote indices referenced for each remote rank 92295fce210SBarry Smith 92395fce210SBarry Smith Level: developer 92495fce210SBarry Smith 925cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGetLeafRanks()` 92695fce210SBarry Smith @*/ 927d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetRootRanks(PetscSF sf, PetscInt *nranks, const PetscMPIInt **ranks, const PetscInt **roffset, const PetscInt **rmine, const PetscInt **rremote) 928d71ae5a4SJacob Faibussowitsch { 92995fce210SBarry Smith PetscFunctionBegin; 93095fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 93128b400f6SJacob Faibussowitsch PetscCheck(sf->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFSetUp() before obtaining ranks"); 932dec1416fSJunchao Zhang if (sf->ops->GetRootRanks) { 9339566063dSJacob Faibussowitsch PetscCall((sf->ops->GetRootRanks)(sf, nranks, ranks, roffset, rmine, rremote)); 934dec1416fSJunchao Zhang } else { 935dec1416fSJunchao Zhang /* The generic implementation */ 93695fce210SBarry Smith if (nranks) *nranks = sf->nranks; 93795fce210SBarry Smith if (ranks) *ranks = sf->ranks; 93895fce210SBarry Smith if (roffset) *roffset = sf->roffset; 93995fce210SBarry Smith if (rmine) *rmine = sf->rmine; 94095fce210SBarry Smith if (rremote) *rremote = sf->rremote; 941dec1416fSJunchao Zhang } 9423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 94395fce210SBarry Smith } 94495fce210SBarry Smith 9458750ddebSJunchao Zhang /*@C 9468750ddebSJunchao Zhang PetscSFGetLeafRanks - Get leaf ranks referencing roots on this process 9478750ddebSJunchao Zhang 9488750ddebSJunchao Zhang Not Collective 9498750ddebSJunchao Zhang 9504165533cSJose E. Roman Input Parameter: 9518750ddebSJunchao Zhang . sf - star forest 9528750ddebSJunchao Zhang 9534165533cSJose E. Roman Output Parameters: 9548750ddebSJunchao Zhang + niranks - number of leaf ranks referencing roots on this process 95520662ed9SBarry Smith . iranks - [`niranks`] array of ranks 95620662ed9SBarry Smith . ioffset - [`niranks`+1] offset in `irootloc` for each rank 95720662ed9SBarry Smith - irootloc - [`ioffset`[`niranks`]] concatenated array holding local indices of roots referenced by each leaf rank 9588750ddebSJunchao Zhang 9598750ddebSJunchao Zhang Level: developer 9608750ddebSJunchao Zhang 961cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGetRootRanks()` 9628750ddebSJunchao Zhang @*/ 963d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetLeafRanks(PetscSF sf, PetscInt *niranks, const PetscMPIInt **iranks, const PetscInt **ioffset, const PetscInt **irootloc) 964d71ae5a4SJacob Faibussowitsch { 9658750ddebSJunchao Zhang PetscFunctionBegin; 9668750ddebSJunchao Zhang PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 96728b400f6SJacob Faibussowitsch PetscCheck(sf->setupcalled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFSetUp() before obtaining ranks"); 9688750ddebSJunchao Zhang if (sf->ops->GetLeafRanks) { 9699566063dSJacob Faibussowitsch PetscCall((sf->ops->GetLeafRanks)(sf, niranks, iranks, ioffset, irootloc)); 9708750ddebSJunchao Zhang } else { 9718750ddebSJunchao Zhang PetscSFType type; 9729566063dSJacob Faibussowitsch PetscCall(PetscSFGetType(sf, &type)); 97398921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "PetscSFGetLeafRanks() is not supported on this StarForest type: %s", type); 9748750ddebSJunchao Zhang } 9753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9768750ddebSJunchao Zhang } 9778750ddebSJunchao Zhang 978d71ae5a4SJacob Faibussowitsch static PetscBool InList(PetscMPIInt needle, PetscMPIInt n, const PetscMPIInt *list) 979d71ae5a4SJacob Faibussowitsch { 980b5a8e515SJed Brown PetscInt i; 981b5a8e515SJed Brown for (i = 0; i < n; i++) { 982b5a8e515SJed Brown if (needle == list[i]) return PETSC_TRUE; 983b5a8e515SJed Brown } 984b5a8e515SJed Brown return PETSC_FALSE; 985b5a8e515SJed Brown } 986b5a8e515SJed Brown 98795fce210SBarry Smith /*@C 988cab54364SBarry Smith PetscSFSetUpRanks - Set up data structures associated with ranks; this is for internal use by `PetscSF` implementations. 98921c688dcSJed Brown 99021c688dcSJed Brown Collective 99121c688dcSJed Brown 9924165533cSJose E. Roman Input Parameters: 993cab54364SBarry Smith + sf - `PetscSF` to set up; `PetscSFSetGraph()` must have been called 994cab54364SBarry Smith - dgroup - `MPI_Group` of ranks to be distinguished (e.g., for self or shared memory exchange) 99521c688dcSJed Brown 99621c688dcSJed Brown Level: developer 99721c688dcSJed Brown 998cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGetRootRanks()` 99921c688dcSJed Brown @*/ 1000d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFSetUpRanks(PetscSF sf, MPI_Group dgroup) 1001d71ae5a4SJacob Faibussowitsch { 1002eec179cfSJacob Faibussowitsch PetscHMapI table; 1003eec179cfSJacob Faibussowitsch PetscHashIter pos; 1004b5a8e515SJed Brown PetscMPIInt size, groupsize, *groupranks; 1005247e8311SStefano Zampini PetscInt *rcount, *ranks; 1006247e8311SStefano Zampini PetscInt i, irank = -1, orank = -1; 100721c688dcSJed Brown 100821c688dcSJed Brown PetscFunctionBegin; 100921c688dcSJed Brown PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 101029046d53SLisandro Dalcin PetscSFCheckGraphSet(sf, 1); 10119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sf), &size)); 1012eec179cfSJacob Faibussowitsch PetscCall(PetscHMapICreateWithSize(10, &table)); 101321c688dcSJed Brown for (i = 0; i < sf->nleaves; i++) { 101421c688dcSJed Brown /* Log 1-based rank */ 1015eec179cfSJacob Faibussowitsch PetscCall(PetscHMapISetWithMode(table, sf->remote[i].rank + 1, 1, ADD_VALUES)); 101621c688dcSJed Brown } 1017eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetSize(table, &sf->nranks)); 10189566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(sf->nranks, &sf->ranks, sf->nranks + 1, &sf->roffset, sf->nleaves, &sf->rmine, sf->nleaves, &sf->rremote)); 10199566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(sf->nranks, &rcount, sf->nranks, &ranks)); 1020eec179cfSJacob Faibussowitsch PetscHashIterBegin(table, pos); 102121c688dcSJed Brown for (i = 0; i < sf->nranks; i++) { 1022eec179cfSJacob Faibussowitsch PetscHashIterGetKey(table, pos, ranks[i]); 1023eec179cfSJacob Faibussowitsch PetscHashIterGetVal(table, pos, rcount[i]); 1024eec179cfSJacob Faibussowitsch PetscHashIterNext(table, pos); 102521c688dcSJed Brown ranks[i]--; /* Convert back to 0-based */ 102621c688dcSJed Brown } 1027eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&table)); 1028b5a8e515SJed Brown 1029b5a8e515SJed Brown /* We expect that dgroup is reliably "small" while nranks could be large */ 1030b5a8e515SJed Brown { 10317fb8a5e4SKarl Rupp MPI_Group group = MPI_GROUP_NULL; 1032b5a8e515SJed Brown PetscMPIInt *dgroupranks; 10339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)sf), &group)); 10349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_size(dgroup, &groupsize)); 10359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(groupsize, &dgroupranks)); 10369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(groupsize, &groupranks)); 1037b5a8e515SJed Brown for (i = 0; i < groupsize; i++) dgroupranks[i] = i; 10389566063dSJacob Faibussowitsch if (groupsize) PetscCallMPI(MPI_Group_translate_ranks(dgroup, groupsize, dgroupranks, group, groupranks)); 10399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&group)); 10409566063dSJacob Faibussowitsch PetscCall(PetscFree(dgroupranks)); 1041b5a8e515SJed Brown } 1042b5a8e515SJed Brown 1043b5a8e515SJed Brown /* Partition ranks[] into distinguished (first sf->ndranks) followed by non-distinguished */ 1044b5a8e515SJed Brown for (sf->ndranks = 0, i = sf->nranks; sf->ndranks < i;) { 1045b5a8e515SJed Brown for (i--; sf->ndranks < i; i--) { /* Scan i backward looking for distinguished rank */ 1046b5a8e515SJed Brown if (InList(ranks[i], groupsize, groupranks)) break; 1047b5a8e515SJed Brown } 1048b5a8e515SJed Brown for (; sf->ndranks <= i; sf->ndranks++) { /* Scan sf->ndranks forward looking for non-distinguished rank */ 1049b5a8e515SJed Brown if (!InList(ranks[sf->ndranks], groupsize, groupranks)) break; 1050b5a8e515SJed Brown } 1051b5a8e515SJed Brown if (sf->ndranks < i) { /* Swap ranks[sf->ndranks] with ranks[i] */ 1052b5a8e515SJed Brown PetscInt tmprank, tmpcount; 1053247e8311SStefano Zampini 1054b5a8e515SJed Brown tmprank = ranks[i]; 1055b5a8e515SJed Brown tmpcount = rcount[i]; 1056b5a8e515SJed Brown ranks[i] = ranks[sf->ndranks]; 1057b5a8e515SJed Brown rcount[i] = rcount[sf->ndranks]; 1058b5a8e515SJed Brown ranks[sf->ndranks] = tmprank; 1059b5a8e515SJed Brown rcount[sf->ndranks] = tmpcount; 1060b5a8e515SJed Brown sf->ndranks++; 1061b5a8e515SJed Brown } 1062b5a8e515SJed Brown } 10639566063dSJacob Faibussowitsch PetscCall(PetscFree(groupranks)); 10649566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithArray(sf->ndranks, ranks, rcount)); 10659566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithArray(sf->nranks - sf->ndranks, ranks + sf->ndranks, rcount + sf->ndranks)); 106621c688dcSJed Brown sf->roffset[0] = 0; 106721c688dcSJed Brown for (i = 0; i < sf->nranks; i++) { 10689566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(ranks[i], sf->ranks + i)); 106921c688dcSJed Brown sf->roffset[i + 1] = sf->roffset[i] + rcount[i]; 107021c688dcSJed Brown rcount[i] = 0; 107121c688dcSJed Brown } 1072247e8311SStefano Zampini for (i = 0, irank = -1, orank = -1; i < sf->nleaves; i++) { 1073247e8311SStefano Zampini /* short circuit */ 1074247e8311SStefano Zampini if (orank != sf->remote[i].rank) { 107521c688dcSJed Brown /* Search for index of iremote[i].rank in sf->ranks */ 10769566063dSJacob Faibussowitsch PetscCall(PetscFindMPIInt(sf->remote[i].rank, sf->ndranks, sf->ranks, &irank)); 1077b5a8e515SJed Brown if (irank < 0) { 10789566063dSJacob Faibussowitsch PetscCall(PetscFindMPIInt(sf->remote[i].rank, sf->nranks - sf->ndranks, sf->ranks + sf->ndranks, &irank)); 1079b5a8e515SJed Brown if (irank >= 0) irank += sf->ndranks; 108021c688dcSJed Brown } 1081247e8311SStefano Zampini orank = sf->remote[i].rank; 1082247e8311SStefano Zampini } 108308401ef6SPierre Jolivet PetscCheck(irank >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find rank %" PetscInt_FMT " in array", sf->remote[i].rank); 108421c688dcSJed Brown sf->rmine[sf->roffset[irank] + rcount[irank]] = sf->mine ? sf->mine[i] : i; 108521c688dcSJed Brown sf->rremote[sf->roffset[irank] + rcount[irank]] = sf->remote[i].index; 108621c688dcSJed Brown rcount[irank]++; 108721c688dcSJed Brown } 10889566063dSJacob Faibussowitsch PetscCall(PetscFree2(rcount, ranks)); 10893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 109021c688dcSJed Brown } 109121c688dcSJed Brown 109221c688dcSJed Brown /*@C 109395fce210SBarry Smith PetscSFGetGroups - gets incoming and outgoing process groups 109495fce210SBarry Smith 109595fce210SBarry Smith Collective 109695fce210SBarry Smith 10974165533cSJose E. Roman Input Parameter: 109895fce210SBarry Smith . sf - star forest 109995fce210SBarry Smith 11004165533cSJose E. Roman Output Parameters: 110195fce210SBarry Smith + incoming - group of origin processes for incoming edges (leaves that reference my roots) 110295fce210SBarry Smith - outgoing - group of destination processes for outgoing edges (roots that I reference) 110395fce210SBarry Smith 110495fce210SBarry Smith Level: developer 110595fce210SBarry Smith 1106cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGetWindow()`, `PetscSFRestoreWindow()` 110795fce210SBarry Smith @*/ 1108d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetGroups(PetscSF sf, MPI_Group *incoming, MPI_Group *outgoing) 1109d71ae5a4SJacob Faibussowitsch { 11107fb8a5e4SKarl Rupp MPI_Group group = MPI_GROUP_NULL; 111195fce210SBarry Smith 111295fce210SBarry Smith PetscFunctionBegin; 111308401ef6SPierre Jolivet PetscCheck(sf->nranks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFSetUpRanks() before obtaining groups"); 111495fce210SBarry Smith if (sf->ingroup == MPI_GROUP_NULL) { 111595fce210SBarry Smith PetscInt i; 111695fce210SBarry Smith const PetscInt *indegree; 111795fce210SBarry Smith PetscMPIInt rank, *outranks, *inranks; 111895fce210SBarry Smith PetscSFNode *remote; 111995fce210SBarry Smith PetscSF bgcount; 112095fce210SBarry Smith 112195fce210SBarry Smith /* Compute the number of incoming ranks */ 11229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(sf->nranks, &remote)); 112395fce210SBarry Smith for (i = 0; i < sf->nranks; i++) { 112495fce210SBarry Smith remote[i].rank = sf->ranks[i]; 112595fce210SBarry Smith remote[i].index = 0; 112695fce210SBarry Smith } 11279566063dSJacob Faibussowitsch PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_CONFONLY, &bgcount)); 11289566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(bgcount, 1, sf->nranks, NULL, PETSC_COPY_VALUES, remote, PETSC_OWN_POINTER)); 11299566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(bgcount, &indegree)); 11309566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(bgcount, &indegree)); 113195fce210SBarry Smith /* Enumerate the incoming ranks */ 11329566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(indegree[0], &inranks, sf->nranks, &outranks)); 11339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank)); 113495fce210SBarry Smith for (i = 0; i < sf->nranks; i++) outranks[i] = rank; 11359566063dSJacob Faibussowitsch PetscCall(PetscSFGatherBegin(bgcount, MPI_INT, outranks, inranks)); 11369566063dSJacob Faibussowitsch PetscCall(PetscSFGatherEnd(bgcount, MPI_INT, outranks, inranks)); 11379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)sf), &group)); 11389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_incl(group, indegree[0], inranks, &sf->ingroup)); 11399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&group)); 11409566063dSJacob Faibussowitsch PetscCall(PetscFree2(inranks, outranks)); 11419566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&bgcount)); 114295fce210SBarry Smith } 114395fce210SBarry Smith *incoming = sf->ingroup; 114495fce210SBarry Smith 114595fce210SBarry Smith if (sf->outgroup == MPI_GROUP_NULL) { 11469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(PetscObjectComm((PetscObject)sf), &group)); 11479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_incl(group, sf->nranks, sf->ranks, &sf->outgroup)); 11489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&group)); 114995fce210SBarry Smith } 115095fce210SBarry Smith *outgoing = sf->outgroup; 11513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115295fce210SBarry Smith } 115395fce210SBarry Smith 115429046d53SLisandro Dalcin /*@ 1155cab54364SBarry Smith PetscSFGetMultiSF - gets the inner `PetscSF` implementing gathers and scatters 115695fce210SBarry Smith 115795fce210SBarry Smith Collective 115895fce210SBarry Smith 11594165533cSJose E. Roman Input Parameter: 116095fce210SBarry Smith . sf - star forest that may contain roots with 0 or with more than 1 vertex 116195fce210SBarry Smith 11624165533cSJose E. Roman Output Parameter: 116395fce210SBarry Smith . multi - star forest with split roots, such that each root has degree exactly 1 116495fce210SBarry Smith 116595fce210SBarry Smith Level: developer 116695fce210SBarry Smith 1167cab54364SBarry Smith Note: 1168cab54364SBarry Smith In most cases, users should use `PetscSFGatherBegin()` and `PetscSFScatterBegin()` instead of manipulating multi 116995fce210SBarry Smith directly. Since multi satisfies the stronger condition that each entry in the global space has exactly one incoming 117095fce210SBarry Smith edge, it is a candidate for future optimization that might involve its removal. 117195fce210SBarry Smith 1172cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFGatherBegin()`, `PetscSFScatterBegin()`, `PetscSFComputeMultiRootOriginalNumbering()` 117395fce210SBarry Smith @*/ 1174d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetMultiSF(PetscSF sf, PetscSF *multi) 1175d71ae5a4SJacob Faibussowitsch { 117695fce210SBarry Smith PetscFunctionBegin; 117795fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 117895fce210SBarry Smith PetscValidPointer(multi, 2); 117995fce210SBarry Smith if (sf->nroots < 0) { /* Graph has not been set yet; why do we need this? */ 11809566063dSJacob Faibussowitsch PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, &sf->multi)); 118195fce210SBarry Smith *multi = sf->multi; 1182013b3241SStefano Zampini sf->multi->multi = sf->multi; 11833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 118495fce210SBarry Smith } 118595fce210SBarry Smith if (!sf->multi) { 118695fce210SBarry Smith const PetscInt *indegree; 11879837ea96SMatthew G. Knepley PetscInt i, *inoffset, *outones, *outoffset, maxlocal; 118895fce210SBarry Smith PetscSFNode *remote; 118929046d53SLisandro Dalcin maxlocal = sf->maxleaf + 1; /* TODO: We should use PetscSFGetLeafRange() */ 11909566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sf, &indegree)); 11919566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sf, &indegree)); 11929566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(sf->nroots + 1, &inoffset, maxlocal, &outones, maxlocal, &outoffset)); 119395fce210SBarry Smith inoffset[0] = 0; 119495fce210SBarry Smith for (i = 0; i < sf->nroots; i++) inoffset[i + 1] = inoffset[i] + indegree[i]; 11959837ea96SMatthew G. Knepley for (i = 0; i < maxlocal; i++) outones[i] = 1; 11969566063dSJacob Faibussowitsch PetscCall(PetscSFFetchAndOpBegin(sf, MPIU_INT, inoffset, outones, outoffset, MPI_SUM)); 11979566063dSJacob Faibussowitsch PetscCall(PetscSFFetchAndOpEnd(sf, MPIU_INT, inoffset, outones, outoffset, MPI_SUM)); 119895fce210SBarry Smith for (i = 0; i < sf->nroots; i++) inoffset[i] -= indegree[i]; /* Undo the increment */ 119976bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { /* Check that the expected number of increments occurred */ 1200ad540459SPierre Jolivet for (i = 0; i < sf->nroots; i++) PetscCheck(inoffset[i] + indegree[i] == inoffset[i + 1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrect result after PetscSFFetchAndOp"); 120176bd3646SJed Brown } 12029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(sf->nleaves, &remote)); 120395fce210SBarry Smith for (i = 0; i < sf->nleaves; i++) { 120495fce210SBarry Smith remote[i].rank = sf->remote[i].rank; 120538e7336fSToby Isaac remote[i].index = outoffset[sf->mine ? sf->mine[i] : i]; 120695fce210SBarry Smith } 12079566063dSJacob Faibussowitsch PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_RANKS, &sf->multi)); 1208013b3241SStefano Zampini sf->multi->multi = sf->multi; 12099566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf->multi, inoffset[sf->nroots], sf->nleaves, sf->mine, PETSC_COPY_VALUES, remote, PETSC_OWN_POINTER)); 121095fce210SBarry Smith if (sf->rankorder) { /* Sort the ranks */ 121195fce210SBarry Smith PetscMPIInt rank; 121295fce210SBarry Smith PetscInt *inranks, *newoffset, *outranks, *newoutoffset, *tmpoffset, maxdegree; 121395fce210SBarry Smith PetscSFNode *newremote; 12149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank)); 121595fce210SBarry Smith for (i = 0, maxdegree = 0; i < sf->nroots; i++) maxdegree = PetscMax(maxdegree, indegree[i]); 12169566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(sf->multi->nroots, &inranks, sf->multi->nroots, &newoffset, maxlocal, &outranks, maxlocal, &newoutoffset, maxdegree, &tmpoffset)); 12179837ea96SMatthew G. Knepley for (i = 0; i < maxlocal; i++) outranks[i] = rank; 12189566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf->multi, MPIU_INT, outranks, inranks, MPI_REPLACE)); 12199566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf->multi, MPIU_INT, outranks, inranks, MPI_REPLACE)); 122095fce210SBarry Smith /* Sort the incoming ranks at each vertex, build the inverse map */ 122195fce210SBarry Smith for (i = 0; i < sf->nroots; i++) { 122295fce210SBarry Smith PetscInt j; 122395fce210SBarry Smith for (j = 0; j < indegree[i]; j++) tmpoffset[j] = j; 12249566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithArray(indegree[i], inranks + inoffset[i], tmpoffset)); 122595fce210SBarry Smith for (j = 0; j < indegree[i]; j++) newoffset[inoffset[i] + tmpoffset[j]] = inoffset[i] + j; 122695fce210SBarry Smith } 12279566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf->multi, MPIU_INT, newoffset, newoutoffset, MPI_REPLACE)); 12289566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf->multi, MPIU_INT, newoffset, newoutoffset, MPI_REPLACE)); 12299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(sf->nleaves, &newremote)); 123095fce210SBarry Smith for (i = 0; i < sf->nleaves; i++) { 123195fce210SBarry Smith newremote[i].rank = sf->remote[i].rank; 123201365b40SToby Isaac newremote[i].index = newoutoffset[sf->mine ? sf->mine[i] : i]; 123395fce210SBarry Smith } 12349566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf->multi, inoffset[sf->nroots], sf->nleaves, sf->mine, PETSC_COPY_VALUES, newremote, PETSC_OWN_POINTER)); 12359566063dSJacob Faibussowitsch PetscCall(PetscFree5(inranks, newoffset, outranks, newoutoffset, tmpoffset)); 123695fce210SBarry Smith } 12379566063dSJacob Faibussowitsch PetscCall(PetscFree3(inoffset, outones, outoffset)); 123895fce210SBarry Smith } 123995fce210SBarry Smith *multi = sf->multi; 12403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 124195fce210SBarry Smith } 124295fce210SBarry Smith 124395fce210SBarry Smith /*@C 124420662ed9SBarry Smith PetscSFCreateEmbeddedRootSF - removes edges from all but the selected roots of a `PetscSF`, does not remap indices 124595fce210SBarry Smith 124695fce210SBarry Smith Collective 124795fce210SBarry Smith 12484165533cSJose E. Roman Input Parameters: 124995fce210SBarry Smith + sf - original star forest 1250ba2a7774SJunchao Zhang . nselected - number of selected roots on this process 1251ba2a7774SJunchao Zhang - selected - indices of the selected roots on this process 125295fce210SBarry Smith 12534165533cSJose E. Roman Output Parameter: 1254cd620004SJunchao Zhang . esf - new star forest 125595fce210SBarry Smith 125695fce210SBarry Smith Level: advanced 125795fce210SBarry Smith 125895fce210SBarry Smith Note: 1259cab54364SBarry Smith To use the new `PetscSF`, it may be necessary to know the indices of the leaves that are still participating. This can 126095fce210SBarry Smith be done by calling PetscSFGetGraph(). 126195fce210SBarry Smith 1262cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 126395fce210SBarry Smith @*/ 1264d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateEmbeddedRootSF(PetscSF sf, PetscInt nselected, const PetscInt *selected, PetscSF *esf) 1265d71ae5a4SJacob Faibussowitsch { 1266cd620004SJunchao Zhang PetscInt i, j, n, nroots, nleaves, esf_nleaves, *new_ilocal, minleaf, maxleaf, maxlocal; 1267cd620004SJunchao Zhang const PetscInt *ilocal; 1268cd620004SJunchao Zhang signed char *rootdata, *leafdata, *leafmem; 1269ba2a7774SJunchao Zhang const PetscSFNode *iremote; 1270f659e5c7SJunchao Zhang PetscSFNode *new_iremote; 1271f659e5c7SJunchao Zhang MPI_Comm comm; 127295fce210SBarry Smith 127395fce210SBarry Smith PetscFunctionBegin; 127495fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 127529046d53SLisandro Dalcin PetscSFCheckGraphSet(sf, 1); 1276dadcf809SJacob Faibussowitsch if (nselected) PetscValidIntPointer(selected, 3); 1277cd620004SJunchao Zhang PetscValidPointer(esf, 4); 12780511a646SMatthew G. Knepley 12799566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 12809566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_EmbedSF, sf, 0, 0, 0)); 12819566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 12829566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, &iremote)); 1283cd620004SJunchao Zhang 128476bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { /* Error out if selected[] has dups or out of range indices */ 1285cd620004SJunchao Zhang PetscBool dups; 12869566063dSJacob Faibussowitsch PetscCall(PetscCheckDupsInt(nselected, selected, &dups)); 128728b400f6SJacob Faibussowitsch PetscCheck(!dups, comm, PETSC_ERR_ARG_WRONG, "selected[] has dups"); 12889371c9d4SSatish Balay for (i = 0; i < nselected; i++) PetscCheck(selected[i] >= 0 && selected[i] < nroots, comm, PETSC_ERR_ARG_OUTOFRANGE, "selected root indice %" PetscInt_FMT " is out of [0,%" PetscInt_FMT ")", selected[i], nroots); 1289cd620004SJunchao Zhang } 1290f659e5c7SJunchao Zhang 1291dbbe0bcdSBarry Smith if (sf->ops->CreateEmbeddedRootSF) PetscUseTypeMethod(sf, CreateEmbeddedRootSF, nselected, selected, esf); 1292dbbe0bcdSBarry Smith else { 1293cd620004SJunchao Zhang /* A generic version of creating embedded sf */ 12949566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(sf, &minleaf, &maxleaf)); 1295cd620004SJunchao Zhang maxlocal = maxleaf - minleaf + 1; 12969566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(nroots, &rootdata, maxlocal, &leafmem)); 1297cd620004SJunchao Zhang leafdata = leafmem - minleaf; 1298cd620004SJunchao Zhang /* Tag selected roots and bcast to leaves */ 1299cd620004SJunchao Zhang for (i = 0; i < nselected; i++) rootdata[selected[i]] = 1; 13009566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPI_SIGNED_CHAR, rootdata, leafdata, MPI_REPLACE)); 13019566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPI_SIGNED_CHAR, rootdata, leafdata, MPI_REPLACE)); 1302ba2a7774SJunchao Zhang 1303cd620004SJunchao Zhang /* Build esf with leaves that are still connected */ 1304cd620004SJunchao Zhang esf_nleaves = 0; 1305cd620004SJunchao Zhang for (i = 0; i < nleaves; i++) { 1306cd620004SJunchao Zhang j = ilocal ? ilocal[i] : i; 1307cd620004SJunchao Zhang /* esf_nleaves += leafdata[j] should work in theory, but failed with SFWindow bugs 1308cd620004SJunchao Zhang with PetscSFBcast. See https://gitlab.com/petsc/petsc/issues/555 1309cd620004SJunchao Zhang */ 1310cd620004SJunchao Zhang esf_nleaves += (leafdata[j] ? 1 : 0); 1311cd620004SJunchao Zhang } 13129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(esf_nleaves, &new_ilocal)); 13139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(esf_nleaves, &new_iremote)); 1314cd620004SJunchao Zhang for (i = n = 0; i < nleaves; i++) { 1315cd620004SJunchao Zhang j = ilocal ? ilocal[i] : i; 1316cd620004SJunchao Zhang if (leafdata[j]) { 1317cd620004SJunchao Zhang new_ilocal[n] = j; 1318cd620004SJunchao Zhang new_iremote[n].rank = iremote[i].rank; 1319cd620004SJunchao Zhang new_iremote[n].index = iremote[i].index; 1320fc1ede2bSMatthew G. Knepley ++n; 132195fce210SBarry Smith } 132295fce210SBarry Smith } 13239566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, esf)); 13249566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(*esf)); 13259566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*esf, nroots, esf_nleaves, new_ilocal, PETSC_OWN_POINTER, new_iremote, PETSC_OWN_POINTER)); 13269566063dSJacob Faibussowitsch PetscCall(PetscFree2(rootdata, leafmem)); 1327f659e5c7SJunchao Zhang } 13289566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_EmbedSF, sf, 0, 0, 0)); 13293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 133095fce210SBarry Smith } 133195fce210SBarry Smith 13322f5fb4c2SMatthew G. Knepley /*@C 133320662ed9SBarry Smith PetscSFCreateEmbeddedLeafSF - removes edges from all but the selected leaves of a `PetscSF`, does not remap indices 13342f5fb4c2SMatthew G. Knepley 13352f5fb4c2SMatthew G. Knepley Collective 13362f5fb4c2SMatthew G. Knepley 13374165533cSJose E. Roman Input Parameters: 13382f5fb4c2SMatthew G. Knepley + sf - original star forest 1339f659e5c7SJunchao Zhang . nselected - number of selected leaves on this process 1340f659e5c7SJunchao Zhang - selected - indices of the selected leaves on this process 13412f5fb4c2SMatthew G. Knepley 13424165533cSJose E. Roman Output Parameter: 13432f5fb4c2SMatthew G. Knepley . newsf - new star forest 13442f5fb4c2SMatthew G. Knepley 13452f5fb4c2SMatthew G. Knepley Level: advanced 13462f5fb4c2SMatthew G. Knepley 1347cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreateEmbeddedRootSF()`, `PetscSFSetGraph()`, `PetscSFGetGraph()` 13482f5fb4c2SMatthew G. Knepley @*/ 1349d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateEmbeddedLeafSF(PetscSF sf, PetscInt nselected, const PetscInt *selected, PetscSF *newsf) 1350d71ae5a4SJacob Faibussowitsch { 1351f659e5c7SJunchao Zhang const PetscSFNode *iremote; 1352f659e5c7SJunchao Zhang PetscSFNode *new_iremote; 1353f659e5c7SJunchao Zhang const PetscInt *ilocal; 1354f659e5c7SJunchao Zhang PetscInt i, nroots, *leaves, *new_ilocal; 1355f659e5c7SJunchao Zhang MPI_Comm comm; 13562f5fb4c2SMatthew G. Knepley 13572f5fb4c2SMatthew G. Knepley PetscFunctionBegin; 13582f5fb4c2SMatthew G. Knepley PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 135929046d53SLisandro Dalcin PetscSFCheckGraphSet(sf, 1); 1360dadcf809SJacob Faibussowitsch if (nselected) PetscValidIntPointer(selected, 3); 13612f5fb4c2SMatthew G. Knepley PetscValidPointer(newsf, 4); 13622f5fb4c2SMatthew G. Knepley 1363f659e5c7SJunchao Zhang /* Uniq selected[] and put results in leaves[] */ 13649566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 13659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nselected, &leaves)); 13669566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(leaves, selected, nselected)); 13679566063dSJacob Faibussowitsch PetscCall(PetscSortedRemoveDupsInt(&nselected, leaves)); 136808401ef6SPierre 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); 1369f659e5c7SJunchao Zhang 1370f659e5c7SJunchao Zhang /* Optimize the routine only when sf is setup and hence we can reuse sf's communication pattern */ 1371dbbe0bcdSBarry Smith if (sf->setupcalled && sf->ops->CreateEmbeddedLeafSF) PetscUseTypeMethod(sf, CreateEmbeddedLeafSF, nselected, leaves, newsf); 1372dbbe0bcdSBarry Smith else { 13739566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, NULL, &ilocal, &iremote)); 13749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nselected, &new_ilocal)); 13759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nselected, &new_iremote)); 1376f659e5c7SJunchao Zhang for (i = 0; i < nselected; ++i) { 1377f659e5c7SJunchao Zhang const PetscInt l = leaves[i]; 1378f659e5c7SJunchao Zhang new_ilocal[i] = ilocal ? ilocal[l] : l; 1379f659e5c7SJunchao Zhang new_iremote[i].rank = iremote[l].rank; 1380f659e5c7SJunchao Zhang new_iremote[i].index = iremote[l].index; 13812f5fb4c2SMatthew G. Knepley } 13829566063dSJacob Faibussowitsch PetscCall(PetscSFDuplicate(sf, PETSCSF_DUPLICATE_CONFONLY, newsf)); 13839566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*newsf, nroots, nselected, new_ilocal, PETSC_OWN_POINTER, new_iremote, PETSC_OWN_POINTER)); 1384f659e5c7SJunchao Zhang } 13859566063dSJacob Faibussowitsch PetscCall(PetscFree(leaves)); 13863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13872f5fb4c2SMatthew G. Knepley } 13882f5fb4c2SMatthew G. Knepley 138995fce210SBarry Smith /*@C 1390cab54364SBarry Smith PetscSFBcastBegin - begin pointwise broadcast with root value being reduced to leaf value, to be concluded with call to `PetscSFBcastEnd()` 13913482bfa8SJunchao Zhang 1392c3339decSBarry Smith Collective 13933482bfa8SJunchao Zhang 13944165533cSJose E. Roman Input Parameters: 13953482bfa8SJunchao Zhang + sf - star forest on which to communicate 13963482bfa8SJunchao Zhang . unit - data type associated with each node 13973482bfa8SJunchao Zhang . rootdata - buffer to broadcast 13983482bfa8SJunchao Zhang - op - operation to use for reduction 13993482bfa8SJunchao Zhang 14004165533cSJose E. Roman Output Parameter: 14013482bfa8SJunchao Zhang . leafdata - buffer to be reduced with values from each leaf's respective root 14023482bfa8SJunchao Zhang 14033482bfa8SJunchao Zhang Level: intermediate 14043482bfa8SJunchao Zhang 140520662ed9SBarry Smith Note: 140620662ed9SBarry Smith When PETSc is configured with device support, it will use its own mechanism to figure out whether the given data pointers 1407da81f932SPierre Jolivet are host pointers or device pointers, which may incur a noticeable cost. If you already knew the info, you should 1408cab54364SBarry Smith use `PetscSFBcastWithMemTypeBegin()` instead. 1409cab54364SBarry Smith 1410cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFBcastEnd()`, `PetscSFBcastWithMemTypeBegin()` 14113482bfa8SJunchao Zhang @*/ 1412d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFBcastBegin(PetscSF sf, MPI_Datatype unit, const void *rootdata, void *leafdata, MPI_Op op) 1413d71ae5a4SJacob Faibussowitsch { 1414eb02082bSJunchao Zhang PetscMemType rootmtype, leafmtype; 14153482bfa8SJunchao Zhang 14163482bfa8SJunchao Zhang PetscFunctionBegin; 14173482bfa8SJunchao Zhang PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 14189566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 14199566063dSJacob Faibussowitsch if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_BcastBegin, sf, 0, 0, 0)); 14209566063dSJacob Faibussowitsch PetscCall(PetscGetMemType(rootdata, &rootmtype)); 14219566063dSJacob Faibussowitsch PetscCall(PetscGetMemType(leafdata, &leafmtype)); 1422dbbe0bcdSBarry Smith PetscUseTypeMethod(sf, BcastBegin, unit, rootmtype, rootdata, leafmtype, leafdata, op); 14239566063dSJacob Faibussowitsch if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_BcastBegin, sf, 0, 0, 0)); 14243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14253482bfa8SJunchao Zhang } 14263482bfa8SJunchao Zhang 14273482bfa8SJunchao Zhang /*@C 142820662ed9SBarry Smith PetscSFBcastWithMemTypeBegin - begin pointwise broadcast with root value being reduced to leaf value with explicit memory types, to be concluded with call 142920662ed9SBarry Smith to `PetscSFBcastEnd()` 1430d0295fc0SJunchao Zhang 1431c3339decSBarry Smith Collective 1432d0295fc0SJunchao Zhang 14334165533cSJose E. Roman Input Parameters: 1434d0295fc0SJunchao Zhang + sf - star forest on which to communicate 1435d0295fc0SJunchao Zhang . unit - data type associated with each node 1436d0295fc0SJunchao Zhang . rootmtype - memory type of rootdata 1437d0295fc0SJunchao Zhang . rootdata - buffer to broadcast 1438d0295fc0SJunchao Zhang . leafmtype - memory type of leafdata 1439d0295fc0SJunchao Zhang - op - operation to use for reduction 1440d0295fc0SJunchao Zhang 14414165533cSJose E. Roman Output Parameter: 1442d0295fc0SJunchao Zhang . leafdata - buffer to be reduced with values from each leaf's respective root 1443d0295fc0SJunchao Zhang 1444d0295fc0SJunchao Zhang Level: intermediate 1445d0295fc0SJunchao Zhang 1446cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFBcastEnd()`, `PetscSFBcastBegin()` 1447d0295fc0SJunchao Zhang @*/ 1448d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFBcastWithMemTypeBegin(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, const void *rootdata, PetscMemType leafmtype, void *leafdata, MPI_Op op) 1449d71ae5a4SJacob Faibussowitsch { 1450d0295fc0SJunchao Zhang PetscFunctionBegin; 1451d0295fc0SJunchao Zhang PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 14529566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 14539566063dSJacob Faibussowitsch if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_BcastBegin, sf, 0, 0, 0)); 1454dbbe0bcdSBarry Smith PetscUseTypeMethod(sf, BcastBegin, unit, rootmtype, rootdata, leafmtype, leafdata, op); 14559566063dSJacob Faibussowitsch if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_BcastBegin, sf, 0, 0, 0)); 14563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1457d0295fc0SJunchao Zhang } 1458d0295fc0SJunchao Zhang 1459d0295fc0SJunchao Zhang /*@C 146020662ed9SBarry Smith PetscSFBcastEnd - end a broadcast and reduce operation started with `PetscSFBcastBegin()` or `PetscSFBcastWithMemTypeBegin()` 14613482bfa8SJunchao Zhang 14623482bfa8SJunchao Zhang Collective 14633482bfa8SJunchao Zhang 14644165533cSJose E. Roman Input Parameters: 14653482bfa8SJunchao Zhang + sf - star forest 14663482bfa8SJunchao Zhang . unit - data type 14673482bfa8SJunchao Zhang . rootdata - buffer to broadcast 14683482bfa8SJunchao Zhang - op - operation to use for reduction 14693482bfa8SJunchao Zhang 14704165533cSJose E. Roman Output Parameter: 14713482bfa8SJunchao Zhang . leafdata - buffer to be reduced with values from each leaf's respective root 14723482bfa8SJunchao Zhang 14733482bfa8SJunchao Zhang Level: intermediate 14743482bfa8SJunchao Zhang 1475cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFReduceEnd()` 14763482bfa8SJunchao Zhang @*/ 1477d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFBcastEnd(PetscSF sf, MPI_Datatype unit, const void *rootdata, void *leafdata, MPI_Op op) 1478d71ae5a4SJacob Faibussowitsch { 14793482bfa8SJunchao Zhang PetscFunctionBegin; 14803482bfa8SJunchao Zhang PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 14819566063dSJacob Faibussowitsch if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_BcastEnd, sf, 0, 0, 0)); 1482dbbe0bcdSBarry Smith PetscUseTypeMethod(sf, BcastEnd, unit, rootdata, leafdata, op); 14839566063dSJacob Faibussowitsch if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_BcastEnd, sf, 0, 0, 0)); 14843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14853482bfa8SJunchao Zhang } 14863482bfa8SJunchao Zhang 14873482bfa8SJunchao Zhang /*@C 1488cab54364SBarry Smith PetscSFReduceBegin - begin reduction of leafdata into rootdata, to be completed with call to `PetscSFReduceEnd()` 148995fce210SBarry Smith 149095fce210SBarry Smith Collective 149195fce210SBarry Smith 14924165533cSJose E. Roman Input Parameters: 149395fce210SBarry Smith + sf - star forest 149495fce210SBarry Smith . unit - data type 149595fce210SBarry Smith . leafdata - values to reduce 149695fce210SBarry Smith - op - reduction operation 149795fce210SBarry Smith 14984165533cSJose E. Roman Output Parameter: 149995fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 150095fce210SBarry Smith 150195fce210SBarry Smith Level: intermediate 150295fce210SBarry Smith 150320662ed9SBarry Smith Note: 150420662ed9SBarry Smith When PETSc is configured with device support, it will use its own mechanism to figure out whether the given data pointers 1505da81f932SPierre Jolivet are host pointers or device pointers, which may incur a noticeable cost. If you already knew the info, you should 1506cab54364SBarry Smith use `PetscSFReduceWithMemTypeBegin()` instead. 1507d0295fc0SJunchao Zhang 150820662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFBcastBegin()`, `PetscSFReduceWithMemTypeBegin()`, `PetscSFReduceEnd()` 150995fce210SBarry Smith @*/ 1510d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFReduceBegin(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *rootdata, MPI_Op op) 1511d71ae5a4SJacob Faibussowitsch { 1512eb02082bSJunchao Zhang PetscMemType rootmtype, leafmtype; 151395fce210SBarry Smith 151495fce210SBarry Smith PetscFunctionBegin; 151595fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 15169566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 15179566063dSJacob Faibussowitsch if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_ReduceBegin, sf, 0, 0, 0)); 15189566063dSJacob Faibussowitsch PetscCall(PetscGetMemType(rootdata, &rootmtype)); 15199566063dSJacob Faibussowitsch PetscCall(PetscGetMemType(leafdata, &leafmtype)); 15209566063dSJacob Faibussowitsch PetscCall((sf->ops->ReduceBegin)(sf, unit, leafmtype, leafdata, rootmtype, rootdata, op)); 15219566063dSJacob Faibussowitsch if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_ReduceBegin, sf, 0, 0, 0)); 15223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 152395fce210SBarry Smith } 152495fce210SBarry Smith 152595fce210SBarry Smith /*@C 1526cab54364SBarry Smith PetscSFReduceWithMemTypeBegin - begin reduction of leafdata into rootdata with explicit memory types, to be completed with call to `PetscSFReduceEnd()` 1527d0295fc0SJunchao Zhang 1528d0295fc0SJunchao Zhang Collective 1529d0295fc0SJunchao Zhang 15304165533cSJose E. Roman Input Parameters: 1531d0295fc0SJunchao Zhang + sf - star forest 1532d0295fc0SJunchao Zhang . unit - data type 1533d0295fc0SJunchao Zhang . leafmtype - memory type of leafdata 1534d0295fc0SJunchao Zhang . leafdata - values to reduce 1535d0295fc0SJunchao Zhang . rootmtype - memory type of rootdata 1536d0295fc0SJunchao Zhang - op - reduction operation 1537d0295fc0SJunchao Zhang 15384165533cSJose E. Roman Output Parameter: 1539d0295fc0SJunchao Zhang . rootdata - result of reduction of values from all leaves of each root 1540d0295fc0SJunchao Zhang 1541d0295fc0SJunchao Zhang Level: intermediate 1542d0295fc0SJunchao Zhang 154320662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFBcastBegin()`, `PetscSFReduceBegin()`, `PetscSFReduceEnd()` 1544d0295fc0SJunchao Zhang @*/ 1545d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFReduceWithMemTypeBegin(PetscSF sf, MPI_Datatype unit, PetscMemType leafmtype, const void *leafdata, PetscMemType rootmtype, void *rootdata, MPI_Op op) 1546d71ae5a4SJacob Faibussowitsch { 1547d0295fc0SJunchao Zhang PetscFunctionBegin; 1548d0295fc0SJunchao Zhang PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 15499566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 15509566063dSJacob Faibussowitsch if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_ReduceBegin, sf, 0, 0, 0)); 15519566063dSJacob Faibussowitsch PetscCall((sf->ops->ReduceBegin)(sf, unit, leafmtype, leafdata, rootmtype, rootdata, op)); 15529566063dSJacob Faibussowitsch if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_ReduceBegin, sf, 0, 0, 0)); 15533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1554d0295fc0SJunchao Zhang } 1555d0295fc0SJunchao Zhang 1556d0295fc0SJunchao Zhang /*@C 155720662ed9SBarry Smith PetscSFReduceEnd - end a reduction operation started with `PetscSFReduceBegin()` or `PetscSFReduceWithMemTypeBegin()` 155895fce210SBarry Smith 155995fce210SBarry Smith Collective 156095fce210SBarry Smith 15614165533cSJose E. Roman Input Parameters: 156295fce210SBarry Smith + sf - star forest 156395fce210SBarry Smith . unit - data type 156495fce210SBarry Smith . leafdata - values to reduce 156595fce210SBarry Smith - op - reduction operation 156695fce210SBarry Smith 15674165533cSJose E. Roman Output Parameter: 156895fce210SBarry Smith . rootdata - result of reduction of values from all leaves of each root 156995fce210SBarry Smith 157095fce210SBarry Smith Level: intermediate 157195fce210SBarry Smith 157220662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFSetGraph()`, `PetscSFBcastEnd()`, `PetscSFReduceBegin()`, `PetscSFReduceWithMemTypeBegin()` 157395fce210SBarry Smith @*/ 1574d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFReduceEnd(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *rootdata, MPI_Op op) 1575d71ae5a4SJacob Faibussowitsch { 157695fce210SBarry Smith PetscFunctionBegin; 157795fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 15789566063dSJacob Faibussowitsch if (!sf->vscat.logging) PetscCall(PetscLogEventBegin(PETSCSF_ReduceEnd, sf, 0, 0, 0)); 1579dbbe0bcdSBarry Smith PetscUseTypeMethod(sf, ReduceEnd, unit, leafdata, rootdata, op); 15809566063dSJacob Faibussowitsch if (!sf->vscat.logging) PetscCall(PetscLogEventEnd(PETSCSF_ReduceEnd, sf, 0, 0, 0)); 15813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 158295fce210SBarry Smith } 158395fce210SBarry Smith 158495fce210SBarry Smith /*@C 1585cab54364SBarry Smith PetscSFFetchAndOpBegin - begin operation that fetches values from root and updates atomically by applying operation using my leaf value, 1586cab54364SBarry Smith to be completed with `PetscSFFetchAndOpEnd()` 1587a1729e3fSJunchao Zhang 1588a1729e3fSJunchao Zhang Collective 1589a1729e3fSJunchao Zhang 15904165533cSJose E. Roman Input Parameters: 1591a1729e3fSJunchao Zhang + sf - star forest 1592a1729e3fSJunchao Zhang . unit - data type 1593a1729e3fSJunchao Zhang . leafdata - leaf values to use in reduction 1594a1729e3fSJunchao Zhang - op - operation to use for reduction 1595a1729e3fSJunchao Zhang 15964165533cSJose E. Roman Output Parameters: 1597a1729e3fSJunchao Zhang + rootdata - root values to be updated, input state is seen by first process to perform an update 1598a1729e3fSJunchao Zhang - leafupdate - state at each leaf's respective root immediately prior to my atomic update 1599a1729e3fSJunchao Zhang 1600a1729e3fSJunchao Zhang Level: advanced 1601a1729e3fSJunchao Zhang 1602a1729e3fSJunchao Zhang Note: 1603a1729e3fSJunchao Zhang The update is only atomic at the granularity provided by the hardware. Different roots referenced by the same process 1604a1729e3fSJunchao Zhang might be updated in a different order. Furthermore, if a composite type is used for the unit datatype, atomicity is 1605a1729e3fSJunchao Zhang not guaranteed across the whole vertex. Therefore, this function is mostly only used with primitive types such as 1606a1729e3fSJunchao Zhang integers. 1607a1729e3fSJunchao Zhang 1608cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeBegin()`, `PetscSFReduceBegin()`, `PetscSFSetGraph()` 1609a1729e3fSJunchao Zhang @*/ 1610d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFFetchAndOpBegin(PetscSF sf, MPI_Datatype unit, void *rootdata, const void *leafdata, void *leafupdate, MPI_Op op) 1611d71ae5a4SJacob Faibussowitsch { 1612eb02082bSJunchao Zhang PetscMemType rootmtype, leafmtype, leafupdatemtype; 1613a1729e3fSJunchao Zhang 1614a1729e3fSJunchao Zhang PetscFunctionBegin; 1615a1729e3fSJunchao Zhang PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 16169566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 16179566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_FetchAndOpBegin, sf, 0, 0, 0)); 16189566063dSJacob Faibussowitsch PetscCall(PetscGetMemType(rootdata, &rootmtype)); 16199566063dSJacob Faibussowitsch PetscCall(PetscGetMemType(leafdata, &leafmtype)); 16209566063dSJacob Faibussowitsch PetscCall(PetscGetMemType(leafupdate, &leafupdatemtype)); 162108401ef6SPierre Jolivet PetscCheck(leafmtype == leafupdatemtype, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for leafdata and leafupdate in different memory types"); 1622dbbe0bcdSBarry Smith PetscUseTypeMethod(sf, FetchAndOpBegin, unit, rootmtype, rootdata, leafmtype, leafdata, leafupdate, op); 16239566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_FetchAndOpBegin, sf, 0, 0, 0)); 16243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1625a1729e3fSJunchao Zhang } 1626a1729e3fSJunchao Zhang 1627a1729e3fSJunchao Zhang /*@C 1628cab54364SBarry Smith PetscSFFetchAndOpWithMemTypeBegin - begin operation with explicit memory types that fetches values from root and updates atomically by 1629cab54364SBarry Smith applying operation using my leaf value, to be completed with `PetscSFFetchAndOpEnd()` 1630d3b3e55cSJunchao Zhang 1631d3b3e55cSJunchao Zhang Collective 1632d3b3e55cSJunchao Zhang 1633d3b3e55cSJunchao Zhang Input Parameters: 1634d3b3e55cSJunchao Zhang + sf - star forest 1635d3b3e55cSJunchao Zhang . unit - data type 1636d3b3e55cSJunchao Zhang . rootmtype - memory type of rootdata 1637d3b3e55cSJunchao Zhang . leafmtype - memory type of leafdata 1638d3b3e55cSJunchao Zhang . leafdata - leaf values to use in reduction 1639d3b3e55cSJunchao Zhang . leafupdatemtype - memory type of leafupdate 1640d3b3e55cSJunchao Zhang - op - operation to use for reduction 1641d3b3e55cSJunchao Zhang 1642d3b3e55cSJunchao Zhang Output Parameters: 1643d3b3e55cSJunchao Zhang + rootdata - root values to be updated, input state is seen by first process to perform an update 1644d3b3e55cSJunchao Zhang - leafupdate - state at each leaf's respective root immediately prior to my atomic update 1645d3b3e55cSJunchao Zhang 1646d3b3e55cSJunchao Zhang Level: advanced 1647d3b3e55cSJunchao Zhang 1648cab54364SBarry Smith Note: 1649cab54364SBarry Smith See `PetscSFFetchAndOpBegin()` for more details. 1650d3b3e55cSJunchao Zhang 165120662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFFetchAndOpBegin()`, `PetscSFComputeDegreeBegin()`, `PetscSFReduceBegin()`, `PetscSFSetGraph()`, `PetscSFFetchAndOpEnd()` 1652d3b3e55cSJunchao Zhang @*/ 1653d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFFetchAndOpWithMemTypeBegin(PetscSF sf, MPI_Datatype unit, PetscMemType rootmtype, void *rootdata, PetscMemType leafmtype, const void *leafdata, PetscMemType leafupdatemtype, void *leafupdate, MPI_Op op) 1654d71ae5a4SJacob Faibussowitsch { 1655d3b3e55cSJunchao Zhang PetscFunctionBegin; 1656d3b3e55cSJunchao Zhang PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 16579566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 16589566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_FetchAndOpBegin, sf, 0, 0, 0)); 165908401ef6SPierre Jolivet PetscCheck(leafmtype == leafupdatemtype, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for leafdata and leafupdate in different memory types"); 1660dbbe0bcdSBarry Smith PetscUseTypeMethod(sf, FetchAndOpBegin, unit, rootmtype, rootdata, leafmtype, leafdata, leafupdate, op); 16619566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_FetchAndOpBegin, sf, 0, 0, 0)); 16623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1663d3b3e55cSJunchao Zhang } 1664d3b3e55cSJunchao Zhang 1665d3b3e55cSJunchao Zhang /*@C 166620662ed9SBarry Smith PetscSFFetchAndOpEnd - end operation started in matching call to `PetscSFFetchAndOpBegin()` or `PetscSFFetchAndOpWithMemTypeBegin()` 166720662ed9SBarry Smith to fetch values from roots and update atomically by applying operation using my leaf value 1668a1729e3fSJunchao Zhang 1669a1729e3fSJunchao Zhang Collective 1670a1729e3fSJunchao Zhang 16714165533cSJose E. Roman Input Parameters: 1672a1729e3fSJunchao Zhang + sf - star forest 1673a1729e3fSJunchao Zhang . unit - data type 1674a1729e3fSJunchao Zhang . leafdata - leaf values to use in reduction 1675a1729e3fSJunchao Zhang - op - operation to use for reduction 1676a1729e3fSJunchao Zhang 16774165533cSJose E. Roman Output Parameters: 1678a1729e3fSJunchao Zhang + rootdata - root values to be updated, input state is seen by first process to perform an update 1679a1729e3fSJunchao Zhang - leafupdate - state at each leaf's respective root immediately prior to my atomic update 1680a1729e3fSJunchao Zhang 1681a1729e3fSJunchao Zhang Level: advanced 1682a1729e3fSJunchao Zhang 168320662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeEnd()`, `PetscSFReduceEnd()`, `PetscSFSetGraph()`, `PetscSFFetchAndOpBegin()`, `PetscSFFetchAndOpWithMemTypeBegin()` 1684a1729e3fSJunchao Zhang @*/ 1685d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFFetchAndOpEnd(PetscSF sf, MPI_Datatype unit, void *rootdata, const void *leafdata, void *leafupdate, MPI_Op op) 1686d71ae5a4SJacob Faibussowitsch { 1687a1729e3fSJunchao Zhang PetscFunctionBegin; 1688a1729e3fSJunchao Zhang PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 16899566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_FetchAndOpEnd, sf, 0, 0, 0)); 1690dbbe0bcdSBarry Smith PetscUseTypeMethod(sf, FetchAndOpEnd, unit, rootdata, leafdata, leafupdate, op); 16919566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_FetchAndOpEnd, sf, 0, 0, 0)); 16923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1693a1729e3fSJunchao Zhang } 1694a1729e3fSJunchao Zhang 1695a1729e3fSJunchao Zhang /*@C 1696cab54364SBarry Smith PetscSFComputeDegreeBegin - begin computation of degree for each root vertex, to be completed with `PetscSFComputeDegreeEnd()` 169795fce210SBarry Smith 169895fce210SBarry Smith Collective 169995fce210SBarry Smith 17004165533cSJose E. Roman Input Parameter: 170195fce210SBarry Smith . sf - star forest 170295fce210SBarry Smith 17034165533cSJose E. Roman Output Parameter: 170495fce210SBarry Smith . degree - degree of each root vertex 170595fce210SBarry Smith 170695fce210SBarry Smith Level: advanced 170795fce210SBarry Smith 1708cab54364SBarry Smith Note: 170920662ed9SBarry Smith The returned array is owned by `PetscSF` and automatically freed by `PetscSFDestroy()`. Hence there is no need to call `PetscFree()` on it. 1710ffe67aa5SVáclav Hapla 1711cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGatherBegin()`, `PetscSFComputeDegreeEnd()` 171295fce210SBarry Smith @*/ 1713d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFComputeDegreeBegin(PetscSF sf, const PetscInt **degree) 1714d71ae5a4SJacob Faibussowitsch { 171595fce210SBarry Smith PetscFunctionBegin; 171695fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 171795fce210SBarry Smith PetscSFCheckGraphSet(sf, 1); 171895fce210SBarry Smith PetscValidPointer(degree, 2); 1719803bd9e8SMatthew G. Knepley if (!sf->degreeknown) { 17205b0d146aSStefano Zampini PetscInt i, nroots = sf->nroots, maxlocal; 172128b400f6SJacob Faibussowitsch PetscCheck(!sf->degree, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Calls to PetscSFComputeDegreeBegin() cannot be nested."); 17225b0d146aSStefano Zampini maxlocal = sf->maxleaf - sf->minleaf + 1; 17239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nroots, &sf->degree)); 17249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PetscMax(maxlocal, 1), &sf->degreetmp)); /* allocate at least one entry, see check in PetscSFComputeDegreeEnd() */ 172529046d53SLisandro Dalcin for (i = 0; i < nroots; i++) sf->degree[i] = 0; 17269837ea96SMatthew G. Knepley for (i = 0; i < maxlocal; i++) sf->degreetmp[i] = 1; 17279566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf, MPIU_INT, sf->degreetmp - sf->minleaf, sf->degree, MPI_SUM)); 172895fce210SBarry Smith } 172995fce210SBarry Smith *degree = NULL; 17303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173195fce210SBarry Smith } 173295fce210SBarry Smith 173395fce210SBarry Smith /*@C 1734cab54364SBarry Smith PetscSFComputeDegreeEnd - complete computation of degree for each root vertex, started with `PetscSFComputeDegreeBegin()` 173595fce210SBarry Smith 173695fce210SBarry Smith Collective 173795fce210SBarry Smith 17384165533cSJose E. Roman Input Parameter: 173995fce210SBarry Smith . sf - star forest 174095fce210SBarry Smith 17414165533cSJose E. Roman Output Parameter: 174295fce210SBarry Smith . degree - degree of each root vertex 174395fce210SBarry Smith 174495fce210SBarry Smith Level: developer 174595fce210SBarry Smith 1746cab54364SBarry Smith Note: 174720662ed9SBarry Smith The returned array is owned by `PetscSF` and automatically freed by `PetscSFDestroy()`. Hence there is no need to call `PetscFree()` on it. 1748ffe67aa5SVáclav Hapla 1749cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFGatherBegin()`, `PetscSFComputeDegreeBegin()` 175095fce210SBarry Smith @*/ 1751d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFComputeDegreeEnd(PetscSF sf, const PetscInt **degree) 1752d71ae5a4SJacob Faibussowitsch { 175395fce210SBarry Smith PetscFunctionBegin; 175495fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 175595fce210SBarry Smith PetscSFCheckGraphSet(sf, 1); 175629046d53SLisandro Dalcin PetscValidPointer(degree, 2); 175795fce210SBarry Smith if (!sf->degreeknown) { 175828b400f6SJacob Faibussowitsch PetscCheck(sf->degreetmp, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PetscSFComputeDegreeBegin() before PetscSFComputeDegreeEnd()"); 17599566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sf, MPIU_INT, sf->degreetmp - sf->minleaf, sf->degree, MPI_SUM)); 17609566063dSJacob Faibussowitsch PetscCall(PetscFree(sf->degreetmp)); 176195fce210SBarry Smith sf->degreeknown = PETSC_TRUE; 176295fce210SBarry Smith } 176395fce210SBarry Smith *degree = sf->degree; 17643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 176595fce210SBarry Smith } 176695fce210SBarry Smith 1767673100f5SVaclav Hapla /*@C 176820662ed9SBarry Smith PetscSFComputeMultiRootOriginalNumbering - Returns original numbering of multi-roots (roots of multi-`PetscSF` returned by `PetscSFGetMultiSF()`). 176966dfcd1aSVaclav Hapla Each multi-root is assigned index of the corresponding original root. 1770673100f5SVaclav Hapla 1771673100f5SVaclav Hapla Collective 1772673100f5SVaclav Hapla 17734165533cSJose E. Roman Input Parameters: 1774673100f5SVaclav Hapla + sf - star forest 1775cab54364SBarry Smith - degree - degree of each root vertex, computed with `PetscSFComputeDegreeBegin()`/`PetscSFComputeDegreeEnd()` 1776673100f5SVaclav Hapla 17774165533cSJose E. Roman Output Parameters: 177820662ed9SBarry Smith + nMultiRoots - (optional) number of multi-roots (roots of multi-`PetscSF`) 177920662ed9SBarry Smith - multiRootsOrigNumbering - original indices of multi-roots; length of this array is `nMultiRoots` 1780673100f5SVaclav Hapla 1781673100f5SVaclav Hapla Level: developer 1782673100f5SVaclav Hapla 1783cab54364SBarry Smith Note: 178420662ed9SBarry Smith The returned array `multiRootsOrigNumbering` is newly allocated and should be destroyed with `PetscFree()` when no longer needed. 1785ffe67aa5SVáclav Hapla 1786cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeBegin()`, `PetscSFComputeDegreeEnd()`, `PetscSFGetMultiSF()` 1787673100f5SVaclav Hapla @*/ 1788d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFComputeMultiRootOriginalNumbering(PetscSF sf, const PetscInt degree[], PetscInt *nMultiRoots, PetscInt *multiRootsOrigNumbering[]) 1789d71ae5a4SJacob Faibussowitsch { 1790673100f5SVaclav Hapla PetscSF msf; 1791673100f5SVaclav Hapla PetscInt i, j, k, nroots, nmroots; 1792673100f5SVaclav Hapla 1793673100f5SVaclav Hapla PetscFunctionBegin; 1794673100f5SVaclav Hapla PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 17959566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 179629328920SVaclav Hapla if (nroots) PetscValidIntPointer(degree, 2); 179766dfcd1aSVaclav Hapla if (nMultiRoots) PetscValidIntPointer(nMultiRoots, 3); 179866dfcd1aSVaclav Hapla PetscValidPointer(multiRootsOrigNumbering, 4); 17999566063dSJacob Faibussowitsch PetscCall(PetscSFGetMultiSF(sf, &msf)); 18009566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(msf, &nmroots, NULL, NULL, NULL)); 18019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nmroots, multiRootsOrigNumbering)); 1802673100f5SVaclav Hapla for (i = 0, j = 0, k = 0; i < nroots; i++) { 1803673100f5SVaclav Hapla if (!degree[i]) continue; 1804ad540459SPierre Jolivet for (j = 0; j < degree[i]; j++, k++) (*multiRootsOrigNumbering)[k] = i; 1805673100f5SVaclav Hapla } 180608401ef6SPierre Jolivet PetscCheck(k == nmroots, PETSC_COMM_SELF, PETSC_ERR_PLIB, "sanity check fail"); 180766dfcd1aSVaclav Hapla if (nMultiRoots) *nMultiRoots = nmroots; 18083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1809673100f5SVaclav Hapla } 1810673100f5SVaclav Hapla 181195fce210SBarry Smith /*@C 1812cab54364SBarry Smith PetscSFGatherBegin - begin pointwise gather of all leaves into multi-roots, to be completed with `PetscSFGatherEnd()` 181395fce210SBarry Smith 181495fce210SBarry Smith Collective 181595fce210SBarry Smith 18164165533cSJose E. Roman Input Parameters: 181795fce210SBarry Smith + sf - star forest 181895fce210SBarry Smith . unit - data type 181995fce210SBarry Smith - leafdata - leaf data to gather to roots 182095fce210SBarry Smith 18214165533cSJose E. Roman Output Parameter: 182295fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 182395fce210SBarry Smith 182495fce210SBarry Smith Level: intermediate 182595fce210SBarry Smith 1826cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeBegin()`, `PetscSFScatterBegin()` 182795fce210SBarry Smith @*/ 1828d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGatherBegin(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *multirootdata) 1829d71ae5a4SJacob Faibussowitsch { 1830a5526d50SJunchao Zhang PetscSF multi = NULL; 183195fce210SBarry Smith 183295fce210SBarry Smith PetscFunctionBegin; 183395fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 18349566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 18359566063dSJacob Faibussowitsch PetscCall(PetscSFGetMultiSF(sf, &multi)); 18369566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(multi, unit, leafdata, multirootdata, MPI_REPLACE)); 18373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 183895fce210SBarry Smith } 183995fce210SBarry Smith 184095fce210SBarry Smith /*@C 1841cab54364SBarry Smith PetscSFGatherEnd - ends pointwise gather operation that was started with `PetscSFGatherBegin()` 184295fce210SBarry Smith 184395fce210SBarry Smith Collective 184495fce210SBarry Smith 18454165533cSJose E. Roman Input Parameters: 184695fce210SBarry Smith + sf - star forest 184795fce210SBarry Smith . unit - data type 184895fce210SBarry Smith - leafdata - leaf data to gather to roots 184995fce210SBarry Smith 18504165533cSJose E. Roman Output Parameter: 185195fce210SBarry Smith . multirootdata - root buffer to gather into, amount of space per root is equal to its degree 185295fce210SBarry Smith 185395fce210SBarry Smith Level: intermediate 185495fce210SBarry Smith 1855cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeEnd()`, `PetscSFScatterEnd()` 185695fce210SBarry Smith @*/ 1857d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGatherEnd(PetscSF sf, MPI_Datatype unit, const void *leafdata, void *multirootdata) 1858d71ae5a4SJacob Faibussowitsch { 1859a5526d50SJunchao Zhang PetscSF multi = NULL; 186095fce210SBarry Smith 186195fce210SBarry Smith PetscFunctionBegin; 186295fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 18639566063dSJacob Faibussowitsch PetscCall(PetscSFGetMultiSF(sf, &multi)); 18649566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(multi, unit, leafdata, multirootdata, MPI_REPLACE)); 18653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 186695fce210SBarry Smith } 186795fce210SBarry Smith 186895fce210SBarry Smith /*@C 1869cab54364SBarry Smith PetscSFScatterBegin - begin pointwise scatter operation from multi-roots to leaves, to be completed with `PetscSFScatterEnd()` 187095fce210SBarry Smith 187195fce210SBarry Smith Collective 187295fce210SBarry Smith 18734165533cSJose E. Roman Input Parameters: 187495fce210SBarry Smith + sf - star forest 187595fce210SBarry Smith . unit - data type 187695fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 187795fce210SBarry Smith 18784165533cSJose E. Roman Output Parameter: 187995fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 188095fce210SBarry Smith 188195fce210SBarry Smith Level: intermediate 188295fce210SBarry Smith 188320662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeBegin()`, `PetscSFScatterEnd()` 188495fce210SBarry Smith @*/ 1885d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFScatterBegin(PetscSF sf, MPI_Datatype unit, const void *multirootdata, void *leafdata) 1886d71ae5a4SJacob Faibussowitsch { 1887a5526d50SJunchao Zhang PetscSF multi = NULL; 188895fce210SBarry Smith 188995fce210SBarry Smith PetscFunctionBegin; 189095fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 18919566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 18929566063dSJacob Faibussowitsch PetscCall(PetscSFGetMultiSF(sf, &multi)); 18939566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(multi, unit, multirootdata, leafdata, MPI_REPLACE)); 18943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 189595fce210SBarry Smith } 189695fce210SBarry Smith 189795fce210SBarry Smith /*@C 1898cab54364SBarry Smith PetscSFScatterEnd - ends pointwise scatter operation that was started with `PetscSFScatterBegin()` 189995fce210SBarry Smith 190095fce210SBarry Smith Collective 190195fce210SBarry Smith 19024165533cSJose E. Roman Input Parameters: 190395fce210SBarry Smith + sf - star forest 190495fce210SBarry Smith . unit - data type 190595fce210SBarry Smith - multirootdata - root buffer to send to each leaf, one unit of data per leaf 190695fce210SBarry Smith 19074165533cSJose E. Roman Output Parameter: 190895fce210SBarry Smith . leafdata - leaf data to be update with personal data from each respective root 190995fce210SBarry Smith 191095fce210SBarry Smith Level: intermediate 191195fce210SBarry Smith 191220662ed9SBarry Smith .seealso: `PetscSF`, `PetscSFComputeDegreeEnd()`, `PetscSFScatterBegin()` 191395fce210SBarry Smith @*/ 1914d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFScatterEnd(PetscSF sf, MPI_Datatype unit, const void *multirootdata, void *leafdata) 1915d71ae5a4SJacob Faibussowitsch { 1916a5526d50SJunchao Zhang PetscSF multi = NULL; 191795fce210SBarry Smith 191895fce210SBarry Smith PetscFunctionBegin; 191995fce210SBarry Smith PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 19209566063dSJacob Faibussowitsch PetscCall(PetscSFGetMultiSF(sf, &multi)); 19219566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(multi, unit, multirootdata, leafdata, MPI_REPLACE)); 19223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 192395fce210SBarry Smith } 1924a7b3aa13SAta Mesgarnejad 1925d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFCheckLeavesUnique_Private(PetscSF sf) 1926d71ae5a4SJacob Faibussowitsch { 1927a072220fSLawrence Mitchell PetscInt i, n, nleaves; 1928a072220fSLawrence Mitchell const PetscInt *ilocal = NULL; 1929a072220fSLawrence Mitchell PetscHSetI seen; 1930a072220fSLawrence Mitchell 1931a072220fSLawrence Mitchell PetscFunctionBegin; 1932b458e8f1SJose E. Roman if (PetscDefined(USE_DEBUG)) { 19339566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, &ilocal, NULL)); 19349566063dSJacob Faibussowitsch PetscCall(PetscHSetICreate(&seen)); 1935a072220fSLawrence Mitchell for (i = 0; i < nleaves; i++) { 1936a072220fSLawrence Mitchell const PetscInt leaf = ilocal ? ilocal[i] : i; 19379566063dSJacob Faibussowitsch PetscCall(PetscHSetIAdd(seen, leaf)); 1938a072220fSLawrence Mitchell } 19399566063dSJacob Faibussowitsch PetscCall(PetscHSetIGetSize(seen, &n)); 194008401ef6SPierre Jolivet PetscCheck(n == nleaves, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Provided leaves have repeated values: all leaves must be unique"); 19419566063dSJacob Faibussowitsch PetscCall(PetscHSetIDestroy(&seen)); 1942b458e8f1SJose E. Roman } 19433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1944a072220fSLawrence Mitchell } 194554729392SStefano Zampini 1946a7b3aa13SAta Mesgarnejad /*@ 1947cab54364SBarry Smith PetscSFCompose - Compose a new `PetscSF` by putting the second `PetscSF` under the first one in a top (roots) down (leaves) view 1948a7b3aa13SAta Mesgarnejad 1949a7b3aa13SAta Mesgarnejad Input Parameters: 1950cab54364SBarry Smith + sfA - The first `PetscSF` 1951cab54364SBarry Smith - sfB - The second `PetscSF` 1952a7b3aa13SAta Mesgarnejad 1953a7b3aa13SAta Mesgarnejad Output Parameters: 1954cab54364SBarry Smith . sfBA - The composite `PetscSF` 1955a7b3aa13SAta Mesgarnejad 1956a7b3aa13SAta Mesgarnejad Level: developer 1957a7b3aa13SAta Mesgarnejad 1958a072220fSLawrence Mitchell Notes: 1959cab54364SBarry Smith Currently, the two `PetscSF`s must be defined on congruent communicators and they must be true star 196054729392SStefano Zampini forests, i.e. the same leaf is not connected with different roots. 196154729392SStefano Zampini 196220662ed9SBarry Smith `sfA`'s leaf space and `sfB`'s root space might be partially overlapped. The composition builds 196320662ed9SBarry Smith a graph with `sfA`'s roots and `sfB`'s leaves only when there is a path between them. Unconnected 196420662ed9SBarry Smith nodes (roots or leaves) are not in `sfBA`. Doing a `PetscSFBcastBegin()`/`PetscSFBcastEnd()` on the new `PetscSF` is equivalent to doing a 196520662ed9SBarry Smith `PetscSFBcastBegin()`/`PetscSFBcastEnd()` on `sfA`, then a `PetscSFBcastBegin()`/`PetscSFBcastEnd()` on `sfB`, on connected nodes. 1966a072220fSLawrence Mitchell 1967db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFComposeInverse()`, `PetscSFGetGraph()`, `PetscSFSetGraph()` 1968a7b3aa13SAta Mesgarnejad @*/ 1969d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCompose(PetscSF sfA, PetscSF sfB, PetscSF *sfBA) 1970d71ae5a4SJacob Faibussowitsch { 1971a7b3aa13SAta Mesgarnejad const PetscSFNode *remotePointsA, *remotePointsB; 1972d41018fbSJunchao Zhang PetscSFNode *remotePointsBA = NULL, *reorderedRemotePointsA = NULL, *leafdataB; 197354729392SStefano Zampini const PetscInt *localPointsA, *localPointsB; 197454729392SStefano Zampini PetscInt *localPointsBA; 197554729392SStefano Zampini PetscInt i, numRootsA, numLeavesA, numRootsB, numLeavesB, minleaf, maxleaf, numLeavesBA; 197654729392SStefano Zampini PetscBool denseB; 1977a7b3aa13SAta Mesgarnejad 1978a7b3aa13SAta Mesgarnejad PetscFunctionBegin; 1979a7b3aa13SAta Mesgarnejad PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1); 198029046d53SLisandro Dalcin PetscSFCheckGraphSet(sfA, 1); 198129046d53SLisandro Dalcin PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 2); 198229046d53SLisandro Dalcin PetscSFCheckGraphSet(sfB, 2); 198354729392SStefano Zampini PetscCheckSameComm(sfA, 1, sfB, 2); 198429046d53SLisandro Dalcin PetscValidPointer(sfBA, 3); 19859566063dSJacob Faibussowitsch PetscCall(PetscSFCheckLeavesUnique_Private(sfA)); 19869566063dSJacob Faibussowitsch PetscCall(PetscSFCheckLeavesUnique_Private(sfB)); 198754729392SStefano Zampini 19889566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA)); 19899566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB)); 199020662ed9SBarry Smith /* Make sure that PetscSFBcast{Begin, End}(sfB, ...) works with root data of size 199120662ed9SBarry Smith numRootsB; otherwise, garbage will be broadcasted. 199220662ed9SBarry Smith Example (comm size = 1): 199320662ed9SBarry Smith sfA: 0 <- (0, 0) 199420662ed9SBarry Smith sfB: 100 <- (0, 0) 199520662ed9SBarry Smith 101 <- (0, 1) 199620662ed9SBarry Smith Here, we have remotePointsA = [(0, 0)], but for remotePointsA to be a valid tartget 199720662ed9SBarry Smith of sfB, it has to be recasted as [(0, 0), (-1, -1)] so that points 100 and 101 would 199820662ed9SBarry Smith receive (0, 0) and (-1, -1), respectively, when PetscSFBcast(sfB, ...) is called on 199920662ed9SBarry Smith remotePointsA; if not recasted, point 101 would receive a garbage value. */ 20009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numRootsB, &reorderedRemotePointsA)); 200154729392SStefano Zampini for (i = 0; i < numRootsB; i++) { 200254729392SStefano Zampini reorderedRemotePointsA[i].rank = -1; 200354729392SStefano Zampini reorderedRemotePointsA[i].index = -1; 200454729392SStefano Zampini } 200554729392SStefano Zampini for (i = 0; i < numLeavesA; i++) { 20060ea77edaSksagiyam PetscInt localp = localPointsA ? localPointsA[i] : i; 20070ea77edaSksagiyam 20080ea77edaSksagiyam if (localp >= numRootsB) continue; 20090ea77edaSksagiyam reorderedRemotePointsA[localp] = remotePointsA[i]; 201054729392SStefano Zampini } 2011d41018fbSJunchao Zhang remotePointsA = reorderedRemotePointsA; 20129566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(sfB, &minleaf, &maxleaf)); 20139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxleaf - minleaf + 1, &leafdataB)); 20140ea77edaSksagiyam for (i = 0; i < maxleaf - minleaf + 1; i++) { 20150ea77edaSksagiyam leafdataB[i].rank = -1; 20160ea77edaSksagiyam leafdataB[i].index = -1; 20170ea77edaSksagiyam } 20189566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sfB, MPIU_2INT, remotePointsA, leafdataB - minleaf, MPI_REPLACE)); 20199566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sfB, MPIU_2INT, remotePointsA, leafdataB - minleaf, MPI_REPLACE)); 20209566063dSJacob Faibussowitsch PetscCall(PetscFree(reorderedRemotePointsA)); 2021d41018fbSJunchao Zhang 202254729392SStefano Zampini denseB = (PetscBool)!localPointsB; 202354729392SStefano Zampini for (i = 0, numLeavesBA = 0; i < numLeavesB; i++) { 202454729392SStefano Zampini if (leafdataB[localPointsB ? localPointsB[i] - minleaf : i].rank == -1) denseB = PETSC_FALSE; 202554729392SStefano Zampini else numLeavesBA++; 202654729392SStefano Zampini } 202754729392SStefano Zampini if (denseB) { 2028d41018fbSJunchao Zhang localPointsBA = NULL; 2029d41018fbSJunchao Zhang remotePointsBA = leafdataB; 2030d41018fbSJunchao Zhang } else { 20319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numLeavesBA, &localPointsBA)); 20329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numLeavesBA, &remotePointsBA)); 203354729392SStefano Zampini for (i = 0, numLeavesBA = 0; i < numLeavesB; i++) { 203454729392SStefano Zampini const PetscInt l = localPointsB ? localPointsB[i] : i; 203554729392SStefano Zampini 203654729392SStefano Zampini if (leafdataB[l - minleaf].rank == -1) continue; 203754729392SStefano Zampini remotePointsBA[numLeavesBA] = leafdataB[l - minleaf]; 203854729392SStefano Zampini localPointsBA[numLeavesBA] = l; 203954729392SStefano Zampini numLeavesBA++; 204054729392SStefano Zampini } 20419566063dSJacob Faibussowitsch PetscCall(PetscFree(leafdataB)); 2042d41018fbSJunchao Zhang } 20439566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfA), sfBA)); 20449566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(*sfBA)); 20459566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sfBA, numRootsA, numLeavesBA, localPointsBA, PETSC_OWN_POINTER, remotePointsBA, PETSC_OWN_POINTER)); 20463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2047a7b3aa13SAta Mesgarnejad } 20481c6ba672SJunchao Zhang 204904c0ada0SJunchao Zhang /*@ 2050cab54364SBarry Smith PetscSFComposeInverse - Compose a new `PetscSF` by putting the inverse of the second `PetscSF` under the first one 205104c0ada0SJunchao Zhang 205204c0ada0SJunchao Zhang Input Parameters: 2053cab54364SBarry Smith + sfA - The first `PetscSF` 2054cab54364SBarry Smith - sfB - The second `PetscSF` 205504c0ada0SJunchao Zhang 205604c0ada0SJunchao Zhang Output Parameters: 2057cab54364SBarry Smith . sfBA - The composite `PetscSF`. 205804c0ada0SJunchao Zhang 205904c0ada0SJunchao Zhang Level: developer 206004c0ada0SJunchao Zhang 206154729392SStefano Zampini Notes: 206220662ed9SBarry Smith Currently, the two `PetscSF`s must be defined on congruent communicators and they must be true star 206354729392SStefano Zampini forests, i.e. the same leaf is not connected with different roots. Even more, all roots of the 206420662ed9SBarry Smith second `PetscSF` must have a degree of 1, i.e., no roots have more than one leaf connected. 206554729392SStefano Zampini 206620662ed9SBarry Smith `sfA`'s leaf space and `sfB`'s leaf space might be partially overlapped. The composition builds 206720662ed9SBarry Smith a graph with `sfA`'s roots and `sfB`'s roots only when there is a path between them. Unconnected 206820662ed9SBarry Smith roots are not in `sfBA`. Doing a `PetscSFBcastBegin()`/`PetscSFBcastEnd()` on the new `PetscSF` is equivalent to doing a `PetscSFBcastBegin()`/`PetscSFBcastEnd()` 206920662ed9SBarry Smith on `sfA`, then 207020662ed9SBarry Smith a `PetscSFReduceBegin()`/`PetscSFReduceEnd()` on `sfB`, on connected roots. 207154729392SStefano Zampini 2072db781477SPatrick Sanan .seealso: `PetscSF`, `PetscSFCompose()`, `PetscSFGetGraph()`, `PetscSFSetGraph()`, `PetscSFCreateInverseSF()` 207304c0ada0SJunchao Zhang @*/ 2074d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFComposeInverse(PetscSF sfA, PetscSF sfB, PetscSF *sfBA) 2075d71ae5a4SJacob Faibussowitsch { 207604c0ada0SJunchao Zhang const PetscSFNode *remotePointsA, *remotePointsB; 207704c0ada0SJunchao Zhang PetscSFNode *remotePointsBA; 207804c0ada0SJunchao Zhang const PetscInt *localPointsA, *localPointsB; 207954729392SStefano Zampini PetscSFNode *reorderedRemotePointsA = NULL; 208054729392SStefano Zampini PetscInt i, numRootsA, numLeavesA, numLeavesBA, numRootsB, numLeavesB, minleaf, maxleaf, *localPointsBA; 20815b0d146aSStefano Zampini MPI_Op op; 20825b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES) 20835b0d146aSStefano Zampini PetscBool iswin; 20845b0d146aSStefano Zampini #endif 208504c0ada0SJunchao Zhang 208604c0ada0SJunchao Zhang PetscFunctionBegin; 208704c0ada0SJunchao Zhang PetscValidHeaderSpecific(sfA, PETSCSF_CLASSID, 1); 208804c0ada0SJunchao Zhang PetscSFCheckGraphSet(sfA, 1); 208904c0ada0SJunchao Zhang PetscValidHeaderSpecific(sfB, PETSCSF_CLASSID, 2); 209004c0ada0SJunchao Zhang PetscSFCheckGraphSet(sfB, 2); 209154729392SStefano Zampini PetscCheckSameComm(sfA, 1, sfB, 2); 209204c0ada0SJunchao Zhang PetscValidPointer(sfBA, 3); 20939566063dSJacob Faibussowitsch PetscCall(PetscSFCheckLeavesUnique_Private(sfA)); 20949566063dSJacob Faibussowitsch PetscCall(PetscSFCheckLeavesUnique_Private(sfB)); 209554729392SStefano Zampini 20969566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfA, &numRootsA, &numLeavesA, &localPointsA, &remotePointsA)); 20979566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfB, &numRootsB, &numLeavesB, &localPointsB, &remotePointsB)); 20985b0d146aSStefano Zampini 20995b0d146aSStefano Zampini /* TODO: Check roots of sfB have degree of 1 */ 21005b0d146aSStefano Zampini /* Once we implement it, we can replace the MPI_MAXLOC 210183df288dSJunchao Zhang with MPI_REPLACE. In that case, MPI_MAXLOC and MPI_REPLACE have the same effect. 21025b0d146aSStefano Zampini We use MPI_MAXLOC only to have a deterministic output from this routine if 21035b0d146aSStefano Zampini the root condition is not meet. 21045b0d146aSStefano Zampini */ 21055b0d146aSStefano Zampini op = MPI_MAXLOC; 21065b0d146aSStefano Zampini #if defined(PETSC_USE_64BIT_INDICES) 21075b0d146aSStefano Zampini /* we accept a non-deterministic output (if any) with PETSCSFWINDOW, since MPI_MAXLOC cannot operate on MPIU_2INT with MPI_Accumulate */ 21089566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)sfB, PETSCSFWINDOW, &iswin)); 210983df288dSJunchao Zhang if (iswin) op = MPI_REPLACE; 21105b0d146aSStefano Zampini #endif 21115b0d146aSStefano Zampini 21129566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(sfB, &minleaf, &maxleaf)); 21139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxleaf - minleaf + 1, &reorderedRemotePointsA)); 211454729392SStefano Zampini for (i = 0; i < maxleaf - minleaf + 1; i++) { 211554729392SStefano Zampini reorderedRemotePointsA[i].rank = -1; 211654729392SStefano Zampini reorderedRemotePointsA[i].index = -1; 211754729392SStefano Zampini } 211854729392SStefano Zampini if (localPointsA) { 211954729392SStefano Zampini for (i = 0; i < numLeavesA; i++) { 212054729392SStefano Zampini if (localPointsA[i] > maxleaf || localPointsA[i] < minleaf) continue; 212154729392SStefano Zampini reorderedRemotePointsA[localPointsA[i] - minleaf] = remotePointsA[i]; 212254729392SStefano Zampini } 212354729392SStefano Zampini } else { 212454729392SStefano Zampini for (i = 0; i < numLeavesA; i++) { 212554729392SStefano Zampini if (i > maxleaf || i < minleaf) continue; 212654729392SStefano Zampini reorderedRemotePointsA[i - minleaf] = remotePointsA[i]; 212754729392SStefano Zampini } 212854729392SStefano Zampini } 212954729392SStefano Zampini 21309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numRootsB, &localPointsBA)); 21319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numRootsB, &remotePointsBA)); 213254729392SStefano Zampini for (i = 0; i < numRootsB; i++) { 213354729392SStefano Zampini remotePointsBA[i].rank = -1; 213454729392SStefano Zampini remotePointsBA[i].index = -1; 213554729392SStefano Zampini } 213654729392SStefano Zampini 21379566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sfB, MPIU_2INT, reorderedRemotePointsA - minleaf, remotePointsBA, op)); 21389566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sfB, MPIU_2INT, reorderedRemotePointsA - minleaf, remotePointsBA, op)); 21399566063dSJacob Faibussowitsch PetscCall(PetscFree(reorderedRemotePointsA)); 214054729392SStefano Zampini for (i = 0, numLeavesBA = 0; i < numRootsB; i++) { 214154729392SStefano Zampini if (remotePointsBA[i].rank == -1) continue; 214254729392SStefano Zampini remotePointsBA[numLeavesBA].rank = remotePointsBA[i].rank; 214354729392SStefano Zampini remotePointsBA[numLeavesBA].index = remotePointsBA[i].index; 214454729392SStefano Zampini localPointsBA[numLeavesBA] = i; 214554729392SStefano Zampini numLeavesBA++; 214654729392SStefano Zampini } 21479566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)sfA), sfBA)); 21489566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(*sfBA)); 21499566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sfBA, numRootsA, numLeavesBA, localPointsBA, PETSC_OWN_POINTER, remotePointsBA, PETSC_OWN_POINTER)); 21503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 215104c0ada0SJunchao Zhang } 215204c0ada0SJunchao Zhang 21531c6ba672SJunchao Zhang /* 2154cab54364SBarry Smith PetscSFCreateLocalSF_Private - Creates a local `PetscSF` that only has intra-process edges of the global `PetscSF` 21551c6ba672SJunchao Zhang 21561c6ba672SJunchao Zhang Input Parameters: 2157cab54364SBarry Smith . sf - The global `PetscSF` 21581c6ba672SJunchao Zhang 21591c6ba672SJunchao Zhang Output Parameters: 2160cab54364SBarry Smith . out - The local `PetscSF` 2161cab54364SBarry Smith 2162cab54364SBarry Smith .seealso: `PetscSF`, `PetscSFCreate()` 21631c6ba672SJunchao Zhang */ 2164d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFCreateLocalSF_Private(PetscSF sf, PetscSF *out) 2165d71ae5a4SJacob Faibussowitsch { 21661c6ba672SJunchao Zhang MPI_Comm comm; 21671c6ba672SJunchao Zhang PetscMPIInt myrank; 21681c6ba672SJunchao Zhang const PetscInt *ilocal; 21691c6ba672SJunchao Zhang const PetscSFNode *iremote; 21701c6ba672SJunchao Zhang PetscInt i, j, nroots, nleaves, lnleaves, *lilocal; 21711c6ba672SJunchao Zhang PetscSFNode *liremote; 21721c6ba672SJunchao Zhang PetscSF lsf; 21731c6ba672SJunchao Zhang 21741c6ba672SJunchao Zhang PetscFunctionBegin; 21751c6ba672SJunchao Zhang PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 2176dbbe0bcdSBarry Smith if (sf->ops->CreateLocalSF) PetscUseTypeMethod(sf, CreateLocalSF, out); 2177dbbe0bcdSBarry Smith else { 21781c6ba672SJunchao Zhang /* Could use PetscSFCreateEmbeddedLeafSF, but since we know the comm is PETSC_COMM_SELF, we can make it fast */ 21799566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 21809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &myrank)); 21811c6ba672SJunchao Zhang 21821c6ba672SJunchao Zhang /* Find out local edges and build a local SF */ 21839566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, &iremote)); 21849371c9d4SSatish Balay for (i = lnleaves = 0; i < nleaves; i++) { 21859371c9d4SSatish Balay if (iremote[i].rank == (PetscInt)myrank) lnleaves++; 21869371c9d4SSatish Balay } 21879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lnleaves, &lilocal)); 21889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lnleaves, &liremote)); 21891c6ba672SJunchao Zhang 21901c6ba672SJunchao Zhang for (i = j = 0; i < nleaves; i++) { 21911c6ba672SJunchao Zhang if (iremote[i].rank == (PetscInt)myrank) { 21921c6ba672SJunchao Zhang lilocal[j] = ilocal ? ilocal[i] : i; /* ilocal=NULL for contiguous storage */ 21931c6ba672SJunchao Zhang liremote[j].rank = 0; /* rank in PETSC_COMM_SELF */ 21941c6ba672SJunchao Zhang liremote[j].index = iremote[i].index; 21951c6ba672SJunchao Zhang j++; 21961c6ba672SJunchao Zhang } 21971c6ba672SJunchao Zhang } 21989566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PETSC_COMM_SELF, &lsf)); 21999566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(lsf)); 22009566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(lsf, nroots, lnleaves, lilocal, PETSC_OWN_POINTER, liremote, PETSC_OWN_POINTER)); 22019566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(lsf)); 22021c6ba672SJunchao Zhang *out = lsf; 22031c6ba672SJunchao Zhang } 22043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22051c6ba672SJunchao Zhang } 2206dd5b3ca6SJunchao Zhang 2207dd5b3ca6SJunchao Zhang /* Similar to PetscSFBcast, but only Bcast to leaves on rank 0 */ 2208d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFBcastToZero_Private(PetscSF sf, MPI_Datatype unit, const void *rootdata, void *leafdata) 2209d71ae5a4SJacob Faibussowitsch { 2210eb02082bSJunchao Zhang PetscMemType rootmtype, leafmtype; 2211dd5b3ca6SJunchao Zhang 2212dd5b3ca6SJunchao Zhang PetscFunctionBegin; 2213dd5b3ca6SJunchao Zhang PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 22149566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sf)); 22159566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PETSCSF_BcastBegin, sf, 0, 0, 0)); 22169566063dSJacob Faibussowitsch PetscCall(PetscGetMemType(rootdata, &rootmtype)); 22179566063dSJacob Faibussowitsch PetscCall(PetscGetMemType(leafdata, &leafmtype)); 2218dbbe0bcdSBarry Smith PetscUseTypeMethod(sf, BcastToZero, unit, rootmtype, rootdata, leafmtype, leafdata); 22199566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PETSCSF_BcastBegin, sf, 0, 0, 0)); 22203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2221dd5b3ca6SJunchao Zhang } 2222dd5b3ca6SJunchao Zhang 2223157edd7aSVaclav Hapla /*@ 2224cab54364SBarry Smith PetscSFConcatenate - concatenate multiple `PetscSF` into one 2225157edd7aSVaclav Hapla 2226157edd7aSVaclav Hapla Input Parameters: 2227157edd7aSVaclav Hapla + comm - the communicator 2228cab54364SBarry Smith . nsfs - the number of input `PetscSF` 2229cab54364SBarry Smith . sfs - the array of input `PetscSF` 22301f40158dSVaclav Hapla . rootMode - the root mode specifying how roots are handled 223120662ed9SBarry Smith - leafOffsets - the array of local leaf offsets, one for each input `PetscSF`, or `NULL` for contiguous storage 2232157edd7aSVaclav Hapla 2233157edd7aSVaclav Hapla Output Parameters: 2234cab54364SBarry Smith . newsf - The resulting `PetscSF` 2235157edd7aSVaclav Hapla 22361f40158dSVaclav Hapla Level: advanced 2237157edd7aSVaclav Hapla 2238157edd7aSVaclav Hapla Notes: 223920662ed9SBarry Smith The communicator of all `PetscSF`s in `sfs` must be comm. 2240157edd7aSVaclav Hapla 224120662ed9SBarry Smith Leaves are always concatenated locally, keeping them ordered by the input `PetscSF` index and original local order. 224220662ed9SBarry Smith 224320662ed9SBarry Smith The offsets in `leafOffsets` are added to the original leaf indices. 224420662ed9SBarry Smith 224520662ed9SBarry Smith If all input SFs use contiguous leaf storage (`ilocal` = `NULL`), `leafOffsets` can be passed as `NULL` as well. 224620662ed9SBarry Smith In this case, `NULL` is also passed as `ilocal` to the resulting `PetscSF`. 224720662ed9SBarry Smith 224820662ed9SBarry Smith If any input `PetscSF` has non-null `ilocal`, `leafOffsets` is needed to distinguish leaves from different input `PetscSF`s. 2249157edd7aSVaclav Hapla In this case, user is responsible to provide correct offsets so that the resulting leaves are unique (otherwise an error occurs). 2250157edd7aSVaclav Hapla 225120662ed9SBarry Smith All root modes retain the essential connectivity condition. 225220662ed9SBarry Smith If two leaves of the same input `PetscSF` are connected (sharing the same root), they are also connected in the output `PetscSF`. 225320662ed9SBarry Smith Parameter `rootMode` controls how the input root spaces are combined. 225420662ed9SBarry Smith For `PETSCSF_CONCATENATE_ROOTMODE_SHARED`, the root space is considered the same for each input `PetscSF` (checked in debug mode) 225520662ed9SBarry Smith and is also the same in the output `PetscSF`. 22561f40158dSVaclav Hapla For `PETSCSF_CONCATENATE_ROOTMODE_LOCAL` and `PETSCSF_CONCATENATE_ROOTMODE_GLOBAL`, the input root spaces are taken as separate and joined. 22571f40158dSVaclav Hapla `PETSCSF_CONCATENATE_ROOTMODE_LOCAL` joins the root spaces locally; 225820662ed9SBarry Smith roots of sfs[0], sfs[1], sfs[2], ... are joined on each rank separately, ordered by input `PetscSF` and original local index, and renumbered contiguously. 22591f40158dSVaclav Hapla `PETSCSF_CONCATENATE_ROOTMODE_GLOBAL` joins the root spaces globally; 226020662ed9SBarry Smith roots of sfs[0], sfs[1], sfs[2, ... are joined globally, ordered by input `PetscSF` index and original global index, and renumbered contiguously; 22611f40158dSVaclav Hapla the original root ranks are ignored. 22621f40158dSVaclav Hapla For both `PETSCSF_CONCATENATE_ROOTMODE_LOCAL` and `PETSCSF_CONCATENATE_ROOTMODE_GLOBAL`, 226320662ed9SBarry Smith the output `PetscSF`'s root layout is such that the local number of roots is a sum of the input `PetscSF`'s local numbers of roots on each rank 226420662ed9SBarry Smith to keep the load balancing. 226520662ed9SBarry Smith However, for `PETSCSF_CONCATENATE_ROOTMODE_GLOBAL`, roots can move to different ranks. 22661f40158dSVaclav Hapla 22671f40158dSVaclav Hapla Example: 22681f40158dSVaclav Hapla We can use src/vec/is/sf/tests/ex18.c to compare the root modes. By running 226920662ed9SBarry Smith .vb 227020662ed9SBarry Smith make -C $PETSC_DIR/src/vec/is/sf/tests ex18 227120662ed9SBarry Smith for m in {local,global,shared}; do 227220662ed9SBarry Smith mpirun -n 2 $PETSC_DIR/src/vec/is/sf/tests/ex18 -nsfs 2 -n 2 -root_mode $m -sf_view 227320662ed9SBarry Smith done 227420662ed9SBarry Smith .ve 227520662ed9SBarry Smith we generate two identical `PetscSF`s sf_0 and sf_1, 227620662ed9SBarry Smith .vb 227720662ed9SBarry Smith PetscSF Object: sf_0 2 MPI processes 227820662ed9SBarry Smith type: basic 227920662ed9SBarry Smith rank #leaves #roots 228020662ed9SBarry Smith [ 0] 4 2 228120662ed9SBarry Smith [ 1] 4 2 228220662ed9SBarry Smith leaves roots roots in global numbering 228320662ed9SBarry Smith ( 0, 0) <- ( 0, 0) = 0 228420662ed9SBarry Smith ( 0, 1) <- ( 0, 1) = 1 228520662ed9SBarry Smith ( 0, 2) <- ( 1, 0) = 2 228620662ed9SBarry Smith ( 0, 3) <- ( 1, 1) = 3 228720662ed9SBarry Smith ( 1, 0) <- ( 0, 0) = 0 228820662ed9SBarry Smith ( 1, 1) <- ( 0, 1) = 1 228920662ed9SBarry Smith ( 1, 2) <- ( 1, 0) = 2 229020662ed9SBarry Smith ( 1, 3) <- ( 1, 1) = 3 229120662ed9SBarry Smith .ve 229220662ed9SBarry Smith and pass them to `PetscSFConcatenate()` along with different choices of `rootMode`, yielding different result_sf: 229320662ed9SBarry Smith .vb 229420662ed9SBarry Smith rootMode = local: 229520662ed9SBarry Smith PetscSF Object: result_sf 2 MPI processes 229620662ed9SBarry Smith type: basic 229720662ed9SBarry Smith rank #leaves #roots 229820662ed9SBarry Smith [ 0] 8 4 229920662ed9SBarry Smith [ 1] 8 4 230020662ed9SBarry Smith leaves roots roots in global numbering 230120662ed9SBarry Smith ( 0, 0) <- ( 0, 0) = 0 230220662ed9SBarry Smith ( 0, 1) <- ( 0, 1) = 1 230320662ed9SBarry Smith ( 0, 2) <- ( 1, 0) = 4 230420662ed9SBarry Smith ( 0, 3) <- ( 1, 1) = 5 230520662ed9SBarry Smith ( 0, 4) <- ( 0, 2) = 2 230620662ed9SBarry Smith ( 0, 5) <- ( 0, 3) = 3 230720662ed9SBarry Smith ( 0, 6) <- ( 1, 2) = 6 230820662ed9SBarry Smith ( 0, 7) <- ( 1, 3) = 7 230920662ed9SBarry Smith ( 1, 0) <- ( 0, 0) = 0 231020662ed9SBarry Smith ( 1, 1) <- ( 0, 1) = 1 231120662ed9SBarry Smith ( 1, 2) <- ( 1, 0) = 4 231220662ed9SBarry Smith ( 1, 3) <- ( 1, 1) = 5 231320662ed9SBarry Smith ( 1, 4) <- ( 0, 2) = 2 231420662ed9SBarry Smith ( 1, 5) <- ( 0, 3) = 3 231520662ed9SBarry Smith ( 1, 6) <- ( 1, 2) = 6 231620662ed9SBarry Smith ( 1, 7) <- ( 1, 3) = 7 231720662ed9SBarry Smith 231820662ed9SBarry Smith rootMode = global: 231920662ed9SBarry Smith PetscSF Object: result_sf 2 MPI processes 232020662ed9SBarry Smith type: basic 232120662ed9SBarry Smith rank #leaves #roots 232220662ed9SBarry Smith [ 0] 8 4 232320662ed9SBarry Smith [ 1] 8 4 232420662ed9SBarry Smith leaves roots roots in global numbering 232520662ed9SBarry Smith ( 0, 0) <- ( 0, 0) = 0 232620662ed9SBarry Smith ( 0, 1) <- ( 0, 1) = 1 232720662ed9SBarry Smith ( 0, 2) <- ( 0, 2) = 2 232820662ed9SBarry Smith ( 0, 3) <- ( 0, 3) = 3 232920662ed9SBarry Smith ( 0, 4) <- ( 1, 0) = 4 233020662ed9SBarry Smith ( 0, 5) <- ( 1, 1) = 5 233120662ed9SBarry Smith ( 0, 6) <- ( 1, 2) = 6 233220662ed9SBarry Smith ( 0, 7) <- ( 1, 3) = 7 233320662ed9SBarry Smith ( 1, 0) <- ( 0, 0) = 0 233420662ed9SBarry Smith ( 1, 1) <- ( 0, 1) = 1 233520662ed9SBarry Smith ( 1, 2) <- ( 0, 2) = 2 233620662ed9SBarry Smith ( 1, 3) <- ( 0, 3) = 3 233720662ed9SBarry Smith ( 1, 4) <- ( 1, 0) = 4 233820662ed9SBarry Smith ( 1, 5) <- ( 1, 1) = 5 233920662ed9SBarry Smith ( 1, 6) <- ( 1, 2) = 6 234020662ed9SBarry Smith ( 1, 7) <- ( 1, 3) = 7 234120662ed9SBarry Smith 234220662ed9SBarry Smith rootMode = shared: 234320662ed9SBarry Smith PetscSF Object: result_sf 2 MPI processes 234420662ed9SBarry Smith type: basic 234520662ed9SBarry Smith rank #leaves #roots 234620662ed9SBarry Smith [ 0] 8 2 234720662ed9SBarry Smith [ 1] 8 2 234820662ed9SBarry Smith leaves roots roots in global numbering 234920662ed9SBarry Smith ( 0, 0) <- ( 0, 0) = 0 235020662ed9SBarry Smith ( 0, 1) <- ( 0, 1) = 1 235120662ed9SBarry Smith ( 0, 2) <- ( 1, 0) = 2 235220662ed9SBarry Smith ( 0, 3) <- ( 1, 1) = 3 235320662ed9SBarry Smith ( 0, 4) <- ( 0, 0) = 0 235420662ed9SBarry Smith ( 0, 5) <- ( 0, 1) = 1 235520662ed9SBarry Smith ( 0, 6) <- ( 1, 0) = 2 235620662ed9SBarry Smith ( 0, 7) <- ( 1, 1) = 3 235720662ed9SBarry Smith ( 1, 0) <- ( 0, 0) = 0 235820662ed9SBarry Smith ( 1, 1) <- ( 0, 1) = 1 235920662ed9SBarry Smith ( 1, 2) <- ( 1, 0) = 2 236020662ed9SBarry Smith ( 1, 3) <- ( 1, 1) = 3 236120662ed9SBarry Smith ( 1, 4) <- ( 0, 0) = 0 236220662ed9SBarry Smith ( 1, 5) <- ( 0, 1) = 1 236320662ed9SBarry Smith ( 1, 6) <- ( 1, 0) = 2 236420662ed9SBarry Smith ( 1, 7) <- ( 1, 1) = 3 236520662ed9SBarry Smith .ve 23661f40158dSVaclav Hapla 23671f40158dSVaclav Hapla .seealso: `PetscSF`, `PetscSFCompose()`, `PetscSFGetGraph()`, `PetscSFSetGraph()`, `PetscSFConcatenateRootMode` 2368157edd7aSVaclav Hapla @*/ 23691f40158dSVaclav Hapla PetscErrorCode PetscSFConcatenate(MPI_Comm comm, PetscInt nsfs, PetscSF sfs[], PetscSFConcatenateRootMode rootMode, PetscInt leafOffsets[], PetscSF *newsf) 2370d71ae5a4SJacob Faibussowitsch { 2371157edd7aSVaclav Hapla PetscInt i, s, nLeaves, nRoots; 2372157edd7aSVaclav Hapla PetscInt *leafArrayOffsets; 2373157edd7aSVaclav Hapla PetscInt *ilocal_new; 2374157edd7aSVaclav Hapla PetscSFNode *iremote_new; 2375157edd7aSVaclav Hapla PetscBool all_ilocal_null = PETSC_FALSE; 23761f40158dSVaclav Hapla PetscLayout glayout = NULL; 23771f40158dSVaclav Hapla PetscInt *gremote = NULL; 23781f40158dSVaclav Hapla PetscMPIInt rank, size; 2379157edd7aSVaclav Hapla 2380157edd7aSVaclav Hapla PetscFunctionBegin; 238112f479c1SVaclav Hapla if (PetscDefined(USE_DEBUG)) { 2382157edd7aSVaclav Hapla PetscSF dummy; /* just to have a PetscObject on comm for input validation */ 2383157edd7aSVaclav Hapla 23849566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, &dummy)); 2385157edd7aSVaclav Hapla PetscValidLogicalCollectiveInt(dummy, nsfs, 2); 2386157edd7aSVaclav Hapla PetscValidPointer(sfs, 3); 2387157edd7aSVaclav Hapla for (i = 0; i < nsfs; i++) { 2388157edd7aSVaclav Hapla PetscValidHeaderSpecific(sfs[i], PETSCSF_CLASSID, 3); 2389157edd7aSVaclav Hapla PetscCheckSameComm(dummy, 1, sfs[i], 3); 2390157edd7aSVaclav Hapla } 23911f40158dSVaclav Hapla PetscValidLogicalCollectiveEnum(dummy, rootMode, 4); 2392157edd7aSVaclav Hapla if (leafOffsets) PetscValidIntPointer(leafOffsets, 5); 2393157edd7aSVaclav Hapla PetscValidPointer(newsf, 6); 23949566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&dummy)); 2395157edd7aSVaclav Hapla } 2396157edd7aSVaclav Hapla if (!nsfs) { 23979566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, newsf)); 23989566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*newsf, 0, 0, NULL, PETSC_OWN_POINTER, NULL, PETSC_OWN_POINTER)); 23993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2400157edd7aSVaclav Hapla } 24019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 24021f40158dSVaclav Hapla PetscCallMPI(MPI_Comm_size(comm, &size)); 2403157edd7aSVaclav Hapla 24041f40158dSVaclav Hapla /* Calculate leaf array offsets */ 24059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsfs + 1, &leafArrayOffsets)); 2406157edd7aSVaclav Hapla leafArrayOffsets[0] = 0; 2407157edd7aSVaclav Hapla for (s = 0; s < nsfs; s++) { 2408157edd7aSVaclav Hapla PetscInt nl; 2409157edd7aSVaclav Hapla 24109566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfs[s], NULL, &nl, NULL, NULL)); 2411157edd7aSVaclav Hapla leafArrayOffsets[s + 1] = leafArrayOffsets[s] + nl; 2412157edd7aSVaclav Hapla } 2413157edd7aSVaclav Hapla nLeaves = leafArrayOffsets[nsfs]; 2414157edd7aSVaclav Hapla 24151f40158dSVaclav Hapla /* Calculate number of roots */ 24161f40158dSVaclav Hapla switch (rootMode) { 24171f40158dSVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_SHARED: { 24181f40158dSVaclav Hapla PetscCall(PetscSFGetGraph(sfs[0], &nRoots, NULL, NULL, NULL)); 24191f40158dSVaclav Hapla if (PetscDefined(USE_DEBUG)) { 24201f40158dSVaclav Hapla for (s = 1; s < nsfs; s++) { 24211f40158dSVaclav Hapla PetscInt nr; 24221f40158dSVaclav Hapla 24231f40158dSVaclav Hapla PetscCall(PetscSFGetGraph(sfs[s], &nr, NULL, NULL, NULL)); 24241f40158dSVaclav Hapla PetscCheck(nr == nRoots, comm, PETSC_ERR_ARG_SIZ, "rootMode = %s but sfs[%" PetscInt_FMT "] has a different number of roots (%" PetscInt_FMT ") than sfs[0] (%" PetscInt_FMT ")", PetscSFConcatenateRootModes[rootMode], s, nr, nRoots); 24251f40158dSVaclav Hapla } 24261f40158dSVaclav Hapla } 24271f40158dSVaclav Hapla } break; 24281f40158dSVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_GLOBAL: { 24291f40158dSVaclav Hapla /* Calculate also global layout in this case */ 24301f40158dSVaclav Hapla PetscInt *nls; 24311f40158dSVaclav Hapla PetscLayout *lts; 24321f40158dSVaclav Hapla PetscInt **inds; 24331f40158dSVaclav Hapla PetscInt j; 24341f40158dSVaclav Hapla PetscInt rootOffset = 0; 24351f40158dSVaclav Hapla 24361f40158dSVaclav Hapla PetscCall(PetscCalloc3(nsfs, <s, nsfs, &nls, nsfs, &inds)); 24371f40158dSVaclav Hapla PetscCall(PetscLayoutCreate(comm, &glayout)); 24381f40158dSVaclav Hapla glayout->bs = 1; 24391f40158dSVaclav Hapla glayout->n = 0; 24401f40158dSVaclav Hapla glayout->N = 0; 24411f40158dSVaclav Hapla for (s = 0; s < nsfs; s++) { 24421f40158dSVaclav Hapla PetscCall(PetscSFGetGraphLayout(sfs[s], <s[s], &nls[s], NULL, &inds[s])); 24431f40158dSVaclav Hapla glayout->n += lts[s]->n; 24441f40158dSVaclav Hapla glayout->N += lts[s]->N; 24451f40158dSVaclav Hapla } 24461f40158dSVaclav Hapla PetscCall(PetscLayoutSetUp(glayout)); 24471f40158dSVaclav Hapla PetscCall(PetscMalloc1(nLeaves, &gremote)); 24481f40158dSVaclav Hapla for (s = 0, j = 0; s < nsfs; s++) { 24491f40158dSVaclav Hapla for (i = 0; i < nls[s]; i++, j++) gremote[j] = inds[s][i] + rootOffset; 24501f40158dSVaclav Hapla rootOffset += lts[s]->N; 24511f40158dSVaclav Hapla PetscCall(PetscLayoutDestroy(<s[s])); 24521f40158dSVaclav Hapla PetscCall(PetscFree(inds[s])); 24531f40158dSVaclav Hapla } 24541f40158dSVaclav Hapla PetscCall(PetscFree3(lts, nls, inds)); 24551f40158dSVaclav Hapla nRoots = glayout->N; 24561f40158dSVaclav Hapla } break; 24571f40158dSVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_LOCAL: 24581f40158dSVaclav Hapla /* nRoots calculated later in this case */ 24591f40158dSVaclav Hapla break; 24601f40158dSVaclav Hapla default: 24611f40158dSVaclav Hapla SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid PetscSFConcatenateRootMode %d", rootMode); 24621f40158dSVaclav Hapla } 24631f40158dSVaclav Hapla 2464157edd7aSVaclav Hapla if (!leafOffsets) { 2465157edd7aSVaclav Hapla all_ilocal_null = PETSC_TRUE; 2466157edd7aSVaclav Hapla for (s = 0; s < nsfs; s++) { 2467157edd7aSVaclav Hapla const PetscInt *ilocal; 2468157edd7aSVaclav Hapla 24699566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfs[s], NULL, NULL, &ilocal, NULL)); 2470157edd7aSVaclav Hapla if (ilocal) { 2471157edd7aSVaclav Hapla all_ilocal_null = PETSC_FALSE; 2472157edd7aSVaclav Hapla break; 2473157edd7aSVaclav Hapla } 2474157edd7aSVaclav Hapla } 2475157edd7aSVaclav 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"); 2476157edd7aSVaclav Hapla } 2477157edd7aSVaclav Hapla 2478157edd7aSVaclav Hapla /* Renumber and concatenate local leaves */ 2479157edd7aSVaclav Hapla ilocal_new = NULL; 2480157edd7aSVaclav Hapla if (!all_ilocal_null) { 24819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nLeaves, &ilocal_new)); 2482157edd7aSVaclav Hapla for (i = 0; i < nLeaves; i++) ilocal_new[i] = -1; 2483157edd7aSVaclav Hapla for (s = 0; s < nsfs; s++) { 2484157edd7aSVaclav Hapla const PetscInt *ilocal; 2485157edd7aSVaclav Hapla PetscInt *ilocal_l = &ilocal_new[leafArrayOffsets[s]]; 2486157edd7aSVaclav Hapla PetscInt i, nleaves_l; 2487157edd7aSVaclav Hapla 24889566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfs[s], NULL, &nleaves_l, &ilocal, NULL)); 2489157edd7aSVaclav Hapla for (i = 0; i < nleaves_l; i++) ilocal_l[i] = (ilocal ? ilocal[i] : i) + leafOffsets[s]; 2490157edd7aSVaclav Hapla } 2491157edd7aSVaclav Hapla } 2492157edd7aSVaclav Hapla 2493157edd7aSVaclav Hapla /* Renumber and concatenate remote roots */ 24941f40158dSVaclav Hapla if (rootMode == PETSCSF_CONCATENATE_ROOTMODE_LOCAL || rootMode == PETSCSF_CONCATENATE_ROOTMODE_SHARED) { 24951f40158dSVaclav Hapla PetscInt rootOffset = 0; 24961f40158dSVaclav Hapla 24979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nLeaves, &iremote_new)); 2498157edd7aSVaclav Hapla for (i = 0; i < nLeaves; i++) { 2499157edd7aSVaclav Hapla iremote_new[i].rank = -1; 2500157edd7aSVaclav Hapla iremote_new[i].index = -1; 2501157edd7aSVaclav Hapla } 2502157edd7aSVaclav Hapla for (s = 0; s < nsfs; s++) { 2503157edd7aSVaclav Hapla PetscInt i, nl, nr; 2504157edd7aSVaclav Hapla PetscSF tmp_sf; 2505157edd7aSVaclav Hapla const PetscSFNode *iremote; 2506157edd7aSVaclav Hapla PetscSFNode *tmp_rootdata; 2507157edd7aSVaclav Hapla PetscSFNode *tmp_leafdata = &iremote_new[leafArrayOffsets[s]]; 2508157edd7aSVaclav Hapla 25099566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfs[s], &nr, &nl, NULL, &iremote)); 25109566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, &tmp_sf)); 2511157edd7aSVaclav Hapla /* create helper SF with contiguous leaves */ 25129566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(tmp_sf, nr, nl, NULL, PETSC_USE_POINTER, (PetscSFNode *)iremote, PETSC_COPY_VALUES)); 25139566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(tmp_sf)); 25149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nr, &tmp_rootdata)); 25151f40158dSVaclav Hapla if (rootMode == PETSCSF_CONCATENATE_ROOTMODE_LOCAL) { 2516157edd7aSVaclav Hapla for (i = 0; i < nr; i++) { 25171f40158dSVaclav Hapla tmp_rootdata[i].index = i + rootOffset; 2518157edd7aSVaclav Hapla tmp_rootdata[i].rank = (PetscInt)rank; 2519157edd7aSVaclav Hapla } 25201f40158dSVaclav Hapla rootOffset += nr; 25211f40158dSVaclav Hapla } else { 25221f40158dSVaclav Hapla for (i = 0; i < nr; i++) { 25231f40158dSVaclav Hapla tmp_rootdata[i].index = i; 25241f40158dSVaclav Hapla tmp_rootdata[i].rank = (PetscInt)rank; 25251f40158dSVaclav Hapla } 25261f40158dSVaclav Hapla } 25279566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(tmp_sf, MPIU_2INT, tmp_rootdata, tmp_leafdata, MPI_REPLACE)); 25289566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(tmp_sf, MPIU_2INT, tmp_rootdata, tmp_leafdata, MPI_REPLACE)); 25299566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&tmp_sf)); 25309566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp_rootdata)); 2531157edd7aSVaclav Hapla } 2532aa624791SPierre Jolivet if (rootMode == PETSCSF_CONCATENATE_ROOTMODE_LOCAL) nRoots = rootOffset; // else nRoots already calculated above 2533157edd7aSVaclav Hapla 2534157edd7aSVaclav Hapla /* Build the new SF */ 25359566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, newsf)); 25369566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*newsf, nRoots, nLeaves, ilocal_new, PETSC_OWN_POINTER, iremote_new, PETSC_OWN_POINTER)); 25371f40158dSVaclav Hapla } else { 25381f40158dSVaclav Hapla /* Build the new SF */ 25391f40158dSVaclav Hapla PetscCall(PetscSFCreate(comm, newsf)); 25401f40158dSVaclav Hapla PetscCall(PetscSFSetGraphLayout(*newsf, glayout, nLeaves, ilocal_new, PETSC_OWN_POINTER, gremote)); 25411f40158dSVaclav Hapla } 25429566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(*newsf)); 25431f40158dSVaclav Hapla PetscCall(PetscSFViewFromOptions(*newsf, NULL, "-sf_concat_view")); 25441f40158dSVaclav Hapla PetscCall(PetscLayoutDestroy(&glayout)); 25451f40158dSVaclav Hapla PetscCall(PetscFree(gremote)); 25469566063dSJacob Faibussowitsch PetscCall(PetscFree(leafArrayOffsets)); 25473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2548157edd7aSVaclav Hapla } 2549