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