1*c4762a1bSJed Brown static char help[] = "Partition a mesh in parallel, perhaps with overlap\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscdmplex.h> 4*c4762a1bSJed Brown #include <petscsf.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown /* Sample usage: 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown Load a file in serial and distribute it on 24 processes: 9*c4762a1bSJed Brown 10*c4762a1bSJed 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 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown Load a file in serial and distribute it on 24 processes using a custom partitioner: 13*c4762a1bSJed Brown 14*c4762a1bSJed 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 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown Load a file in serial, distribute it, and then redistribute it on 24 processes using two different partitioners: 17*c4762a1bSJed Brown 18*c4762a1bSJed 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 19*c4762a1bSJed Brown 20*c4762a1bSJed 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: 21*c4762a1bSJed Brown 22*c4762a1bSJed 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 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown */ 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown enum {STAGE_LOAD, STAGE_DISTRIBUTE, STAGE_REFINE, STAGE_REDISTRIBUTE}; 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown typedef struct { 29*c4762a1bSJed Brown /* Domain and mesh definition */ 30*c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 31*c4762a1bSJed Brown PetscBool cellSimplex; /* Use simplices or hexes */ 32*c4762a1bSJed Brown PetscInt cells[3]; /* The initial domain division */ 33*c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */ 34*c4762a1bSJed Brown PetscInt overlap; /* The cell overlap to use during partitioning */ 35*c4762a1bSJed Brown PetscBool testPartition; /* Use a fixed partitioning for testing */ 36*c4762a1bSJed Brown PetscBool testRedundant; /* Use a redundant partitioning for testing */ 37*c4762a1bSJed Brown PetscBool loadBalance; /* Load balance via a second distribute step */ 38*c4762a1bSJed Brown PetscBool partitionBalance; /* Balance shared point partition */ 39*c4762a1bSJed Brown PetscLogStage stages[4]; 40*c4762a1bSJed Brown } AppCtx; 41*c4762a1bSJed Brown 42*c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 43*c4762a1bSJed Brown { 44*c4762a1bSJed Brown PetscInt n; 45*c4762a1bSJed Brown PetscErrorCode ierr; 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown PetscFunctionBegin; 48*c4762a1bSJed Brown options->dim = 2; 49*c4762a1bSJed Brown options->cellSimplex = PETSC_TRUE; 50*c4762a1bSJed Brown options->cells[0] = 2; 51*c4762a1bSJed Brown options->cells[1] = 2; 52*c4762a1bSJed Brown options->cells[2] = 2; 53*c4762a1bSJed Brown options->filename[0] = '\0'; 54*c4762a1bSJed Brown options->overlap = 0; 55*c4762a1bSJed Brown options->testPartition = PETSC_FALSE; 56*c4762a1bSJed Brown options->testRedundant = PETSC_FALSE; 57*c4762a1bSJed Brown options->loadBalance = PETSC_FALSE; 58*c4762a1bSJed Brown options->partitionBalance = PETSC_FALSE; 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex12.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex12.c", options->cellSimplex, &options->cellSimplex, NULL);CHKERRQ(ierr); 63*c4762a1bSJed Brown n = 3; 64*c4762a1bSJed Brown ierr = PetscOptionsIntArray("-cells", "The initial mesh division", "ex12.c", options->cells, &n, NULL);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = PetscOptionsString("-filename", "The mesh file", "ex12.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex12.c", options->overlap, &options->overlap, NULL,0);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex12.c", options->testPartition, &options->testPartition, NULL);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = PetscOptionsBool("-test_redundant", "Use a redundant partition for testing", "ex12.c", options->testRedundant, &options->testRedundant, NULL);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = PetscOptionsBool("-load_balance", "Perform parallel load balancing in a second distribution step", "ex12.c", options->loadBalance, &options->loadBalance, NULL);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = PetscOptionsBool("-partition_balance", "Balance the ownership of shared points", "ex12.c", options->partitionBalance, &options->partitionBalance, NULL);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = PetscOptionsEnd(); 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshLoad", &options->stages[STAGE_LOAD]);CHKERRQ(ierr); 74*c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE]);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshRefine", &options->stages[STAGE_REFINE]);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = PetscLogStageRegister("MeshRedistribute", &options->stages[STAGE_REDISTRIBUTE]);CHKERRQ(ierr); 77*c4762a1bSJed Brown PetscFunctionReturn(0); 78*c4762a1bSJed Brown } 79*c4762a1bSJed Brown 80*c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 81*c4762a1bSJed Brown { 82*c4762a1bSJed Brown DM pdm = NULL; 83*c4762a1bSJed Brown PetscInt dim = user->dim; 84*c4762a1bSJed Brown PetscBool cellSimplex = user->cellSimplex; 85*c4762a1bSJed Brown const char *filename = user->filename; 86*c4762a1bSJed Brown PetscInt triSizes_n2[2] = {4, 4}; 87*c4762a1bSJed Brown PetscInt triPoints_n2[8] = {0, 1, 4, 6, 2, 3, 5, 7}; 88*c4762a1bSJed Brown PetscInt triSizes_n3[3] = {3, 2, 3}; 89*c4762a1bSJed Brown PetscInt triPoints_n3[8] = {3, 5, 6, 1, 7, 0, 2, 4}; 90*c4762a1bSJed Brown PetscInt triSizes_n4[4] = {2, 2, 2, 2}; 91*c4762a1bSJed Brown PetscInt triPoints_n4[8] = {0, 7, 1, 5, 2, 3, 4, 6}; 92*c4762a1bSJed Brown PetscInt triSizes_n8[8] = {1, 1, 1, 1, 1, 1, 1, 1}; 93*c4762a1bSJed Brown PetscInt triPoints_n8[8] = {0, 1, 2, 3, 4, 5, 6, 7}; 94*c4762a1bSJed Brown PetscInt quadSizes[2] = {2, 2}; 95*c4762a1bSJed Brown PetscInt quadPoints[4] = {2, 3, 0, 1}; 96*c4762a1bSJed Brown PetscInt overlap = user->overlap >= 0 ? user->overlap : 0; 97*c4762a1bSJed Brown size_t len; 98*c4762a1bSJed Brown PetscMPIInt rank, size; 99*c4762a1bSJed Brown PetscErrorCode ierr; 100*c4762a1bSJed Brown 101*c4762a1bSJed Brown PetscFunctionBegin; 102*c4762a1bSJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = PetscStrlen(filename, &len);CHKERRQ(ierr); 105*c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_LOAD]);CHKERRQ(ierr); 106*c4762a1bSJed Brown if (len) {ierr = DMPlexCreateFromFile(comm, filename, PETSC_TRUE, dm);CHKERRQ(ierr);} 107*c4762a1bSJed Brown else {ierr = DMPlexCreateBoxMesh(comm, dim, cellSimplex, user->cells, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);} 108*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-orig_dm_view");CHKERRQ(ierr); 109*c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 110*c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_DISTRIBUTE]);CHKERRQ(ierr); 111*c4762a1bSJed Brown if (!user->testRedundant) { 112*c4762a1bSJed Brown PetscPartitioner part; 113*c4762a1bSJed Brown 114*c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 116*c4762a1bSJed Brown ierr = DMPlexSetPartitionBalance(*dm, user->partitionBalance);CHKERRQ(ierr); 117*c4762a1bSJed Brown if (user->testPartition) { 118*c4762a1bSJed Brown const PetscInt *sizes = NULL; 119*c4762a1bSJed Brown const PetscInt *points = NULL; 120*c4762a1bSJed Brown 121*c4762a1bSJed Brown if (!rank) { 122*c4762a1bSJed Brown if (dim == 2 && cellSimplex && size == 2) { 123*c4762a1bSJed Brown sizes = triSizes_n2; points = triPoints_n2; 124*c4762a1bSJed Brown } else if (dim == 2 && cellSimplex && size == 3) { 125*c4762a1bSJed Brown sizes = triSizes_n3; points = triPoints_n3; 126*c4762a1bSJed Brown } else if (dim == 2 && cellSimplex && size == 4) { 127*c4762a1bSJed Brown sizes = triSizes_n4; points = triPoints_n4; 128*c4762a1bSJed Brown } else if (dim == 2 && cellSimplex && size == 8) { 129*c4762a1bSJed Brown sizes = triSizes_n8; points = triPoints_n8; 130*c4762a1bSJed Brown } else if (dim == 2 && !cellSimplex && size == 2) { 131*c4762a1bSJed Brown sizes = quadSizes; points = quadPoints; 132*c4762a1bSJed Brown } 133*c4762a1bSJed Brown } 134*c4762a1bSJed Brown ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr); 135*c4762a1bSJed Brown ierr = PetscPartitionerShellSetPartition(part, size, sizes, points);CHKERRQ(ierr); 136*c4762a1bSJed Brown } 137*c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, overlap, NULL, &pdm);CHKERRQ(ierr); 138*c4762a1bSJed Brown } else { 139*c4762a1bSJed Brown PetscSF sf; 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown ierr = DMPlexGetRedundantDM(*dm, &sf, &pdm);CHKERRQ(ierr); 142*c4762a1bSJed Brown if (sf) { 143*c4762a1bSJed Brown DM test; 144*c4762a1bSJed Brown 145*c4762a1bSJed Brown ierr = DMPlexCreate(comm,&test);CHKERRQ(ierr); 146*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)test, "Test SF-migrated Redundant Mesh");CHKERRQ(ierr); 147*c4762a1bSJed Brown ierr = DMPlexMigrate(*dm, sf, test);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = DMViewFromOptions(test, NULL, "-redundant_migrated_dm_view");CHKERRQ(ierr); 149*c4762a1bSJed Brown ierr = DMDestroy(&test);CHKERRQ(ierr); 150*c4762a1bSJed Brown } 151*c4762a1bSJed Brown ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 152*c4762a1bSJed Brown } 153*c4762a1bSJed Brown if (pdm) { 154*c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 155*c4762a1bSJed Brown *dm = pdm; 156*c4762a1bSJed Brown } 157*c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 158*c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 159*c4762a1bSJed Brown if (user->loadBalance) { 160*c4762a1bSJed Brown PetscPartitioner part; 161*c4762a1bSJed Brown 162*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-prelb_dm_view");CHKERRQ(ierr); 163*c4762a1bSJed Brown ierr = DMPlexSetOptionsPrefix(*dm, "lb_");CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_REDISTRIBUTE]);CHKERRQ(ierr); 165*c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject) part, "lb_");CHKERRQ(ierr); 167*c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 168*c4762a1bSJed Brown if (user->testPartition) { 169*c4762a1bSJed Brown PetscInt reSizes_n2[2] = {2, 2}; 170*c4762a1bSJed Brown PetscInt rePoints_n2[4] = {2, 3, 0, 1}; 171*c4762a1bSJed Brown if (rank) {rePoints_n2[0] = 1; rePoints_n2[1] = 2, rePoints_n2[2] = 0, rePoints_n2[3] = 3;} 172*c4762a1bSJed Brown 173*c4762a1bSJed Brown ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr); 174*c4762a1bSJed Brown ierr = PetscPartitionerShellSetPartition(part, size, reSizes_n2, rePoints_n2);CHKERRQ(ierr); 175*c4762a1bSJed Brown } 176*c4762a1bSJed Brown ierr = DMPlexSetPartitionBalance(*dm, user->partitionBalance);CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, overlap, NULL, &pdm);CHKERRQ(ierr); 178*c4762a1bSJed Brown if (pdm) { 179*c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 180*c4762a1bSJed Brown *dm = pdm; 181*c4762a1bSJed Brown } 182*c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 183*c4762a1bSJed Brown } 184*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *dm, cellSimplex ? "Simplicial Mesh" : "Tensor Product Mesh");CHKERRQ(ierr); 185*c4762a1bSJed Brown ierr = PetscLogStagePush(user->stages[STAGE_REFINE]);CHKERRQ(ierr); 186*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 187*c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 188*c4762a1bSJed Brown PetscFunctionReturn(0); 189*c4762a1bSJed Brown } 190*c4762a1bSJed Brown 191*c4762a1bSJed Brown 192*c4762a1bSJed Brown int main(int argc, char **argv) 193*c4762a1bSJed Brown { 194*c4762a1bSJed Brown DM dm; 195*c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 196*c4762a1bSJed Brown PetscErrorCode ierr; 197*c4762a1bSJed Brown 198*c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 199*c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 200*c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 201*c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 202*c4762a1bSJed Brown ierr = PetscFinalize(); 203*c4762a1bSJed Brown return ierr; 204*c4762a1bSJed Brown } 205*c4762a1bSJed Brown 206*c4762a1bSJed Brown /*TEST 207*c4762a1bSJed Brown # Parallel, no overlap tests 0-2 208*c4762a1bSJed Brown test: 209*c4762a1bSJed Brown suffix: 0 210*c4762a1bSJed Brown requires: triangle 211*c4762a1bSJed Brown args: -dm_view ascii:mesh.tex:ascii_latex 212*c4762a1bSJed Brown test: 213*c4762a1bSJed Brown suffix: 1 214*c4762a1bSJed Brown requires: triangle 215*c4762a1bSJed Brown nsize: 3 216*c4762a1bSJed Brown args: -test_partition -dm_view ascii::ascii_info_detail 217*c4762a1bSJed Brown test: 218*c4762a1bSJed Brown suffix: 2 219*c4762a1bSJed Brown requires: triangle 220*c4762a1bSJed Brown nsize: 8 221*c4762a1bSJed Brown args: -test_partition -dm_view ascii::ascii_info_detail 222*c4762a1bSJed Brown # Parallel, level-1 overlap tests 3-4 223*c4762a1bSJed Brown test: 224*c4762a1bSJed Brown suffix: 3 225*c4762a1bSJed Brown requires: triangle 226*c4762a1bSJed Brown nsize: 3 227*c4762a1bSJed Brown args: -test_partition -overlap 1 -dm_view ascii::ascii_info_detail 228*c4762a1bSJed Brown test: 229*c4762a1bSJed Brown suffix: 4 230*c4762a1bSJed Brown requires: triangle 231*c4762a1bSJed Brown nsize: 8 232*c4762a1bSJed Brown args: -test_partition -overlap 1 -dm_view ascii::ascii_info_detail 233*c4762a1bSJed Brown # Parallel, level-2 overlap test 5 234*c4762a1bSJed Brown test: 235*c4762a1bSJed Brown suffix: 5 236*c4762a1bSJed Brown requires: triangle 237*c4762a1bSJed Brown nsize: 8 238*c4762a1bSJed Brown args: -test_partition -overlap 2 -dm_view ascii::ascii_info_detail 239*c4762a1bSJed Brown # Parallel load balancing, test 6-7 240*c4762a1bSJed Brown test: 241*c4762a1bSJed Brown suffix: 6 242*c4762a1bSJed Brown requires: triangle 243*c4762a1bSJed Brown nsize: 2 244*c4762a1bSJed Brown args: -test_partition -overlap 1 -dm_view ascii::ascii_info_detail 245*c4762a1bSJed Brown test: 246*c4762a1bSJed Brown suffix: 7 247*c4762a1bSJed Brown requires: triangle 248*c4762a1bSJed Brown nsize: 2 249*c4762a1bSJed Brown args: -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail 250*c4762a1bSJed Brown # Parallel redundant copying, test 8 251*c4762a1bSJed Brown test: 252*c4762a1bSJed Brown suffix: 8 253*c4762a1bSJed Brown requires: triangle 254*c4762a1bSJed Brown nsize: 2 255*c4762a1bSJed Brown args: -test_redundant -redundant_migrated_dm_view ascii::ascii_info_detail -dm_view ascii::ascii_info_detail 256*c4762a1bSJed Brown test: 257*c4762a1bSJed Brown suffix: lb_0 258*c4762a1bSJed Brown requires: parmetis 259*c4762a1bSJed Brown nsize: 4 260*c4762a1bSJed Brown args: -cell_simplex 0 -cells 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 261*c4762a1bSJed Brown 262*c4762a1bSJed Brown # Same tests as above, but with balancing of the shared point partition 263*c4762a1bSJed Brown test: 264*c4762a1bSJed Brown suffix: 9 265*c4762a1bSJed Brown requires: triangle 266*c4762a1bSJed Brown args: -dm_view ascii:mesh.tex:ascii_latex -partition_balance 267*c4762a1bSJed Brown test: 268*c4762a1bSJed Brown suffix: 10 269*c4762a1bSJed Brown requires: triangle 270*c4762a1bSJed Brown nsize: 3 271*c4762a1bSJed Brown args: -test_partition -dm_view ascii::ascii_info_detail -partition_balance 272*c4762a1bSJed Brown test: 273*c4762a1bSJed Brown suffix: 11 274*c4762a1bSJed Brown requires: triangle 275*c4762a1bSJed Brown nsize: 8 276*c4762a1bSJed Brown args: -test_partition -dm_view ascii::ascii_info_detail -partition_balance 277*c4762a1bSJed Brown # Parallel, level-1 overlap tests 3-4 278*c4762a1bSJed Brown test: 279*c4762a1bSJed Brown suffix: 12 280*c4762a1bSJed Brown requires: triangle 281*c4762a1bSJed Brown nsize: 3 282*c4762a1bSJed Brown args: -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance 283*c4762a1bSJed Brown test: 284*c4762a1bSJed Brown suffix: 13 285*c4762a1bSJed Brown requires: triangle 286*c4762a1bSJed Brown nsize: 8 287*c4762a1bSJed Brown args: -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance 288*c4762a1bSJed Brown # Parallel, level-2 overlap test 5 289*c4762a1bSJed Brown test: 290*c4762a1bSJed Brown suffix: 14 291*c4762a1bSJed Brown requires: triangle 292*c4762a1bSJed Brown nsize: 8 293*c4762a1bSJed Brown args: -test_partition -overlap 2 -dm_view ascii::ascii_info_detail -partition_balance 294*c4762a1bSJed Brown # Parallel load balancing, test 6-7 295*c4762a1bSJed Brown test: 296*c4762a1bSJed Brown suffix: 15 297*c4762a1bSJed Brown requires: triangle 298*c4762a1bSJed Brown nsize: 2 299*c4762a1bSJed Brown args: -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance 300*c4762a1bSJed Brown test: 301*c4762a1bSJed Brown suffix: 16 302*c4762a1bSJed Brown requires: triangle 303*c4762a1bSJed Brown nsize: 2 304*c4762a1bSJed Brown args: -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail -partition_balance 305*c4762a1bSJed Brown # Parallel redundant copying, test 8 306*c4762a1bSJed Brown test: 307*c4762a1bSJed Brown suffix: 17 308*c4762a1bSJed Brown requires: triangle 309*c4762a1bSJed Brown nsize: 2 310*c4762a1bSJed Brown args: -test_redundant -dm_view ascii::ascii_info_detail -partition_balance 311*c4762a1bSJed Brown test: 312*c4762a1bSJed Brown suffix: lb_1 313*c4762a1bSJed Brown requires: parmetis 314*c4762a1bSJed Brown nsize: 4 315*c4762a1bSJed Brown args: -cell_simplex 0 -cells 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 316*c4762a1bSJed Brown TEST*/ 317