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 1030266083SMatthew G. Knepley make -f ./gmakefile test search="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, distribute it, and then redistribute it on 24 processes using two different partitioners: 13c4762a1bSJed Brown 1430266083SMatthew G. Knepley make -f ./gmakefile test search="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 15c4762a1bSJed Brown 16c4762a1bSJed 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: 17c4762a1bSJed Brown 1830266083SMatthew G. Knepley make -f ./gmakefile test search="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 19c4762a1bSJed Brown 20c4762a1bSJed Brown */ 21c4762a1bSJed Brown 229371c9d4SSatish Balay enum { 239371c9d4SSatish Balay STAGE_LOAD, 249371c9d4SSatish Balay STAGE_DISTRIBUTE, 259371c9d4SSatish Balay STAGE_REFINE, 269371c9d4SSatish Balay STAGE_REDISTRIBUTE 279371c9d4SSatish Balay }; 28c4762a1bSJed Brown 29c4762a1bSJed Brown typedef struct { 30c4762a1bSJed Brown /* Domain and mesh definition */ 31c4762a1bSJed Brown PetscInt overlap; /* The cell overlap to use during partitioning */ 32*3b9d9b65SStefano Zampini PetscBool testSection; /* Use a PetscSection to specify cell weights */ 33c4762a1bSJed Brown PetscBool testPartition; /* Use a fixed partitioning for testing */ 34c4762a1bSJed Brown PetscBool testRedundant; /* Use a redundant partitioning for testing */ 35c4762a1bSJed Brown PetscBool loadBalance; /* Load balance via a second distribute step */ 36c4762a1bSJed Brown PetscBool partitionBalance; /* Balance shared point partition */ 37c4762a1bSJed Brown PetscLogStage stages[4]; 38c4762a1bSJed Brown } AppCtx; 39c4762a1bSJed Brown 40d71ae5a4SJacob Faibussowitsch PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 41d71ae5a4SJacob Faibussowitsch { 42c4762a1bSJed Brown PetscFunctionBegin; 43c4762a1bSJed Brown options->overlap = 0; 44c4762a1bSJed Brown options->testPartition = PETSC_FALSE; 45*3b9d9b65SStefano Zampini options->testSection = PETSC_FALSE; 46c4762a1bSJed Brown options->testRedundant = PETSC_FALSE; 47c4762a1bSJed Brown options->loadBalance = PETSC_FALSE; 48c4762a1bSJed Brown options->partitionBalance = PETSC_FALSE; 49c4762a1bSJed Brown 50d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 519566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex12.c", options->overlap, &options->overlap, NULL, 0)); 529566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex12.c", options->testPartition, &options->testPartition, NULL)); 53*3b9d9b65SStefano Zampini PetscCall(PetscOptionsBool("-test_section", "Use a PetscSection for cell weights", "ex12.c", options->testSection, &options->testSection, NULL)); 549566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_redundant", "Use a redundant partition for testing", "ex12.c", options->testRedundant, &options->testRedundant, NULL)); 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-load_balance", "Perform parallel load balancing in a second distribution step", "ex12.c", options->loadBalance, &options->loadBalance, NULL)); 569566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-partition_balance", "Balance the ownership of shared points", "ex12.c", options->partitionBalance, &options->partitionBalance, NULL)); 57d0609cedSBarry Smith PetscOptionsEnd(); 58c4762a1bSJed Brown 599566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("MeshLoad", &options->stages[STAGE_LOAD])); 609566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("MeshDistribute", &options->stages[STAGE_DISTRIBUTE])); 619566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("MeshRefine", &options->stages[STAGE_REFINE])); 629566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("MeshRedistribute", &options->stages[STAGE_REDISTRIBUTE])); 633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 67d71ae5a4SJacob Faibussowitsch { 68c4762a1bSJed Brown DM pdm = NULL; 69c4762a1bSJed Brown PetscInt triSizes_n2[2] = {4, 4}; 70c4762a1bSJed Brown PetscInt triPoints_n2[8] = {0, 1, 4, 6, 2, 3, 5, 7}; 71c4762a1bSJed Brown PetscInt triSizes_n3[3] = {3, 2, 3}; 72c4762a1bSJed Brown PetscInt triPoints_n3[8] = {3, 5, 6, 1, 7, 0, 2, 4}; 73c4762a1bSJed Brown PetscInt triSizes_n4[4] = {2, 2, 2, 2}; 74c4762a1bSJed Brown PetscInt triPoints_n4[8] = {0, 7, 1, 5, 2, 3, 4, 6}; 75c4762a1bSJed Brown PetscInt triSizes_n8[8] = {1, 1, 1, 1, 1, 1, 1, 1}; 76c4762a1bSJed Brown PetscInt triPoints_n8[8] = {0, 1, 2, 3, 4, 5, 6, 7}; 77c4762a1bSJed Brown PetscInt quadSizes[2] = {2, 2}; 78c4762a1bSJed Brown PetscInt quadPoints[4] = {2, 3, 0, 1}; 79c4762a1bSJed Brown PetscInt overlap = user->overlap >= 0 ? user->overlap : 0; 8030602db0SMatthew G. Knepley PetscInt dim; 8130602db0SMatthew G. Knepley PetscBool simplex; 82c4762a1bSJed Brown PetscMPIInt rank, size; 83c4762a1bSJed Brown 84c4762a1bSJed Brown PetscFunctionBegin; 859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 879566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[STAGE_LOAD])); 889566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 899566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 909566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 919566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 92*3b9d9b65SStefano Zampini PetscCall(DMGetDimension(*dm, &dim)); 93*3b9d9b65SStefano Zampini if (user->testSection) { 94*3b9d9b65SStefano Zampini PetscInt numComp[] = {1}; 95*3b9d9b65SStefano Zampini PetscInt numDof[] = {0, 0, 0, 1}; 96*3b9d9b65SStefano Zampini PetscSection section; 97*3b9d9b65SStefano Zampini 98*3b9d9b65SStefano Zampini PetscCall(DMSetNumFields(*dm, 1)); 99*3b9d9b65SStefano Zampini PetscCall(DMPlexCreateSection(*dm, NULL, numComp, numDof + 3 - dim, 0, NULL, NULL, NULL, NULL, §ion)); 100*3b9d9b65SStefano Zampini PetscCall(DMSetLocalSection(*dm, section)); 101*3b9d9b65SStefano Zampini PetscCall(PetscSectionViewFromOptions(section, NULL, "-cell_section_view")); 102*3b9d9b65SStefano Zampini PetscCall(PetscSectionDestroy(§ion)); 103*3b9d9b65SStefano Zampini } 1049566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 1059566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 1069566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(*dm, &simplex)); 1079566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[STAGE_DISTRIBUTE])); 108c4762a1bSJed Brown if (!user->testRedundant) { 109c4762a1bSJed Brown PetscPartitioner part; 110c4762a1bSJed Brown 1119566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(*dm, &part)); 1129566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 113*3b9d9b65SStefano Zampini PetscCall(PetscPartitionerViewFromOptions(part, NULL, "-view_partitioner_pre")); 1149566063dSJacob Faibussowitsch PetscCall(DMPlexSetPartitionBalance(*dm, user->partitionBalance)); 115c4762a1bSJed Brown if (user->testPartition) { 116c4762a1bSJed Brown const PetscInt *sizes = NULL; 117c4762a1bSJed Brown const PetscInt *points = NULL; 118c4762a1bSJed Brown 119dd400576SPatrick Sanan if (rank == 0) { 12030602db0SMatthew G. Knepley if (dim == 2 && simplex && size == 2) { 1219371c9d4SSatish Balay sizes = triSizes_n2; 1229371c9d4SSatish Balay points = triPoints_n2; 12330602db0SMatthew G. Knepley } else if (dim == 2 && simplex && size == 3) { 1249371c9d4SSatish Balay sizes = triSizes_n3; 1259371c9d4SSatish Balay points = triPoints_n3; 12630602db0SMatthew G. Knepley } else if (dim == 2 && simplex && size == 4) { 1279371c9d4SSatish Balay sizes = triSizes_n4; 1289371c9d4SSatish Balay points = triPoints_n4; 12930602db0SMatthew G. Knepley } else if (dim == 2 && simplex && size == 8) { 1309371c9d4SSatish Balay sizes = triSizes_n8; 1319371c9d4SSatish Balay points = triPoints_n8; 13230602db0SMatthew G. Knepley } else if (dim == 2 && !simplex && size == 2) { 1339371c9d4SSatish Balay sizes = quadSizes; 1349371c9d4SSatish Balay points = quadPoints; 135c4762a1bSJed Brown } 136c4762a1bSJed Brown } 1379566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL)); 1389566063dSJacob Faibussowitsch PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points)); 139c4762a1bSJed Brown } 1409566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(*dm, overlap, NULL, &pdm)); 141*3b9d9b65SStefano Zampini PetscCall(PetscPartitionerViewFromOptions(part, NULL, "-view_partitioner")); 142c4762a1bSJed Brown } else { 143c4762a1bSJed Brown PetscSF sf; 144c4762a1bSJed Brown 1459566063dSJacob Faibussowitsch PetscCall(DMPlexGetRedundantDM(*dm, &sf, &pdm)); 146c4762a1bSJed Brown if (sf) { 147c4762a1bSJed Brown DM test; 148c4762a1bSJed Brown 1499566063dSJacob Faibussowitsch PetscCall(DMPlexCreate(comm, &test)); 1509566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)test, "Test SF-migrated Redundant Mesh")); 1519566063dSJacob Faibussowitsch PetscCall(DMPlexMigrate(*dm, sf, test)); 1529566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(test, NULL, "-redundant_migrated_dm_view")); 1539566063dSJacob Faibussowitsch PetscCall(DMDestroy(&test)); 154c4762a1bSJed Brown } 1559566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown if (pdm) { 1589566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 159c4762a1bSJed Brown *dm = pdm; 160c4762a1bSJed Brown } 1619566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 1629566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 163c4762a1bSJed Brown if (user->loadBalance) { 164c4762a1bSJed Brown PetscPartitioner part; 165c4762a1bSJed Brown 1669566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-prelb_dm_view")); 1679566063dSJacob Faibussowitsch PetscCall(DMPlexSetOptionsPrefix(*dm, "lb_")); 1689566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[STAGE_REDISTRIBUTE])); 1699566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(*dm, &part)); 170*3b9d9b65SStefano Zampini PetscCall(PetscPartitionerViewFromOptions(part, NULL, "-view_partitioner_pre")); 1719566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part, "lb_")); 1729566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 173c4762a1bSJed Brown if (user->testPartition) { 174c4762a1bSJed Brown PetscInt reSizes_n2[2] = {2, 2}; 175c4762a1bSJed Brown PetscInt rePoints_n2[4] = {2, 3, 0, 1}; 1769371c9d4SSatish Balay if (rank) { 1779371c9d4SSatish Balay rePoints_n2[0] = 1; 1789371c9d4SSatish Balay rePoints_n2[1] = 2, rePoints_n2[2] = 0, rePoints_n2[3] = 3; 1799371c9d4SSatish Balay } 180c4762a1bSJed Brown 1819566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL)); 1829566063dSJacob Faibussowitsch PetscCall(PetscPartitionerShellSetPartition(part, size, reSizes_n2, rePoints_n2)); 183c4762a1bSJed Brown } 1849566063dSJacob Faibussowitsch PetscCall(DMPlexSetPartitionBalance(*dm, user->partitionBalance)); 1859566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(*dm, overlap, NULL, &pdm)); 186*3b9d9b65SStefano Zampini PetscCall(PetscPartitionerViewFromOptions(part, NULL, "-view_partitioner")); 187c4762a1bSJed Brown if (pdm) { 1889566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 189c4762a1bSJed Brown *dm = pdm; 190c4762a1bSJed Brown } 1919566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 192c4762a1bSJed Brown } 1939566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[STAGE_REFINE])); 1949566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 1959566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 200d71ae5a4SJacob Faibussowitsch { 201c4762a1bSJed Brown DM dm; 202c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 203c4762a1bSJed Brown 204327415f7SBarry Smith PetscFunctionBeginUser; 2059566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2069566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 2079566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2089566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2099566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 210b122ec5aSJacob Faibussowitsch return 0; 211c4762a1bSJed Brown } 212c4762a1bSJed Brown 213c4762a1bSJed Brown /*TEST 214c4762a1bSJed Brown # Parallel, no overlap tests 0-2 215c4762a1bSJed Brown test: 216c4762a1bSJed Brown suffix: 0 217c4762a1bSJed Brown requires: triangle 218*3b9d9b65SStefano Zampini args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex -petscpartitioner_type {{simple multistage}} 219e0008caeSPierre Jolivet output_file: output/empty.out 220c4762a1bSJed Brown test: 221c4762a1bSJed Brown suffix: 1 222c4762a1bSJed Brown requires: triangle 223c4762a1bSJed Brown nsize: 3 22430602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail 225c4762a1bSJed Brown test: 226*3b9d9b65SStefano Zampini suffix: 1_ms_default 227*3b9d9b65SStefano Zampini requires: triangle 228*3b9d9b65SStefano Zampini nsize: 3 229*3b9d9b65SStefano Zampini args: -dm_coord_space 0 -petscpartitioner_type multistage -petscpartitioner_multistage_levels_petscpartitioner_type simple 230*3b9d9b65SStefano Zampini output_file: output/empty.out 231*3b9d9b65SStefano Zampini test: 232*3b9d9b65SStefano Zampini suffix: 1_ms_cellsection 233*3b9d9b65SStefano Zampini requires: triangle 234*3b9d9b65SStefano Zampini nsize: 3 235*3b9d9b65SStefano Zampini args: -dm_coord_space 0 -test_section -dm_view ascii::ascii_info_detail -petscpartitioner_type multistage -petscpartitioner_multistage_levels_petscpartitioner_type simple -petscpartitioner_multistage_node_size 2 -view_partitioner ::ascii_info_detail 236*3b9d9b65SStefano Zampini test: 237c4762a1bSJed Brown suffix: 2 238c4762a1bSJed Brown requires: triangle 239c4762a1bSJed Brown nsize: 8 24030602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail 241*3b9d9b65SStefano Zampini test: 242*3b9d9b65SStefano Zampini suffix: 2_ms 243*3b9d9b65SStefano Zampini requires: triangle 244*3b9d9b65SStefano Zampini nsize: 8 245*3b9d9b65SStefano Zampini args: -dm_coord_space 0 -dm_view ascii::ascii_info_detail \ 246*3b9d9b65SStefano Zampini -petscpartitioner_type multistage \ 247*3b9d9b65SStefano Zampini -petscpartitioner_multistage_strategy node -petscpartitioner_multistage_node_size {{1 2 3 4}separate output} -petscpartitioner_multistage_node_interleaved {{0 1}separate output} \ 248*3b9d9b65SStefano Zampini -petscpartitioner_multistage_levels_petscpartitioner_type simple -petscpartitioner_view ascii::ascii_info_detail 249c4762a1bSJed Brown # Parallel, level-1 overlap tests 3-4 250c4762a1bSJed Brown test: 251c4762a1bSJed Brown suffix: 3 252c4762a1bSJed Brown requires: triangle 253c4762a1bSJed Brown nsize: 3 25430602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail 255c4762a1bSJed Brown test: 256c4762a1bSJed Brown suffix: 4 257c4762a1bSJed Brown requires: triangle 258c4762a1bSJed Brown nsize: 8 25930602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail 260c4762a1bSJed Brown # Parallel, level-2 overlap test 5 261c4762a1bSJed Brown test: 262c4762a1bSJed Brown suffix: 5 263c4762a1bSJed Brown requires: triangle 264c4762a1bSJed Brown nsize: 8 26530602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail 266*3b9d9b65SStefano Zampini test: 267*3b9d9b65SStefano Zampini suffix: 5_ms 268*3b9d9b65SStefano Zampini requires: triangle 269*3b9d9b65SStefano Zampini nsize: 8 270*3b9d9b65SStefano Zampini args: -dm_coord_space 0 -overlap 2 -dm_view ascii::ascii_info_detail \ 271*3b9d9b65SStefano Zampini -petscpartitioner_type multistage \ 272*3b9d9b65SStefano Zampini -petscpartitioner_multistage_strategy msection -petscpartitioner_multistage_msection {{2 3}separate output} \ 273*3b9d9b65SStefano Zampini -petscpartitioner_multistage_levels_petscpartitioner_type simple -petscpartitioner_view ascii::ascii_info_detail 274c4762a1bSJed Brown # Parallel load balancing, test 6-7 275c4762a1bSJed Brown test: 276c4762a1bSJed Brown suffix: 6 277c4762a1bSJed Brown requires: triangle 278c4762a1bSJed Brown nsize: 2 27930602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail 280c4762a1bSJed Brown test: 281c4762a1bSJed Brown suffix: 7 282c4762a1bSJed Brown requires: triangle 283c4762a1bSJed Brown nsize: 2 28430602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail 285c4762a1bSJed Brown # Parallel redundant copying, test 8 286c4762a1bSJed Brown test: 287c4762a1bSJed Brown suffix: 8 288c4762a1bSJed Brown requires: triangle 289c4762a1bSJed Brown nsize: 2 29030602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_redundant -redundant_migrated_dm_view ascii::ascii_info_detail -dm_view ascii::ascii_info_detail 291c4762a1bSJed Brown test: 292c4762a1bSJed Brown suffix: lb_0 293c4762a1bSJed Brown requires: parmetis 294c4762a1bSJed Brown nsize: 4 29530602db0SMatthew 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 296*3b9d9b65SStefano Zampini test: 297*3b9d9b65SStefano Zampini suffix: lb_0_ms 298*3b9d9b65SStefano Zampini requires: parmetis 299*3b9d9b65SStefano Zampini nsize: 4 300*3b9d9b65SStefano Zampini args: -dm_coord_space 0 -dm_plex_simplex 0 -dm_plex_box_faces 4,4 -petscpartitioner_type shell -petscpartitioner_shell_random -lb_petscpartitioner_type multistage -load_balance -lb_petscpartitioner_view ::ascii_info_detail -prelb_dm_view ::load_balance -dm_view ::load_balance -lb_petscpartitioner_multistage_levels_petscpartitioner_type parmetis -lb_petscpartitioner_multistage_node_size 2 301c4762a1bSJed Brown 302c4762a1bSJed Brown # Same tests as above, but with balancing of the shared point partition 303c4762a1bSJed Brown test: 304c4762a1bSJed Brown suffix: 9 305c4762a1bSJed Brown requires: triangle 306*3b9d9b65SStefano Zampini args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex -partition_balance -petscpartitioner_type {{simple multistage}} 307e0008caeSPierre Jolivet output_file: output/empty.out 308c4762a1bSJed Brown test: 309c4762a1bSJed Brown suffix: 10 310c4762a1bSJed Brown requires: triangle 311c4762a1bSJed Brown nsize: 3 31230602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance 313c4762a1bSJed Brown test: 314c4762a1bSJed Brown suffix: 11 315c4762a1bSJed Brown requires: triangle 316c4762a1bSJed Brown nsize: 8 31730602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance 318*3b9d9b65SStefano Zampini test: 319*3b9d9b65SStefano Zampini suffix: 11_ms 320*3b9d9b65SStefano Zampini requires: triangle 321*3b9d9b65SStefano Zampini nsize: 8 322*3b9d9b65SStefano Zampini args: -dm_coord_space 0 -dm_view ascii::ascii_info_detail -partition_balance \ 323*3b9d9b65SStefano Zampini -petscpartitioner_type multistage \ 324*3b9d9b65SStefano Zampini -petscpartitioner_multistage_strategy node -petscpartitioner_multistage_node_size {{1 2 3 4}separate output} -petscpartitioner_multistage_node_interleaved {{0 1}separate output} \ 325*3b9d9b65SStefano Zampini -petscpartitioner_multistage_levels_petscpartitioner_type simple -petscpartitioner_view ascii::ascii_info_detail 326c4762a1bSJed Brown # Parallel, level-1 overlap tests 3-4 327c4762a1bSJed Brown test: 328c4762a1bSJed Brown suffix: 12 329c4762a1bSJed Brown requires: triangle 330c4762a1bSJed Brown nsize: 3 33130602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown suffix: 13 334c4762a1bSJed Brown requires: triangle 335c4762a1bSJed Brown nsize: 8 33630602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance 337c4762a1bSJed Brown # Parallel, level-2 overlap test 5 338c4762a1bSJed Brown test: 339c4762a1bSJed Brown suffix: 14 340c4762a1bSJed Brown requires: triangle 341c4762a1bSJed Brown nsize: 8 34230602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail -partition_balance 343c4762a1bSJed Brown # Parallel load balancing, test 6-7 344c4762a1bSJed Brown test: 345c4762a1bSJed Brown suffix: 15 346c4762a1bSJed Brown requires: triangle 347c4762a1bSJed Brown nsize: 2 34830602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance 349c4762a1bSJed Brown test: 350c4762a1bSJed Brown suffix: 16 351c4762a1bSJed Brown requires: triangle 352c4762a1bSJed Brown nsize: 2 35330602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail -partition_balance 354c4762a1bSJed Brown # Parallel redundant copying, test 8 355c4762a1bSJed Brown test: 356c4762a1bSJed Brown suffix: 17 357c4762a1bSJed Brown requires: triangle 358c4762a1bSJed Brown nsize: 2 35930602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_redundant -dm_view ascii::ascii_info_detail -partition_balance 360c4762a1bSJed Brown test: 361c4762a1bSJed Brown suffix: lb_1 362c4762a1bSJed Brown requires: parmetis 363c4762a1bSJed Brown nsize: 4 36430602db0SMatthew 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 365c4762a1bSJed Brown TEST*/ 366