1f210f596SVaclav Hapla 2f210f596SVaclav Hapla static char help[] = "Test PetscSFConcatenate()\n\n"; 3f210f596SVaclav Hapla 4f210f596SVaclav Hapla #include <petscsf.h> 5f210f596SVaclav Hapla 6f210f596SVaclav Hapla typedef struct { 7f210f596SVaclav Hapla MPI_Comm comm; 8f210f596SVaclav Hapla PetscMPIInt rank, size; 9*c1730cc1SVaclav Hapla PetscInt leaveStep, nsfs, n; 10*c1730cc1SVaclav Hapla PetscBool sparseLeaves; 11*c1730cc1SVaclav Hapla PetscBool compare; 12*c1730cc1SVaclav Hapla PetscBool irregular; 13*c1730cc1SVaclav Hapla PetscSFConcatenateRootMode rootMode; 14*c1730cc1SVaclav Hapla PetscViewer viewer; 15f210f596SVaclav Hapla } AppCtx; 16f210f596SVaclav Hapla 17d71ae5a4SJacob Faibussowitsch static PetscErrorCode GetOptions(MPI_Comm comm, AppCtx *ctx) 18d71ae5a4SJacob Faibussowitsch { 19*c1730cc1SVaclav Hapla PetscViewerFormat format; 20*c1730cc1SVaclav Hapla 21f210f596SVaclav Hapla PetscFunctionBegin; 22f210f596SVaclav Hapla ctx->comm = comm; 23f210f596SVaclav Hapla ctx->nsfs = 3; 24*c1730cc1SVaclav Hapla ctx->n = 1; 25f210f596SVaclav Hapla ctx->leaveStep = 1; 26f210f596SVaclav Hapla ctx->sparseLeaves = PETSC_FALSE; 27*c1730cc1SVaclav Hapla ctx->compare = PETSC_FALSE; 28*c1730cc1SVaclav Hapla ctx->irregular = PETSC_FALSE; 29*c1730cc1SVaclav Hapla ctx->rootMode = PETSCSF_CONCATENATE_ROOTMODE_LOCAL; 30*c1730cc1SVaclav Hapla ctx->viewer = NULL; 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nsfs", &ctx->nsfs, NULL)); 32*c1730cc1SVaclav Hapla PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &ctx->n, NULL)); 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-leave_step", &ctx->leaveStep, NULL)); 34*c1730cc1SVaclav Hapla PetscCall(PetscOptionsGetBool(NULL, NULL, "-irregular", &ctx->irregular, NULL)); 35*c1730cc1SVaclav Hapla PetscCall(PetscOptionsGetBool(NULL, NULL, "-compare_to_reference", &ctx->compare, NULL)); 36*c1730cc1SVaclav Hapla PetscCall(PetscOptionsGetEnum(NULL, NULL, "-root_mode", PetscSFConcatenateRootModes, (PetscEnum *)&ctx->rootMode, NULL)); 37*c1730cc1SVaclav Hapla PetscCall(PetscOptionsGetViewer(comm, NULL, NULL, "-sf_view", &ctx->viewer, &format, NULL)); 38*c1730cc1SVaclav Hapla if (ctx->viewer) PetscCall(PetscViewerPushFormat(ctx->viewer, format)); 39f210f596SVaclav Hapla ctx->sparseLeaves = (PetscBool)(ctx->leaveStep != 1); 409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &ctx->size)); 419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &ctx->rank)); 42f210f596SVaclav Hapla PetscFunctionReturn(0); 43f210f596SVaclav Hapla } 44f210f596SVaclav Hapla 45d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFCheckEqual_Private(PetscSF sf0, PetscSF sf1) 46d71ae5a4SJacob Faibussowitsch { 47f210f596SVaclav Hapla PetscInt nRoot, nLeave; 48f210f596SVaclav Hapla Vec vecRoot0, vecLeave0, vecRoot1, vecLeave1; 49f210f596SVaclav Hapla MPI_Comm comm; 50f210f596SVaclav Hapla PetscBool flg; 51f210f596SVaclav Hapla 52f210f596SVaclav Hapla PetscFunctionBegin; 539566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf0, &comm)); 549566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf0, &nRoot, NULL, NULL, NULL)); 559566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(sf0, NULL, &nLeave)); 56f210f596SVaclav Hapla nLeave++; 579566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(comm, nRoot, PETSC_DECIDE, &vecRoot0)); 589566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(comm, nLeave, PETSC_DECIDE, &vecLeave0)); 599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(vecRoot0, &vecRoot1)); 609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(vecLeave0, &vecLeave1)); 61f210f596SVaclav Hapla { 62f210f596SVaclav Hapla PetscRandom rand; 63f210f596SVaclav Hapla 649566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rand)); 659566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 669566063dSJacob Faibussowitsch PetscCall(VecSetRandom(vecRoot0, rand)); 679566063dSJacob Faibussowitsch PetscCall(VecSetRandom(vecLeave0, rand)); 689566063dSJacob Faibussowitsch PetscCall(VecCopy(vecRoot0, vecRoot1)); 699566063dSJacob Faibussowitsch PetscCall(VecCopy(vecLeave0, vecLeave1)); 709566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 71f210f596SVaclav Hapla } 72f210f596SVaclav Hapla 739566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD)); 749566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD)); 759566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD)); 769566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD)); 779566063dSJacob Faibussowitsch PetscCall(VecEqual(vecLeave0, vecLeave1, &flg)); 78f210f596SVaclav Hapla PetscCheck(flg, comm, PETSC_ERR_PLIB, "leave vectors differ"); 79f210f596SVaclav Hapla 809566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE)); 819566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE)); 829566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE)); 839566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE)); 849566063dSJacob Faibussowitsch PetscCall(VecEqual(vecRoot0, vecRoot1, &flg)); 85f210f596SVaclav Hapla PetscCheck(flg, comm, PETSC_ERR_PLIB, "root vectors differ"); 86f210f596SVaclav Hapla 879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecRoot0)); 889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecRoot1)); 899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecLeave0)); 909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecLeave1)); 91f210f596SVaclav Hapla PetscFunctionReturn(0); 92f210f596SVaclav Hapla } 93f210f596SVaclav Hapla 94*c1730cc1SVaclav Hapla static PetscErrorCode PetscSFViewCustom(PetscSF sf, PetscViewer viewer) 95d71ae5a4SJacob Faibussowitsch { 96*c1730cc1SVaclav Hapla PetscMPIInt rank; 97*c1730cc1SVaclav Hapla PetscInt i, nroots, nleaves, nranks; 98*c1730cc1SVaclav Hapla const PetscInt *ilocal; 99*c1730cc1SVaclav Hapla const PetscSFNode *iremote; 100*c1730cc1SVaclav Hapla PetscLayout rootLayout; 101*c1730cc1SVaclav Hapla PetscInt *gremote; 102*c1730cc1SVaclav Hapla 103*c1730cc1SVaclav Hapla PetscFunctionBegin; 104*c1730cc1SVaclav Hapla PetscCall(PetscSFSetUp(sf)); 105*c1730cc1SVaclav Hapla PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, &ilocal, &iremote)); 106*c1730cc1SVaclav Hapla PetscCall(PetscSFGetRootRanks(sf, &nranks, NULL, NULL, NULL, NULL)); 107*c1730cc1SVaclav Hapla PetscCall(PetscSFGetGraphLayout(sf, &rootLayout, NULL, NULL, &gremote)); 108*c1730cc1SVaclav Hapla PetscCheck(nroots == rootLayout->n, PetscObjectComm((PetscObject)sf), PETSC_ERR_PLIB, "Assertion failed: nroots == rootLayout->n"); 109*c1730cc1SVaclav Hapla PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)sf, viewer)); 110*c1730cc1SVaclav Hapla PetscCall(PetscViewerASCIIPushTab(viewer)); 111*c1730cc1SVaclav Hapla PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank)); 112*c1730cc1SVaclav Hapla PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 113*c1730cc1SVaclav Hapla if (rank == 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "rank #leaves #roots\n")); 114*c1730cc1SVaclav Hapla PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%2d] %7" PetscInt_FMT " %6" PetscInt_FMT "\n", rank, nleaves, nroots)); 115*c1730cc1SVaclav Hapla PetscCall(PetscViewerFlush(viewer)); 116*c1730cc1SVaclav Hapla if (rank == 0) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "leaves roots roots in global numbering\n")); 117*c1730cc1SVaclav Hapla for (i = 0; i < nleaves; i++) 118*c1730cc1SVaclav Hapla PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "(%2d, %2" PetscInt_FMT ") <- (%2" PetscInt_FMT ", %2" PetscInt_FMT ") = %3" PetscInt_FMT "\n", rank, ilocal ? ilocal[i] : i, iremote[i].rank, iremote[i].index, gremote[i])); 119*c1730cc1SVaclav Hapla PetscCall(PetscViewerFlush(viewer)); 120*c1730cc1SVaclav Hapla PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 121*c1730cc1SVaclav Hapla PetscCall(PetscViewerASCIIPopTab(viewer)); 122*c1730cc1SVaclav Hapla PetscCall(PetscLayoutDestroy(&rootLayout)); 123*c1730cc1SVaclav Hapla PetscCall(PetscFree(gremote)); 124*c1730cc1SVaclav Hapla PetscFunctionReturn(0); 125*c1730cc1SVaclav Hapla } 126*c1730cc1SVaclav Hapla 127*c1730cc1SVaclav Hapla PetscErrorCode CreateReferenceSF_Regular(AppCtx *ctx, PetscSF *refSF) 128*c1730cc1SVaclav Hapla { 129*c1730cc1SVaclav Hapla PetscInt j; 130f210f596SVaclav Hapla PetscInt *ilocal = NULL; 131*c1730cc1SVaclav Hapla PetscInt nLeaves = ctx->nsfs * ctx->n * ctx->size; 132*c1730cc1SVaclav Hapla PetscInt nroots = ctx->n * ctx->nsfs; 133f210f596SVaclav Hapla PetscSF sf; 134f210f596SVaclav Hapla 135f210f596SVaclav Hapla PetscFunctionBegin; 136f210f596SVaclav Hapla ilocal = NULL; 13748a46eb9SPierre Jolivet if (ctx->sparseLeaves) PetscCall(PetscCalloc1(nLeaves + 1, &ilocal)); 1389566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(ctx->comm, &sf)); 139*c1730cc1SVaclav Hapla for (j = 0; j < nLeaves; j++) { 140*c1730cc1SVaclav Hapla if (ctx->sparseLeaves) ilocal[j + 1] = ilocal[j] + ctx->leaveStep; 141*c1730cc1SVaclav Hapla } 142*c1730cc1SVaclav Hapla switch (ctx->rootMode) { 143*c1730cc1SVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_SHARED: 144*c1730cc1SVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_LOCAL: { 145*c1730cc1SVaclav Hapla PetscInt i, k; 146*c1730cc1SVaclav Hapla PetscMPIInt r; 147*c1730cc1SVaclav Hapla PetscSFNode *iremote; 148*c1730cc1SVaclav Hapla 149*c1730cc1SVaclav Hapla PetscCall(PetscCalloc1(nLeaves, &iremote)); 150f210f596SVaclav Hapla for (i = 0, j = 0; i < ctx->nsfs; i++) { 151f210f596SVaclav Hapla for (r = 0; r < ctx->size; r++) { 152*c1730cc1SVaclav Hapla for (k = 0; k < ctx->n; k++, j++) { 153f210f596SVaclav Hapla iremote[j].rank = r; 154*c1730cc1SVaclav Hapla iremote[j].index = k + i * ctx->n; 155f210f596SVaclav Hapla } 156f210f596SVaclav Hapla } 157f210f596SVaclav Hapla } 1589566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 159*c1730cc1SVaclav Hapla } break; 160*c1730cc1SVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_GLOBAL: { 161*c1730cc1SVaclav Hapla PetscLayout map = NULL; 162*c1730cc1SVaclav Hapla PetscInt *gremote; 163*c1730cc1SVaclav Hapla 164*c1730cc1SVaclav Hapla PetscCall(PetscLayoutCreateFromSizes(ctx->comm, nroots, PETSC_DECIDE, 1, &map)); 165*c1730cc1SVaclav Hapla PetscCall(PetscMalloc1(nLeaves, &gremote)); 166*c1730cc1SVaclav Hapla for (j = 0; j < nLeaves; j++) gremote[j] = j; 167*c1730cc1SVaclav Hapla PetscCall(PetscSFSetGraphLayout(sf, map, nLeaves, ilocal, PETSC_OWN_POINTER, gremote)); 168*c1730cc1SVaclav Hapla PetscCall(PetscFree(gremote)); 169*c1730cc1SVaclav Hapla PetscCall(PetscLayoutDestroy(&map)); 170*c1730cc1SVaclav Hapla } break; 171*c1730cc1SVaclav Hapla default: 172*c1730cc1SVaclav Hapla SETERRQ(ctx->comm, PETSC_ERR_SUP, "unsupported rootmode %d", ctx->rootMode); 173*c1730cc1SVaclav Hapla } 174*c1730cc1SVaclav Hapla PetscCall(PetscObjectSetName((PetscObject)sf, "reference_sf")); 175*c1730cc1SVaclav Hapla if (ctx->viewer) PetscCall(PetscSFViewCustom(sf, ctx->viewer)); 176f210f596SVaclav Hapla *refSF = sf; 177f210f596SVaclav Hapla PetscFunctionReturn(0); 178f210f596SVaclav Hapla } 179f210f596SVaclav Hapla 180*c1730cc1SVaclav Hapla PetscErrorCode CreateSFs_Irregular(AppCtx *ctx, PetscSF *newSFs[], PetscInt *leafOffsets[]) 181d71ae5a4SJacob Faibussowitsch { 182f210f596SVaclav Hapla PetscInt i; 183f210f596SVaclav Hapla PetscInt *lOffsets = NULL; 184f210f596SVaclav Hapla PetscSF *sfs; 185*c1730cc1SVaclav Hapla PetscInt nLeaves = ctx->n * ctx->size + (ctx->size - 1) * ctx->size / 2; 186*c1730cc1SVaclav Hapla PetscInt nroots = ctx->n + ctx->rank + ctx->nsfs - 1 + ctx->size - 1; 187f210f596SVaclav Hapla 188f210f596SVaclav Hapla PetscFunctionBegin; 18948a46eb9SPierre Jolivet if (ctx->sparseLeaves) PetscCall(PetscCalloc1(ctx->nsfs + 1, &lOffsets)); 1909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ctx->nsfs, &sfs)); 191f210f596SVaclav Hapla for (i = 0; i < ctx->nsfs; i++) { 192*c1730cc1SVaclav Hapla PetscSF sf; 193f210f596SVaclav Hapla PetscInt j, k; 194f210f596SVaclav Hapla PetscMPIInt r; 195f210f596SVaclav Hapla PetscInt *ilocal = NULL; 196f210f596SVaclav Hapla PetscSFNode *iremote; 197*c1730cc1SVaclav Hapla char name[32]; 198f210f596SVaclav Hapla 19948a46eb9SPierre Jolivet if (ctx->sparseLeaves) PetscCall(PetscCalloc1(nLeaves + 1, &ilocal)); 2009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nLeaves, &iremote)); 201*c1730cc1SVaclav Hapla for (r = ctx->size - 1, j = 0; r >= 0; r--) { 202*c1730cc1SVaclav Hapla for (k = 0; k < ctx->n + r; k++, j++) { 203ad540459SPierre Jolivet if (ctx->sparseLeaves) ilocal[j + 1] = ilocal[j] + ctx->leaveStep; 204f210f596SVaclav Hapla iremote[j].rank = r; 205*c1730cc1SVaclav Hapla iremote[j].index = k + i + ctx->rank; 206f210f596SVaclav Hapla } 207f210f596SVaclav Hapla } 208f210f596SVaclav Hapla if (ctx->sparseLeaves) lOffsets[i + 1] = lOffsets[i] + ilocal[j]; 209f210f596SVaclav Hapla 210*c1730cc1SVaclav Hapla PetscCall(PetscSFCreate(ctx->comm, &sf)); 211*c1730cc1SVaclav Hapla PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 212*c1730cc1SVaclav Hapla PetscCall(PetscSNPrintf(name, sizeof(name), "sf_%" PetscInt_FMT, i)); 213*c1730cc1SVaclav Hapla PetscCall(PetscObjectSetName((PetscObject)sf, name)); 214*c1730cc1SVaclav Hapla if (ctx->viewer) PetscCall(PetscSFViewCustom(sf, ctx->viewer)); 215*c1730cc1SVaclav Hapla sfs[i] = sf; 216*c1730cc1SVaclav Hapla } 217*c1730cc1SVaclav Hapla *newSFs = sfs; 218*c1730cc1SVaclav Hapla *leafOffsets = lOffsets; 219*c1730cc1SVaclav Hapla PetscFunctionReturn(0); 220*c1730cc1SVaclav Hapla } 221*c1730cc1SVaclav Hapla 222*c1730cc1SVaclav Hapla PetscErrorCode CreateSFs_Regular(AppCtx *ctx, PetscSF *newSFs[], PetscInt *leafOffsets[]) 223*c1730cc1SVaclav Hapla { 224*c1730cc1SVaclav Hapla PetscInt i; 225*c1730cc1SVaclav Hapla PetscInt *lOffsets = NULL; 226*c1730cc1SVaclav Hapla PetscInt nLeaves = ctx->n * ctx->size; 227*c1730cc1SVaclav Hapla PetscSF *sfs; 228*c1730cc1SVaclav Hapla PetscSFConcatenateRootMode mode = ctx->compare ? ctx->rootMode : PETSCSF_CONCATENATE_ROOTMODE_LOCAL; 229*c1730cc1SVaclav Hapla 230*c1730cc1SVaclav Hapla PetscFunctionBegin; 231*c1730cc1SVaclav Hapla if (ctx->sparseLeaves) PetscCall(PetscCalloc1(ctx->nsfs + 1, &lOffsets)); 232*c1730cc1SVaclav Hapla PetscCall(PetscCalloc1(ctx->nsfs, &sfs)); 233*c1730cc1SVaclav Hapla for (i = 0; i < ctx->nsfs; i++) { 234*c1730cc1SVaclav Hapla PetscSF sf; 235*c1730cc1SVaclav Hapla PetscInt j; 236*c1730cc1SVaclav Hapla PetscInt *ilocal = NULL; 237*c1730cc1SVaclav Hapla char name[32]; 238*c1730cc1SVaclav Hapla 239*c1730cc1SVaclav Hapla PetscCall(PetscSFCreate(ctx->comm, &sf)); 240*c1730cc1SVaclav Hapla if (ctx->sparseLeaves) { 241*c1730cc1SVaclav Hapla PetscCall(PetscCalloc1(nLeaves + 1, &ilocal)); 242*c1730cc1SVaclav Hapla for (j = 0; j < nLeaves; j++) { ilocal[j + 1] = ilocal[j] + ctx->leaveStep; } 243*c1730cc1SVaclav Hapla lOffsets[i + 1] = lOffsets[i] + ilocal[nLeaves]; 244*c1730cc1SVaclav Hapla } 245*c1730cc1SVaclav Hapla switch (mode) { 246*c1730cc1SVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_LOCAL: { 247*c1730cc1SVaclav Hapla PetscInt k, nroots = ctx->n; 248*c1730cc1SVaclav Hapla PetscMPIInt r; 249*c1730cc1SVaclav Hapla PetscSFNode *iremote; 250*c1730cc1SVaclav Hapla 251*c1730cc1SVaclav Hapla PetscCall(PetscMalloc1(nLeaves, &iremote)); 252*c1730cc1SVaclav Hapla for (r = 0, j = 0; r < ctx->size; r++) { 253*c1730cc1SVaclav Hapla for (k = 0; k < ctx->n; k++, j++) { 254*c1730cc1SVaclav Hapla iremote[j].rank = r; 255*c1730cc1SVaclav Hapla iremote[j].index = k; 256*c1730cc1SVaclav Hapla } 257*c1730cc1SVaclav Hapla } 258*c1730cc1SVaclav Hapla PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 259*c1730cc1SVaclav Hapla } break; 260*c1730cc1SVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_SHARED: { 261*c1730cc1SVaclav Hapla PetscInt k, nroots = ctx->n * ctx->nsfs; 262*c1730cc1SVaclav Hapla PetscMPIInt r; 263*c1730cc1SVaclav Hapla PetscSFNode *iremote; 264*c1730cc1SVaclav Hapla 265*c1730cc1SVaclav Hapla PetscCall(PetscMalloc1(nLeaves, &iremote)); 266*c1730cc1SVaclav Hapla for (r = 0, j = 0; r < ctx->size; r++) { 267*c1730cc1SVaclav Hapla for (k = 0; k < ctx->n; k++, j++) { 268*c1730cc1SVaclav Hapla iremote[j].rank = r; 269*c1730cc1SVaclav Hapla iremote[j].index = k + i * ctx->n; 270*c1730cc1SVaclav Hapla } 271*c1730cc1SVaclav Hapla } 272*c1730cc1SVaclav Hapla PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 273*c1730cc1SVaclav Hapla } break; 274*c1730cc1SVaclav Hapla case PETSCSF_CONCATENATE_ROOTMODE_GLOBAL: { 275*c1730cc1SVaclav Hapla PetscInt nroots = ctx->n; 276*c1730cc1SVaclav Hapla PetscLayout map = NULL; 277*c1730cc1SVaclav Hapla PetscInt *gremote; 278*c1730cc1SVaclav Hapla 279*c1730cc1SVaclav Hapla PetscCall(PetscLayoutCreateFromSizes(ctx->comm, nroots, PETSC_DECIDE, 1, &map)); 280*c1730cc1SVaclav Hapla PetscCall(PetscMalloc1(nLeaves, &gremote)); 281*c1730cc1SVaclav Hapla for (j = 0; j < nLeaves; j++) gremote[j] = j; 282*c1730cc1SVaclav Hapla PetscCall(PetscSFSetGraphLayout(sf, map, nLeaves, ilocal, PETSC_OWN_POINTER, gremote)); 283*c1730cc1SVaclav Hapla PetscCall(PetscFree(gremote)); 284*c1730cc1SVaclav Hapla PetscCall(PetscLayoutDestroy(&map)); 285*c1730cc1SVaclav Hapla } break; 286*c1730cc1SVaclav Hapla default: 287*c1730cc1SVaclav Hapla SETERRQ(ctx->comm, PETSC_ERR_SUP, "unsupported rootmode %d", ctx->rootMode); 288*c1730cc1SVaclav Hapla } 289*c1730cc1SVaclav Hapla PetscCall(PetscSNPrintf(name, sizeof(name), "sf_%" PetscInt_FMT, i)); 290*c1730cc1SVaclav Hapla PetscCall(PetscObjectSetName((PetscObject)sf, name)); 291*c1730cc1SVaclav Hapla if (ctx->viewer) PetscCall(PetscSFViewCustom(sf, ctx->viewer)); 292*c1730cc1SVaclav Hapla sfs[i] = sf; 293f210f596SVaclav Hapla } 294f210f596SVaclav Hapla *newSFs = sfs; 295f210f596SVaclav Hapla *leafOffsets = lOffsets; 296f210f596SVaclav Hapla PetscFunctionReturn(0); 297f210f596SVaclav Hapla } 298f210f596SVaclav Hapla 299d71ae5a4SJacob Faibussowitsch PetscErrorCode DestroySFs(AppCtx *ctx, PetscSF *sfs[]) 300d71ae5a4SJacob Faibussowitsch { 301f210f596SVaclav Hapla PetscInt i; 302f210f596SVaclav Hapla 303f210f596SVaclav Hapla PetscFunctionBegin; 30448a46eb9SPierre Jolivet for (i = 0; i < ctx->nsfs; i++) PetscCall(PetscSFDestroy(&(*sfs)[i])); 3059566063dSJacob Faibussowitsch PetscCall(PetscFree(*sfs)); 306f210f596SVaclav Hapla PetscFunctionReturn(0); 307f210f596SVaclav Hapla } 308f210f596SVaclav Hapla 309d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 310d71ae5a4SJacob Faibussowitsch { 311*c1730cc1SVaclav Hapla AppCtx ctx_; 312*c1730cc1SVaclav Hapla AppCtx *ctx = &ctx_; 313*c1730cc1SVaclav Hapla PetscSF sf; 314*c1730cc1SVaclav Hapla PetscSF *sfs = NULL; 315*c1730cc1SVaclav Hapla PetscInt *leafOffsets = NULL; 316f210f596SVaclav Hapla MPI_Comm comm; 317f210f596SVaclav Hapla 318327415f7SBarry Smith PetscFunctionBeginUser; 3199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 320f210f596SVaclav Hapla comm = PETSC_COMM_WORLD; 321*c1730cc1SVaclav Hapla PetscCall(GetOptions(comm, ctx)); 322f210f596SVaclav Hapla 323*c1730cc1SVaclav Hapla if (ctx->irregular) { 324*c1730cc1SVaclav Hapla PetscCall(CreateSFs_Irregular(ctx, &sfs, &leafOffsets)); 325*c1730cc1SVaclav Hapla } else { 326*c1730cc1SVaclav Hapla PetscCall(CreateSFs_Regular(ctx, &sfs, &leafOffsets)); 327*c1730cc1SVaclav Hapla } 328*c1730cc1SVaclav Hapla PetscCall(PetscSFConcatenate(comm, ctx->nsfs, sfs, ctx->rootMode, leafOffsets, &sf)); 329*c1730cc1SVaclav Hapla PetscCall(PetscObjectSetName((PetscObject)sf, "result_sf")); 330*c1730cc1SVaclav Hapla if (ctx->viewer) { 331*c1730cc1SVaclav Hapla PetscCall(PetscPrintf(comm, "rootMode = %s:\n", PetscSFConcatenateRootModes[ctx->rootMode])); 332*c1730cc1SVaclav Hapla PetscCall(PetscSFViewCustom(sf, ctx->viewer)); 333*c1730cc1SVaclav Hapla } 334*c1730cc1SVaclav Hapla if (ctx->compare) { 335*c1730cc1SVaclav Hapla PetscSF sfRef; 336*c1730cc1SVaclav Hapla 337*c1730cc1SVaclav Hapla PetscAssert(!ctx->irregular, comm, PETSC_ERR_SUP, "Combination -compare_to_reference true -irregular true not implemented"); 338*c1730cc1SVaclav Hapla PetscCall(CreateReferenceSF_Regular(ctx, &sfRef)); 3399566063dSJacob Faibussowitsch PetscCall(PetscSFCheckEqual_Private(sf, sfRef)); 340*c1730cc1SVaclav Hapla PetscCall(PetscSFDestroy(&sfRef)); 341*c1730cc1SVaclav Hapla } 342*c1730cc1SVaclav Hapla PetscCall(DestroySFs(ctx, &sfs)); 3439566063dSJacob Faibussowitsch PetscCall(PetscFree(leafOffsets)); 3449566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 345*c1730cc1SVaclav Hapla if (ctx->viewer) { 346*c1730cc1SVaclav Hapla PetscCall(PetscViewerPopFormat(ctx->viewer)); 347*c1730cc1SVaclav Hapla PetscCall(PetscViewerDestroy(&ctx->viewer)); 348*c1730cc1SVaclav Hapla } 3499566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 350b122ec5aSJacob Faibussowitsch return 0; 351f210f596SVaclav Hapla } 352f210f596SVaclav Hapla 353f210f596SVaclav Hapla /*TEST 354f210f596SVaclav Hapla test: 355f210f596SVaclav Hapla nsize: {{1 3}} 356*c1730cc1SVaclav Hapla args: -compare_to_reference -nsfs {{1 3}} -n {{0 1 5}} -leave_step {{1 3}} -root_mode {{local shared global}} 357*c1730cc1SVaclav Hapla 358*c1730cc1SVaclav Hapla test: 359*c1730cc1SVaclav Hapla suffix: 2 360*c1730cc1SVaclav Hapla nsize: 2 361*c1730cc1SVaclav Hapla args: -irregular {{false true}separate output} -sf_view -nsfs 3 -n 1 -leave_step {{1 3}separate output} -root_mode {{local shared global}separate output} 362f210f596SVaclav Hapla TEST*/ 363