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