1c4762a1bSJed Brown static char help[] = "Test that shared points on interface of partitions can be rebalanced.\n\n"; 2c4762a1bSJed Brown static char FILENAME[] = "ex31.c"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdmplex.h> 5c4762a1bSJed Brown #include <petscviewerhdf5.h> 630602db0SMatthew G. Knepley #include <petscsf.h> 7c4762a1bSJed Brown 8c4762a1bSJed Brown typedef struct { 9c4762a1bSJed Brown PetscBool parallel; /* Use ParMetis or Metis */ 10c4762a1bSJed Brown PetscBool useInitialGuess; /* Only active when in parallel, uses RefineKway of ParMetis */ 11c4762a1bSJed Brown PetscInt entityDepth; /* depth of the entities to rebalance ( 0 => vertices) */ 12c4762a1bSJed Brown } AppCtx; 13c4762a1bSJed Brown 14c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 15c4762a1bSJed Brown { 16c4762a1bSJed Brown PetscFunctionBegin; 17c4762a1bSJed Brown options->parallel = PETSC_FALSE; 18c4762a1bSJed Brown options->useInitialGuess = PETSC_FALSE; 19c4762a1bSJed Brown options->entityDepth = 0; 2030602db0SMatthew G. Knepley 21d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX"); 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-entity_depth", "Depth of the entities to rebalance (0 => vertices)", FILENAME, options->entityDepth, &options->entityDepth, NULL,0)); 239566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-parallel", "Use ParMetis instead of Metis", FILENAME, options->parallel, &options->parallel, NULL)); 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_initial_guess", "Use RefineKway function of ParMetis", FILENAME, options->useInitialGuess, &options->useInitialGuess, NULL)); 25d0609cedSBarry Smith PetscOptionsEnd(); 26c4762a1bSJed Brown PetscFunctionReturn(0); 27c4762a1bSJed Brown } 28c4762a1bSJed Brown 29c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 30c4762a1bSJed Brown { 31c4762a1bSJed Brown PetscFunctionBegin; 329566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 339566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 349566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 359566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 36c4762a1bSJed Brown PetscFunctionReturn(0); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown 39c4762a1bSJed Brown int main(int argc, char **argv) 40c4762a1bSJed Brown { 41c4762a1bSJed Brown MPI_Comm comm; 42c4762a1bSJed Brown DM dm, dmdist; 43c4762a1bSJed Brown PetscPartitioner part; 44c4762a1bSJed Brown AppCtx user; 45c4762a1bSJed Brown IS is=NULL; 46c4762a1bSJed Brown PetscSection s=NULL, gsection=NULL; 47c4762a1bSJed Brown PetscMPIInt size; 48c4762a1bSJed Brown PetscSF sf; 49c4762a1bSJed Brown PetscInt pStart, pEnd, p, minBefore, maxBefore, minAfter, maxAfter, gSizeBefore, gSizeAfter; 50c4762a1bSJed Brown PetscBool success; 51c4762a1bSJed Brown 529566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 53c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 559566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 569566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &user, &dm)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* partition dm using PETSCPARTITIONERPARMETIS */ 599566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(dm, &part)); 609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part,"p_")); 619566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERPARMETIS)); 629566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 639566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &s)); 649566063dSJacob Faibussowitsch PetscCall(PetscPartitionerDMPlexPartition(part, dm, NULL, s, &is)); 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(dm, 0, NULL, &dmdist)); 67c4762a1bSJed Brown if (dmdist) { 689566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 69c4762a1bSJed Brown dm = dmdist; 70c4762a1bSJed Brown } 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* cleanup */ 739566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&s)); 749566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* We make a PetscSection with a DOF on every mesh entity of depth 77c4762a1bSJed Brown * user.entityDepth, then make a global section and look at its storage size. 78c4762a1bSJed Brown * We do the same thing after the rebalancing and then assert that the size 79c4762a1bSJed Brown * remains the same. We also make sure that the balance has improved at least 80c4762a1bSJed Brown * a little bit compared to the initial decomposition. */ 81c4762a1bSJed Brown 82c4762a1bSJed Brown if (size>1) { 839566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &s)); 849566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(s, 1)); 859566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(s, 0, 1)); 869566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, user.entityDepth, &pStart, &pEnd)); 879566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(s, pStart, pEnd)); 88c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 899566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(s, p, 1)); 909566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(s, p, 0, 1)); 91c4762a1bSJed Brown } 929566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(s)); 939566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(dm, &sf)); 949566063dSJacob Faibussowitsch PetscCall(PetscSectionCreateGlobalSection(s, sf, PETSC_FALSE, PETSC_FALSE, &gsection)); 959566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(gsection, &gSizeBefore)); 96c4762a1bSJed Brown minBefore = gSizeBefore; 97c4762a1bSJed Brown maxBefore = gSizeBefore; 989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &gSizeBefore, 1, MPIU_INT, MPI_SUM, comm)); 999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &minBefore, 1, MPIU_INT, MPI_MIN, comm)); 1009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &maxBefore, 1, MPIU_INT, MPI_MAX, comm)); 1019566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&gsection)); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 1049566063dSJacob Faibussowitsch PetscCall(DMPlexRebalanceSharedPoints(dm, user.entityDepth, user.useInitialGuess, user.parallel, &success)); 105c4762a1bSJed Brown 106c4762a1bSJed Brown if (size>1) { 1079566063dSJacob Faibussowitsch PetscCall(PetscSectionCreateGlobalSection(s, sf, PETSC_FALSE, PETSC_FALSE, &gsection)); 1089566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(gsection, &gSizeAfter)); 109c4762a1bSJed Brown minAfter = gSizeAfter; 110c4762a1bSJed Brown maxAfter = gSizeAfter; 1119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &gSizeAfter, 1, MPIU_INT, MPI_SUM, comm)); 1129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &minAfter, 1, MPIU_INT, MPI_MIN, comm)); 1139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &maxAfter, 1, MPIU_INT, MPI_MAX, comm)); 11408401ef6SPierre Jolivet PetscCheck(gSizeAfter == gSizeBefore,comm, PETSC_ERR_PLIB, "Global section has not the same size before and after."); 115*1dca8a05SBarry Smith PetscCheck(minAfter >= minBefore && maxAfter <= maxBefore && (minAfter > minBefore || maxAfter < maxBefore),comm, PETSC_ERR_PLIB, "DMPlexRebalanceSharedPoints did not improve mesh point balance."); 1169566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&gsection)); 1179566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&s)); 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 1209566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1219566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 122b122ec5aSJacob Faibussowitsch return 0; 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125c4762a1bSJed Brown /*TEST 126c4762a1bSJed Brown 12730602db0SMatthew G. Knepley testset: 12830602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 12930602db0SMatthew G. Knepley 130c4762a1bSJed Brown test: 131c4762a1bSJed Brown # rebalance a mesh 132c4762a1bSJed Brown suffix: 0 133c4762a1bSJed Brown nsize: {{2 3 4}} 134c4762a1bSJed Brown requires: parmetis 13530602db0SMatthew G. Knepley args: -dm_plex_box_faces {{2,3,4 5,4,3 7,11,5}} -entity_depth {{0 1}} -parallel {{FALSE TRUE}} -use_initial_guess FALSE 136c4762a1bSJed Brown 137c4762a1bSJed Brown test: 138c4762a1bSJed Brown # rebalance a mesh but use the initial guess (uses a random algorithm and gives different results on different machines, so just check that it runs). 139c4762a1bSJed Brown suffix: 1 140c4762a1bSJed Brown nsize: {{2 3 4}} 141c4762a1bSJed Brown requires: parmetis 14230602db0SMatthew G. Knepley args: -dm_plex_box_faces {{2,3,4 5,4,3 7,11,5}} -entity_depth {{0 1}} -parallel TRUE -use_initial_guess TRUE 143c4762a1bSJed Brown 144c4762a1bSJed Brown test: 145c4762a1bSJed Brown # no-op in serial 146c4762a1bSJed Brown suffix: 2 147c4762a1bSJed Brown nsize: {{1}} 148c4762a1bSJed Brown requires: parmetis 14930602db0SMatthew G. Knepley args: -dm_plex_box_faces 2,3,4 -entity_depth 0 -parallel FALSE -use_initial_guess FALSE 150c4762a1bSJed Brown 151c4762a1bSJed Brown TEST*/ 152