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); 245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-entity_depth", "Depth of the entities to rebalance (0 => vertices)", FILENAME, options->entityDepth, &options->entityDepth, NULL,0)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-parallel", "Use ParMetis instead of Metis", FILENAME, options->parallel, &options->parallel, NULL)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-use_initial_guess", "Use RefineKway function of ParMetis", FILENAME, options->useInitialGuess, &options->useInitialGuess, NULL)); 271e1ea65dSPierre 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 PetscFunctionBegin; 345f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 38c4762a1bSJed Brown PetscFunctionReturn(0); 39c4762a1bSJed Brown } 40c4762a1bSJed Brown 41c4762a1bSJed Brown int main(int argc, char **argv) 42c4762a1bSJed Brown { 43c4762a1bSJed Brown MPI_Comm comm; 44c4762a1bSJed Brown DM dm, dmdist; 45c4762a1bSJed Brown PetscPartitioner part; 46c4762a1bSJed Brown AppCtx user; 47c4762a1bSJed Brown IS is=NULL; 48c4762a1bSJed Brown PetscSection s=NULL, gsection=NULL; 49c4762a1bSJed Brown PetscMPIInt size; 50c4762a1bSJed Brown PetscSF sf; 51c4762a1bSJed Brown PetscInt pStart, pEnd, p, minBefore, maxBefore, minAfter, maxAfter, gSizeBefore, gSizeAfter; 52c4762a1bSJed Brown PetscBool success; 53c4762a1bSJed Brown 54*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 55c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 565f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm, &size)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(comm, &user)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(comm, &user, &dm)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* partition dm using PETSCPARTITIONERPARMETIS */ 615f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetPartitioner(dm, &part)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)part,"p_")); 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPartitionerSetType(part, PETSCPARTITIONERPARMETIS)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPartitionerSetFromOptions(part)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionCreate(comm, &s)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPartitionerDMPlexPartition(part, dm, NULL, s, &is)); 67c4762a1bSJed Brown 685f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistribute(dm, 0, NULL, &dmdist)); 69c4762a1bSJed Brown if (dmdist) { 705f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 71c4762a1bSJed Brown dm = dmdist; 72c4762a1bSJed Brown } 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* cleanup */ 755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(&s)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* We make a PetscSection with a DOF on every mesh entity of depth 79c4762a1bSJed Brown * user.entityDepth, then make a global section and look at its storage size. 80c4762a1bSJed Brown * We do the same thing after the rebalancing and then assert that the size 81c4762a1bSJed Brown * remains the same. We also make sure that the balance has improved at least 82c4762a1bSJed Brown * a little bit compared to the initial decomposition. */ 83c4762a1bSJed Brown 84c4762a1bSJed Brown if (size>1) { 855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionCreate(comm, &s)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetNumFields(s, 1)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldComponents(s, 0, 1)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(dm, user.entityDepth, &pStart, &pEnd)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetChart(s, pStart, pEnd)); 90c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetDof(s, p, 1)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldDof(s, p, 0, 1)); 93c4762a1bSJed Brown } 945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetUp(s)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetPointSF(dm, &sf)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionCreateGlobalSection(s, sf, PETSC_FALSE, PETSC_FALSE, &gsection)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetStorageSize(gsection, &gSizeBefore)); 98c4762a1bSJed Brown minBefore = gSizeBefore; 99c4762a1bSJed Brown maxBefore = gSizeBefore; 1005f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(MPI_IN_PLACE, &gSizeBefore, 1, MPIU_INT, MPI_SUM, comm)); 1015f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(MPI_IN_PLACE, &minBefore, 1, MPIU_INT, MPI_MIN, comm)); 1025f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(MPI_IN_PLACE, &maxBefore, 1, MPIU_INT, MPI_MAX, comm)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(&gsection)); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 1065f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRebalanceSharedPoints(dm, user.entityDepth, user.useInitialGuess, user.parallel, &success)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown if (size>1) { 1095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionCreateGlobalSection(s, sf, PETSC_FALSE, PETSC_FALSE, &gsection)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetStorageSize(gsection, &gSizeAfter)); 111c4762a1bSJed Brown minAfter = gSizeAfter; 112c4762a1bSJed Brown maxAfter = gSizeAfter; 1135f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(MPI_IN_PLACE, &gSizeAfter, 1, MPIU_INT, MPI_SUM, comm)); 1145f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(MPI_IN_PLACE, &minAfter, 1, MPIU_INT, MPI_MIN, comm)); 1155f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(MPI_IN_PLACE, &maxAfter, 1, MPIU_INT, MPI_MAX, comm)); 1162c71b3e2SJacob Faibussowitsch PetscCheckFalse(gSizeAfter != gSizeBefore,comm, PETSC_ERR_PLIB, "Global section has not the same size before and after."); 1172c71b3e2SJacob Faibussowitsch PetscCheckFalse(!(minAfter >= minBefore && maxAfter <= maxBefore && (minAfter > minBefore || maxAfter < maxBefore)),comm, PETSC_ERR_PLIB, "DMPlexRebalanceSharedPoints did not improve mesh point balance."); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(&gsection)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(&s)); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 1225f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 123*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 124*b122ec5aSJacob Faibussowitsch return 0; 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown /*TEST 128c4762a1bSJed Brown 12930602db0SMatthew G. Knepley testset: 13030602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 13130602db0SMatthew G. Knepley 132c4762a1bSJed Brown test: 133c4762a1bSJed Brown # rebalance a mesh 134c4762a1bSJed Brown suffix: 0 135c4762a1bSJed Brown nsize: {{2 3 4}} 136c4762a1bSJed Brown requires: parmetis 13730602db0SMatthew 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 138c4762a1bSJed Brown 139c4762a1bSJed Brown test: 140c4762a1bSJed 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). 141c4762a1bSJed Brown suffix: 1 142c4762a1bSJed Brown nsize: {{2 3 4}} 143c4762a1bSJed Brown requires: parmetis 14430602db0SMatthew 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 145c4762a1bSJed Brown 146c4762a1bSJed Brown test: 147c4762a1bSJed Brown # no-op in serial 148c4762a1bSJed Brown suffix: 2 149c4762a1bSJed Brown nsize: {{1}} 150c4762a1bSJed Brown requires: parmetis 15130602db0SMatthew G. Knepley args: -dm_plex_box_faces 2,3,4 -entity_depth 0 -parallel FALSE -use_initial_guess FALSE 152c4762a1bSJed Brown 153c4762a1bSJed Brown TEST*/ 154