1c4762a1bSJed Brown static char help[] = "Partition a mesh in parallel, perhaps with overlap\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown #include <petscsf.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown /* Sample usage: 7c4762a1bSJed Brown 8c4762a1bSJed Brown Load a file in serial and distribute it on 24 processes: 9c4762a1bSJed Brown 10c4762a1bSJed Brown make -f ./gmakefile test globsearch="dm_impls_plex_tests-ex12_0" EXTRA_OPTIONS="-filename $PETSC_DIR/share/petsc/datafiles/meshes/squaremotor-30.exo -orig_dm_view -dm_view" NP=24 11c4762a1bSJed Brown 12c4762a1bSJed Brown Load a file in serial and distribute it on 24 processes using a custom partitioner: 13c4762a1bSJed Brown 14c4762a1bSJed Brown make -f ./gmakefile test globsearch="dm_impls_plex_tests-ex12_0" EXTRA_OPTIONS="-filename $PETSC_DIR/share/petsc/datafiles/meshes/cylinder.med -petscpartitioner_type simple -orig_dm_view -dm_view" NP=24 15c4762a1bSJed Brown 16c4762a1bSJed Brown Load a file in serial, distribute it, and then redistribute it on 24 processes using two different partitioners: 17c4762a1bSJed Brown 18c4762a1bSJed Brown make -f ./gmakefile test globsearch="dm_impls_plex_tests-ex12_0" EXTRA_OPTIONS="-filename $PETSC_DIR/share/petsc/datafiles/meshes/squaremotor-30.exo -petscpartitioner_type simple -load_balance -lb_petscpartitioner_type parmetis -orig_dm_view -dm_view" NP=24 19c4762a1bSJed Brown 20c4762a1bSJed Brown Load a file in serial, distribute it randomly, refine it in parallel, and then redistribute it on 24 processes using two different partitioners, and view to VTK: 21c4762a1bSJed Brown 22c4762a1bSJed Brown make -f ./gmakefile test globsearch="dm_impls_plex_tests-ex12_0" EXTRA_OPTIONS="-filename $PETSC_DIR/share/petsc/datafiles/meshes/squaremotor-30.exo -petscpartitioner_type shell -petscpartitioner_shell_random -dm_refine 1 -load_balance -lb_petscpartitioner_type parmetis -prelb_dm_view vtk:$PWD/prelb.vtk -dm_view vtk:$PWD/balance.vtk -dm_partition_view" NP=24 23c4762a1bSJed Brown 24c4762a1bSJed Brown */ 25c4762a1bSJed Brown 26c4762a1bSJed Brown enum {STAGE_LOAD, STAGE_DISTRIBUTE, STAGE_REFINE, STAGE_REDISTRIBUTE}; 27c4762a1bSJed Brown 28c4762a1bSJed Brown typedef struct { 29c4762a1bSJed Brown /* Domain and mesh definition */ 30c4762a1bSJed Brown PetscInt overlap; /* The cell overlap to use during partitioning */ 31c4762a1bSJed Brown PetscBool testPartition; /* Use a fixed partitioning for testing */ 32c4762a1bSJed Brown PetscBool testRedundant; /* Use a redundant partitioning for testing */ 33c4762a1bSJed Brown PetscBool loadBalance; /* Load balance via a second distribute step */ 34c4762a1bSJed Brown PetscBool partitionBalance; /* Balance shared point partition */ 35c4762a1bSJed Brown PetscLogStage stages[4]; 36c4762a1bSJed Brown } AppCtx; 37c4762a1bSJed Brown 38c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 39c4762a1bSJed Brown { 40c4762a1bSJed Brown PetscErrorCode ierr; 41c4762a1bSJed Brown 42c4762a1bSJed Brown PetscFunctionBegin; 43c4762a1bSJed Brown options->overlap = 0; 44c4762a1bSJed Brown options->testPartition = PETSC_FALSE; 45c4762a1bSJed Brown options->testRedundant = PETSC_FALSE; 46c4762a1bSJed Brown options->loadBalance = PETSC_FALSE; 47c4762a1bSJed Brown options->partitionBalance = PETSC_FALSE; 48c4762a1bSJed Brown 49c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex12.c", options->overlap, &options->overlap, NULL,0);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex12.c", options->testPartition, &options->testPartition, NULL);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = PetscOptionsBool("-test_redundant", "Use a redundant partition for testing", "ex12.c", options->testRedundant, &options->testRedundant, NULL);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = PetscOptionsBool("-load_balance", "Perform parallel load balancing in a second distribution step", "ex12.c", options->loadBalance, &options->loadBalance, NULL);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = PetscOptionsBool("-partition_balance", "Balance the ownership of shared points", "ex12.c", options->partitionBalance, &options->partitionBalance, NULL);CHKERRQ(ierr); 55*1e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 56c4762a1bSJed Brown 57c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshLoad", &options->stages[STAGE_LOAD]);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]);CHKERRQ(ierr); 59c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshRefine", &options->stages[STAGE_REFINE]);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshRedistribute", &options->stages[STAGE_REDISTRIBUTE]);CHKERRQ(ierr); 61c4762a1bSJed Brown PetscFunctionReturn(0); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 65c4762a1bSJed Brown { 66c4762a1bSJed Brown DM pdm = NULL; 67c4762a1bSJed Brown PetscInt triSizes_n2[2] = {4, 4}; 68c4762a1bSJed Brown PetscInt triPoints_n2[8] = {0, 1, 4, 6, 2, 3, 5, 7}; 69c4762a1bSJed Brown PetscInt triSizes_n3[3] = {3, 2, 3}; 70c4762a1bSJed Brown PetscInt triPoints_n3[8] = {3, 5, 6, 1, 7, 0, 2, 4}; 71c4762a1bSJed Brown PetscInt triSizes_n4[4] = {2, 2, 2, 2}; 72c4762a1bSJed Brown PetscInt triPoints_n4[8] = {0, 7, 1, 5, 2, 3, 4, 6}; 73c4762a1bSJed Brown PetscInt triSizes_n8[8] = {1, 1, 1, 1, 1, 1, 1, 1}; 74c4762a1bSJed Brown PetscInt triPoints_n8[8] = {0, 1, 2, 3, 4, 5, 6, 7}; 75c4762a1bSJed Brown PetscInt quadSizes[2] = {2, 2}; 76c4762a1bSJed Brown PetscInt quadPoints[4] = {2, 3, 0, 1}; 77c4762a1bSJed Brown PetscInt overlap = user->overlap >= 0 ? user->overlap : 0; 7830602db0SMatthew G. Knepley PetscInt dim; 7930602db0SMatthew G. Knepley PetscBool simplex; 80c4762a1bSJed Brown PetscMPIInt rank, size; 81c4762a1bSJed Brown PetscErrorCode ierr; 82c4762a1bSJed Brown 83c4762a1bSJed Brown PetscFunctionBegin; 84ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 85ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 86c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_LOAD]);CHKERRQ(ierr); 8730602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 8830602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 8930602db0SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-orig_dm_view");CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 9230602db0SMatthew G. Knepley ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr); 9330602db0SMatthew G. Knepley ierr = DMPlexIsSimplex(*dm, &simplex);CHKERRQ(ierr); 94c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]);CHKERRQ(ierr); 95c4762a1bSJed Brown if (!user->testRedundant) { 96c4762a1bSJed Brown PetscPartitioner part; 97c4762a1bSJed Brown 98c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = DMPlexSetPartitionBalance(*dm, user->partitionBalance);CHKERRQ(ierr); 101c4762a1bSJed Brown if (user->testPartition) { 102c4762a1bSJed Brown const PetscInt *sizes = NULL; 103c4762a1bSJed Brown const PetscInt *points = NULL; 104c4762a1bSJed Brown 105c4762a1bSJed Brown if (!rank) { 10630602db0SMatthew G. Knepley if (dim == 2 && simplex && size == 2) { 107c4762a1bSJed Brown sizes = triSizes_n2; points = triPoints_n2; 10830602db0SMatthew G. Knepley } else if (dim == 2 && simplex && size == 3) { 109c4762a1bSJed Brown sizes = triSizes_n3; points = triPoints_n3; 11030602db0SMatthew G. Knepley } else if (dim == 2 && simplex && size == 4) { 111c4762a1bSJed Brown sizes = triSizes_n4; points = triPoints_n4; 11230602db0SMatthew G. Knepley } else if (dim == 2 && simplex && size == 8) { 113c4762a1bSJed Brown sizes = triSizes_n8; points = triPoints_n8; 11430602db0SMatthew G. Knepley } else if (dim == 2 && !simplex && size == 2) { 115c4762a1bSJed Brown sizes = quadSizes; points = quadPoints; 116c4762a1bSJed Brown } 117c4762a1bSJed Brown } 118c4762a1bSJed Brown ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = PetscPartitionerShellSetPartition(part, size, sizes, points);CHKERRQ(ierr); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, overlap, NULL, &pdm);CHKERRQ(ierr); 122c4762a1bSJed Brown } else { 123c4762a1bSJed Brown PetscSF sf; 124c4762a1bSJed Brown 125c4762a1bSJed Brown ierr = DMPlexGetRedundantDM(*dm, &sf, &pdm);CHKERRQ(ierr); 126c4762a1bSJed Brown if (sf) { 127c4762a1bSJed Brown DM test; 128c4762a1bSJed Brown 129c4762a1bSJed Brown ierr = DMPlexCreate(comm,&test);CHKERRQ(ierr); 130c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)test, "Test SF-migrated Redundant Mesh");CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = DMPlexMigrate(*dm, sf, test);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = DMViewFromOptions(test, NULL, "-redundant_migrated_dm_view");CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = DMDestroy(&test);CHKERRQ(ierr); 134c4762a1bSJed Brown } 135c4762a1bSJed Brown ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown if (pdm) { 138c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 139c4762a1bSJed Brown *dm = pdm; 140c4762a1bSJed Brown } 141c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 143c4762a1bSJed Brown if (user->loadBalance) { 144c4762a1bSJed Brown PetscPartitioner part; 145c4762a1bSJed Brown 146c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-prelb_dm_view");CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = DMPlexSetOptionsPrefix(*dm, "lb_");CHKERRQ(ierr); 148c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_REDISTRIBUTE]);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject) part, "lb_");CHKERRQ(ierr); 151c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 152c4762a1bSJed Brown if (user->testPartition) { 153c4762a1bSJed Brown PetscInt reSizes_n2[2] = {2, 2}; 154c4762a1bSJed Brown PetscInt rePoints_n2[4] = {2, 3, 0, 1}; 155c4762a1bSJed Brown if (rank) {rePoints_n2[0] = 1; rePoints_n2[1] = 2, rePoints_n2[2] = 0, rePoints_n2[3] = 3;} 156c4762a1bSJed Brown 157c4762a1bSJed Brown ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = PetscPartitionerShellSetPartition(part, size, reSizes_n2, rePoints_n2);CHKERRQ(ierr); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown ierr = DMPlexSetPartitionBalance(*dm, user->partitionBalance);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, overlap, NULL, &pdm);CHKERRQ(ierr); 162c4762a1bSJed Brown if (pdm) { 163c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 164c4762a1bSJed Brown *dm = pdm; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_REFINE]);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 171c4762a1bSJed Brown PetscFunctionReturn(0); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown 174c4762a1bSJed Brown int main(int argc, char **argv) 175c4762a1bSJed Brown { 176c4762a1bSJed Brown DM dm; 177c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 178c4762a1bSJed Brown PetscErrorCode ierr; 179c4762a1bSJed Brown 180c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 181c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 182c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 183c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 184c4762a1bSJed Brown ierr = PetscFinalize(); 185c4762a1bSJed Brown return ierr; 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 188c4762a1bSJed Brown /*TEST 189c4762a1bSJed Brown # Parallel, no overlap tests 0-2 190c4762a1bSJed Brown test: 191c4762a1bSJed Brown suffix: 0 192c4762a1bSJed Brown requires: triangle 19330602db0SMatthew G. Knepley args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex 194c4762a1bSJed Brown test: 195c4762a1bSJed Brown suffix: 1 196c4762a1bSJed Brown requires: triangle 197c4762a1bSJed Brown nsize: 3 19830602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail 199c4762a1bSJed Brown test: 200c4762a1bSJed Brown suffix: 2 201c4762a1bSJed Brown requires: triangle 202c4762a1bSJed Brown nsize: 8 20330602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail 204c4762a1bSJed Brown # Parallel, level-1 overlap tests 3-4 205c4762a1bSJed Brown test: 206c4762a1bSJed Brown suffix: 3 207c4762a1bSJed Brown requires: triangle 208c4762a1bSJed Brown nsize: 3 20930602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail 210c4762a1bSJed Brown test: 211c4762a1bSJed Brown suffix: 4 212c4762a1bSJed Brown requires: triangle 213c4762a1bSJed Brown nsize: 8 21430602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail 215c4762a1bSJed Brown # Parallel, level-2 overlap test 5 216c4762a1bSJed Brown test: 217c4762a1bSJed Brown suffix: 5 218c4762a1bSJed Brown requires: triangle 219c4762a1bSJed Brown nsize: 8 22030602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail 221c4762a1bSJed Brown # Parallel load balancing, test 6-7 222c4762a1bSJed Brown test: 223c4762a1bSJed Brown suffix: 6 224c4762a1bSJed Brown requires: triangle 225c4762a1bSJed Brown nsize: 2 22630602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail 227c4762a1bSJed Brown test: 228c4762a1bSJed Brown suffix: 7 229c4762a1bSJed Brown requires: triangle 230c4762a1bSJed Brown nsize: 2 23130602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail 232c4762a1bSJed Brown # Parallel redundant copying, test 8 233c4762a1bSJed Brown test: 234c4762a1bSJed Brown suffix: 8 235c4762a1bSJed Brown requires: triangle 236c4762a1bSJed Brown nsize: 2 23730602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_redundant -redundant_migrated_dm_view ascii::ascii_info_detail -dm_view ascii::ascii_info_detail 238c4762a1bSJed Brown test: 239c4762a1bSJed Brown suffix: lb_0 240c4762a1bSJed Brown requires: parmetis 241c4762a1bSJed Brown nsize: 4 24230602db0SMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_simplex 0 -dm_plex_box_faces 4,4 -petscpartitioner_type shell -petscpartitioner_shell_random -lb_petscpartitioner_type parmetis -load_balance -lb_petscpartitioner_view -prelb_dm_view ::load_balance -dm_view ::load_balance 243c4762a1bSJed Brown 244c4762a1bSJed Brown # Same tests as above, but with balancing of the shared point partition 245c4762a1bSJed Brown test: 246c4762a1bSJed Brown suffix: 9 247c4762a1bSJed Brown requires: triangle 24830602db0SMatthew G. Knepley args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex -partition_balance 249c4762a1bSJed Brown test: 250c4762a1bSJed Brown suffix: 10 251c4762a1bSJed Brown requires: triangle 252c4762a1bSJed Brown nsize: 3 25330602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance 254c4762a1bSJed Brown test: 255c4762a1bSJed Brown suffix: 11 256c4762a1bSJed Brown requires: triangle 257c4762a1bSJed Brown nsize: 8 25830602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance 259c4762a1bSJed Brown # Parallel, level-1 overlap tests 3-4 260c4762a1bSJed Brown test: 261c4762a1bSJed Brown suffix: 12 262c4762a1bSJed Brown requires: triangle 263c4762a1bSJed Brown nsize: 3 26430602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance 265c4762a1bSJed Brown test: 266c4762a1bSJed Brown suffix: 13 267c4762a1bSJed Brown requires: triangle 268c4762a1bSJed Brown nsize: 8 26930602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance 270c4762a1bSJed Brown # Parallel, level-2 overlap test 5 271c4762a1bSJed Brown test: 272c4762a1bSJed Brown suffix: 14 273c4762a1bSJed Brown requires: triangle 274c4762a1bSJed Brown nsize: 8 27530602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail -partition_balance 276c4762a1bSJed Brown # Parallel load balancing, test 6-7 277c4762a1bSJed Brown test: 278c4762a1bSJed Brown suffix: 15 279c4762a1bSJed Brown requires: triangle 280c4762a1bSJed Brown nsize: 2 28130602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance 282c4762a1bSJed Brown test: 283c4762a1bSJed Brown suffix: 16 284c4762a1bSJed Brown requires: triangle 285c4762a1bSJed Brown nsize: 2 28630602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail -partition_balance 287c4762a1bSJed Brown # Parallel redundant copying, test 8 288c4762a1bSJed Brown test: 289c4762a1bSJed Brown suffix: 17 290c4762a1bSJed Brown requires: triangle 291c4762a1bSJed Brown nsize: 2 29230602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_redundant -dm_view ascii::ascii_info_detail -partition_balance 293c4762a1bSJed Brown test: 294c4762a1bSJed Brown suffix: lb_1 295c4762a1bSJed Brown requires: parmetis 296c4762a1bSJed Brown nsize: 4 29730602db0SMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_simplex 0 -dm_plex_box_faces 4,4 -petscpartitioner_type shell -petscpartitioner_shell_random -lb_petscpartitioner_type parmetis -load_balance -lb_petscpartitioner_view -partition_balance -prelb_dm_view ::load_balance -dm_view ::load_balance 298c4762a1bSJed Brown TEST*/ 299