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 PetscErrorCode ierr; 17c4762a1bSJed Brown 18c4762a1bSJed Brown PetscFunctionBegin; 19c4762a1bSJed Brown options->parallel = PETSC_FALSE; 20c4762a1bSJed Brown options->useInitialGuess = PETSC_FALSE; 21c4762a1bSJed Brown options->entityDepth = 0; 2230602db0SMatthew G. Knepley 23c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr); 24c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-entity_depth", "Depth of the entities to rebalance (0 => vertices)", FILENAME, options->entityDepth, &options->entityDepth, NULL,0);CHKERRQ(ierr); 25c4762a1bSJed Brown ierr = PetscOptionsBool("-parallel", "Use ParMetis instead of Metis", FILENAME, options->parallel, &options->parallel, NULL);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = PetscOptionsBool("-use_initial_guess", "Use RefineKway function of ParMetis", FILENAME, options->useInitialGuess, &options->useInitialGuess, NULL);CHKERRQ(ierr); 27*1e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 28c4762a1bSJed Brown PetscFunctionReturn(0); 29c4762a1bSJed Brown } 30c4762a1bSJed Brown 31c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 32c4762a1bSJed Brown { 33c4762a1bSJed Brown PetscErrorCode ierr; 34c4762a1bSJed Brown 35c4762a1bSJed Brown PetscFunctionBegin; 3630602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 3730602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 3830602db0SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 3930602db0SMatthew G. Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 40c4762a1bSJed Brown PetscFunctionReturn(0); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43c4762a1bSJed Brown int main(int argc, char **argv) 44c4762a1bSJed Brown { 45c4762a1bSJed Brown MPI_Comm comm; 46c4762a1bSJed Brown DM dm, dmdist; 47c4762a1bSJed Brown PetscPartitioner part; 48c4762a1bSJed Brown AppCtx user; 49c4762a1bSJed Brown IS is=NULL; 50c4762a1bSJed Brown PetscSection s=NULL, gsection=NULL; 51c4762a1bSJed Brown PetscErrorCode ierr; 52c4762a1bSJed Brown PetscMPIInt size; 53c4762a1bSJed Brown PetscSF sf; 54c4762a1bSJed Brown PetscInt pStart, pEnd, p, minBefore, maxBefore, minAfter, maxAfter, gSizeBefore, gSizeAfter; 55c4762a1bSJed Brown PetscBool success; 56c4762a1bSJed Brown 57c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 58c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 59ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 60c4762a1bSJed Brown ierr = ProcessOptions(comm, &user);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = CreateMesh(comm, &user, &dm);CHKERRQ(ierr); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* partition dm using PETSCPARTITIONERPARMETIS */ 64c4762a1bSJed Brown ierr = DMPlexGetPartitioner(dm, &part);CHKERRQ(ierr); 65c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)part,"p_");CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = PetscPartitionerSetType(part, PETSCPARTITIONERPARMETIS);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = PetscSectionCreate(comm, &s);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = PetscPartitionerDMPlexPartition(part, dm, NULL, s, &is);CHKERRQ(ierr); 70c4762a1bSJed Brown 71c4762a1bSJed Brown ierr = DMPlexDistribute(dm, 0, NULL, &dmdist);CHKERRQ(ierr); 72c4762a1bSJed Brown if (dmdist) { 73c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 74c4762a1bSJed Brown dm = dmdist; 75c4762a1bSJed Brown } 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* cleanup */ 78c4762a1bSJed Brown ierr = PetscSectionDestroy(&s);CHKERRQ(ierr); 79c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* We make a PetscSection with a DOF on every mesh entity of depth 82c4762a1bSJed Brown * user.entityDepth, then make a global section and look at its storage size. 83c4762a1bSJed Brown * We do the same thing after the rebalancing and then assert that the size 84c4762a1bSJed Brown * remains the same. We also make sure that the balance has improved at least 85c4762a1bSJed Brown * a little bit compared to the initial decomposition. */ 86c4762a1bSJed Brown 87c4762a1bSJed Brown if (size>1) { 88c4762a1bSJed Brown ierr = PetscSectionCreate(comm, &s);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = PetscSectionSetNumFields(s, 1);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = PetscSectionSetFieldComponents(s, 0, 1);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = DMPlexGetDepthStratum(dm, user.entityDepth, &pStart, &pEnd);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = PetscSectionSetChart(s, pStart, pEnd);CHKERRQ(ierr); 93c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 94c4762a1bSJed Brown ierr = PetscSectionSetDof(s, p, 1);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = PetscSectionSetFieldDof(s, p, 0, 1);CHKERRQ(ierr); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown ierr = PetscSectionSetUp(s);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = PetscSectionCreateGlobalSection(s, sf, PETSC_FALSE, PETSC_FALSE, &gsection);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = PetscSectionGetStorageSize(gsection, &gSizeBefore);CHKERRQ(ierr); 101c4762a1bSJed Brown minBefore = gSizeBefore; 102c4762a1bSJed Brown maxBefore = gSizeBefore; 103ffc4695bSBarry Smith ierr = MPI_Allreduce(MPI_IN_PLACE, &gSizeBefore, 1, MPIU_INT, MPI_SUM, comm);CHKERRMPI(ierr); 104ffc4695bSBarry Smith ierr = MPI_Allreduce(MPI_IN_PLACE, &minBefore, 1, MPIU_INT, MPI_MIN, comm);CHKERRMPI(ierr); 105ffc4695bSBarry Smith ierr = MPI_Allreduce(MPI_IN_PLACE, &maxBefore, 1, MPIU_INT, MPI_MAX, comm);CHKERRMPI(ierr); 106c4762a1bSJed Brown ierr = PetscSectionDestroy(&gsection);CHKERRQ(ierr); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109c4762a1bSJed Brown ierr = DMPlexRebalanceSharedPoints(dm, user.entityDepth, user.useInitialGuess, user.parallel, &success);CHKERRQ(ierr); 110c4762a1bSJed Brown 111c4762a1bSJed Brown if (size>1) { 112c4762a1bSJed Brown ierr = PetscSectionCreateGlobalSection(s, sf, PETSC_FALSE, PETSC_FALSE, &gsection);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = PetscSectionGetStorageSize(gsection, &gSizeAfter);CHKERRQ(ierr); 114c4762a1bSJed Brown minAfter = gSizeAfter; 115c4762a1bSJed Brown maxAfter = gSizeAfter; 116ffc4695bSBarry Smith ierr = MPI_Allreduce(MPI_IN_PLACE, &gSizeAfter, 1, MPIU_INT, MPI_SUM, comm);CHKERRMPI(ierr); 117ffc4695bSBarry Smith ierr = MPI_Allreduce(MPI_IN_PLACE, &minAfter, 1, MPIU_INT, MPI_MIN, comm);CHKERRMPI(ierr); 118ffc4695bSBarry Smith ierr = MPI_Allreduce(MPI_IN_PLACE, &maxAfter, 1, MPIU_INT, MPI_MAX, comm);CHKERRMPI(ierr); 119c4762a1bSJed Brown if (gSizeAfter != gSizeBefore) SETERRQ(comm, PETSC_ERR_PLIB, "Global section has not the same size before and after."); 120c4762a1bSJed Brown if (!(minAfter >= minBefore && maxAfter <= maxBefore && (minAfter > minBefore || maxAfter < maxBefore))) SETERRQ(comm, PETSC_ERR_PLIB, "DMPlexRebalanceSharedPoints did not improve mesh point balance."); 121c4762a1bSJed Brown ierr = PetscSectionDestroy(&gsection);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = PetscSectionDestroy(&s);CHKERRQ(ierr); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = PetscFinalize(); 127c4762a1bSJed Brown return ierr; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown /*TEST 131c4762a1bSJed Brown 13230602db0SMatthew G. Knepley testset: 13330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 13430602db0SMatthew G. Knepley 135c4762a1bSJed Brown test: 136c4762a1bSJed Brown # rebalance a mesh 137c4762a1bSJed Brown suffix: 0 138c4762a1bSJed Brown nsize: {{2 3 4}} 139c4762a1bSJed Brown requires: parmetis 14030602db0SMatthew 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 141c4762a1bSJed Brown 142c4762a1bSJed Brown test: 143c4762a1bSJed 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). 144c4762a1bSJed Brown suffix: 1 145c4762a1bSJed Brown nsize: {{2 3 4}} 146c4762a1bSJed Brown requires: parmetis 14730602db0SMatthew 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 148c4762a1bSJed Brown 149c4762a1bSJed Brown test: 150c4762a1bSJed Brown # no-op in serial 151c4762a1bSJed Brown suffix: 2 152c4762a1bSJed Brown nsize: {{1}} 153c4762a1bSJed Brown requires: parmetis 15430602db0SMatthew G. Knepley args: -dm_plex_box_faces 2,3,4 -entity_depth 0 -parallel FALSE -use_initial_guess FALSE 155c4762a1bSJed Brown 156c4762a1bSJed Brown TEST*/ 157