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 26*9371c9d4SSatish Balay enum { 27*9371c9d4SSatish Balay STAGE_LOAD, 28*9371c9d4SSatish Balay STAGE_DISTRIBUTE, 29*9371c9d4SSatish Balay STAGE_REFINE, 30*9371c9d4SSatish Balay STAGE_REDISTRIBUTE 31*9371c9d4SSatish Balay }; 32c4762a1bSJed Brown 33c4762a1bSJed Brown typedef struct { 34c4762a1bSJed Brown /* Domain and mesh definition */ 35c4762a1bSJed Brown PetscInt overlap; /* The cell overlap to use during partitioning */ 36c4762a1bSJed Brown PetscBool testPartition; /* Use a fixed partitioning for testing */ 37c4762a1bSJed Brown PetscBool testRedundant; /* Use a redundant partitioning for testing */ 38c4762a1bSJed Brown PetscBool loadBalance; /* Load balance via a second distribute step */ 39c4762a1bSJed Brown PetscBool partitionBalance; /* Balance shared point partition */ 40c4762a1bSJed Brown PetscLogStage stages[4]; 41c4762a1bSJed Brown } AppCtx; 42c4762a1bSJed Brown 43*9371c9d4SSatish Balay PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 44c4762a1bSJed Brown PetscFunctionBegin; 45c4762a1bSJed Brown options->overlap = 0; 46c4762a1bSJed Brown options->testPartition = PETSC_FALSE; 47c4762a1bSJed Brown options->testRedundant = PETSC_FALSE; 48c4762a1bSJed Brown options->loadBalance = PETSC_FALSE; 49c4762a1bSJed Brown options->partitionBalance = PETSC_FALSE; 50c4762a1bSJed Brown 51d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 529566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-overlap", "The cell overlap for partitioning", "ex12.c", options->overlap, &options->overlap, NULL, 0)); 539566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex12.c", options->testPartition, &options->testPartition, 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])); 63c4762a1bSJed Brown PetscFunctionReturn(0); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66*9371c9d4SSatish Balay PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 67c4762a1bSJed Brown DM pdm = NULL; 68c4762a1bSJed Brown PetscInt triSizes_n2[2] = {4, 4}; 69c4762a1bSJed Brown PetscInt triPoints_n2[8] = {0, 1, 4, 6, 2, 3, 5, 7}; 70c4762a1bSJed Brown PetscInt triSizes_n3[3] = {3, 2, 3}; 71c4762a1bSJed Brown PetscInt triPoints_n3[8] = {3, 5, 6, 1, 7, 0, 2, 4}; 72c4762a1bSJed Brown PetscInt triSizes_n4[4] = {2, 2, 2, 2}; 73c4762a1bSJed Brown PetscInt triPoints_n4[8] = {0, 7, 1, 5, 2, 3, 4, 6}; 74c4762a1bSJed Brown PetscInt triSizes_n8[8] = {1, 1, 1, 1, 1, 1, 1, 1}; 75c4762a1bSJed Brown PetscInt triPoints_n8[8] = {0, 1, 2, 3, 4, 5, 6, 7}; 76c4762a1bSJed Brown PetscInt quadSizes[2] = {2, 2}; 77c4762a1bSJed Brown PetscInt quadPoints[4] = {2, 3, 0, 1}; 78c4762a1bSJed Brown PetscInt overlap = user->overlap >= 0 ? user->overlap : 0; 7930602db0SMatthew G. Knepley PetscInt dim; 8030602db0SMatthew G. Knepley PetscBool simplex; 81c4762a1bSJed Brown PetscMPIInt rank, size; 82c4762a1bSJed Brown 83c4762a1bSJed Brown PetscFunctionBegin; 849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 869566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[STAGE_LOAD])); 879566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 889566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 899566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 909566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 919566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 929566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 939566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 949566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(*dm, &simplex)); 959566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[STAGE_DISTRIBUTE])); 96c4762a1bSJed Brown if (!user->testRedundant) { 97c4762a1bSJed Brown PetscPartitioner part; 98c4762a1bSJed Brown 999566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(*dm, &part)); 1009566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 1019566063dSJacob Faibussowitsch PetscCall(DMPlexSetPartitionBalance(*dm, user->partitionBalance)); 102c4762a1bSJed Brown if (user->testPartition) { 103c4762a1bSJed Brown const PetscInt *sizes = NULL; 104c4762a1bSJed Brown const PetscInt *points = NULL; 105c4762a1bSJed Brown 106dd400576SPatrick Sanan if (rank == 0) { 10730602db0SMatthew G. Knepley if (dim == 2 && simplex && size == 2) { 108*9371c9d4SSatish Balay sizes = triSizes_n2; 109*9371c9d4SSatish Balay points = triPoints_n2; 11030602db0SMatthew G. Knepley } else if (dim == 2 && simplex && size == 3) { 111*9371c9d4SSatish Balay sizes = triSizes_n3; 112*9371c9d4SSatish Balay points = triPoints_n3; 11330602db0SMatthew G. Knepley } else if (dim == 2 && simplex && size == 4) { 114*9371c9d4SSatish Balay sizes = triSizes_n4; 115*9371c9d4SSatish Balay points = triPoints_n4; 11630602db0SMatthew G. Knepley } else if (dim == 2 && simplex && size == 8) { 117*9371c9d4SSatish Balay sizes = triSizes_n8; 118*9371c9d4SSatish Balay points = triPoints_n8; 11930602db0SMatthew G. Knepley } else if (dim == 2 && !simplex && size == 2) { 120*9371c9d4SSatish Balay sizes = quadSizes; 121*9371c9d4SSatish Balay points = quadPoints; 122c4762a1bSJed Brown } 123c4762a1bSJed Brown } 1249566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL)); 1259566063dSJacob Faibussowitsch PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points)); 126c4762a1bSJed Brown } 1279566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(*dm, overlap, NULL, &pdm)); 128c4762a1bSJed Brown } else { 129c4762a1bSJed Brown PetscSF sf; 130c4762a1bSJed Brown 1319566063dSJacob Faibussowitsch PetscCall(DMPlexGetRedundantDM(*dm, &sf, &pdm)); 132c4762a1bSJed Brown if (sf) { 133c4762a1bSJed Brown DM test; 134c4762a1bSJed Brown 1359566063dSJacob Faibussowitsch PetscCall(DMPlexCreate(comm, &test)); 1369566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)test, "Test SF-migrated Redundant Mesh")); 1379566063dSJacob Faibussowitsch PetscCall(DMPlexMigrate(*dm, sf, test)); 1389566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(test, NULL, "-redundant_migrated_dm_view")); 1399566063dSJacob Faibussowitsch PetscCall(DMDestroy(&test)); 140c4762a1bSJed Brown } 1419566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown if (pdm) { 1449566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 145c4762a1bSJed Brown *dm = pdm; 146c4762a1bSJed Brown } 1479566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 1489566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 149c4762a1bSJed Brown if (user->loadBalance) { 150c4762a1bSJed Brown PetscPartitioner part; 151c4762a1bSJed Brown 1529566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-prelb_dm_view")); 1539566063dSJacob Faibussowitsch PetscCall(DMPlexSetOptionsPrefix(*dm, "lb_")); 1549566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[STAGE_REDISTRIBUTE])); 1559566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(*dm, &part)); 1569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part, "lb_")); 1579566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 158c4762a1bSJed Brown if (user->testPartition) { 159c4762a1bSJed Brown PetscInt reSizes_n2[2] = {2, 2}; 160c4762a1bSJed Brown PetscInt rePoints_n2[4] = {2, 3, 0, 1}; 161*9371c9d4SSatish Balay if (rank) { 162*9371c9d4SSatish Balay rePoints_n2[0] = 1; 163*9371c9d4SSatish Balay rePoints_n2[1] = 2, rePoints_n2[2] = 0, rePoints_n2[3] = 3; 164*9371c9d4SSatish Balay } 165c4762a1bSJed Brown 1669566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL)); 1679566063dSJacob Faibussowitsch PetscCall(PetscPartitionerShellSetPartition(part, size, reSizes_n2, rePoints_n2)); 168c4762a1bSJed Brown } 1699566063dSJacob Faibussowitsch PetscCall(DMPlexSetPartitionBalance(*dm, user->partitionBalance)); 1709566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(*dm, overlap, NULL, &pdm)); 171c4762a1bSJed Brown if (pdm) { 1729566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 173c4762a1bSJed Brown *dm = pdm; 174c4762a1bSJed Brown } 1759566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 176c4762a1bSJed Brown } 1779566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(user->stages[STAGE_REFINE])); 1789566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 1799566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 180c4762a1bSJed Brown PetscFunctionReturn(0); 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183*9371c9d4SSatish Balay int main(int argc, char **argv) { 184c4762a1bSJed Brown DM dm; 185c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 186c4762a1bSJed Brown 187327415f7SBarry Smith PetscFunctionBeginUser; 1889566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1899566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 1909566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 1919566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1929566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 193b122ec5aSJacob Faibussowitsch return 0; 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196c4762a1bSJed Brown /*TEST 197c4762a1bSJed Brown # Parallel, no overlap tests 0-2 198c4762a1bSJed Brown test: 199c4762a1bSJed Brown suffix: 0 200c4762a1bSJed Brown requires: triangle 20130602db0SMatthew G. Knepley args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex 202c4762a1bSJed Brown test: 203c4762a1bSJed Brown suffix: 1 204c4762a1bSJed Brown requires: triangle 205c4762a1bSJed Brown nsize: 3 20630602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail 207c4762a1bSJed Brown test: 208c4762a1bSJed Brown suffix: 2 209c4762a1bSJed Brown requires: triangle 210c4762a1bSJed Brown nsize: 8 21130602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail 212c4762a1bSJed Brown # Parallel, level-1 overlap tests 3-4 213c4762a1bSJed Brown test: 214c4762a1bSJed Brown suffix: 3 215c4762a1bSJed Brown requires: triangle 216c4762a1bSJed Brown nsize: 3 21730602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail 218c4762a1bSJed Brown test: 219c4762a1bSJed Brown suffix: 4 220c4762a1bSJed Brown requires: triangle 221c4762a1bSJed Brown nsize: 8 22230602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail 223c4762a1bSJed Brown # Parallel, level-2 overlap test 5 224c4762a1bSJed Brown test: 225c4762a1bSJed Brown suffix: 5 226c4762a1bSJed Brown requires: triangle 227c4762a1bSJed Brown nsize: 8 22830602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail 229c4762a1bSJed Brown # Parallel load balancing, test 6-7 230c4762a1bSJed Brown test: 231c4762a1bSJed Brown suffix: 6 232c4762a1bSJed Brown requires: triangle 233c4762a1bSJed Brown nsize: 2 23430602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail 235c4762a1bSJed Brown test: 236c4762a1bSJed Brown suffix: 7 237c4762a1bSJed Brown requires: triangle 238c4762a1bSJed Brown nsize: 2 23930602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail 240c4762a1bSJed Brown # Parallel redundant copying, test 8 241c4762a1bSJed Brown test: 242c4762a1bSJed Brown suffix: 8 243c4762a1bSJed Brown requires: triangle 244c4762a1bSJed Brown nsize: 2 24530602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_redundant -redundant_migrated_dm_view ascii::ascii_info_detail -dm_view ascii::ascii_info_detail 246c4762a1bSJed Brown test: 247c4762a1bSJed Brown suffix: lb_0 248c4762a1bSJed Brown requires: parmetis 249c4762a1bSJed Brown nsize: 4 25030602db0SMatthew 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 251c4762a1bSJed Brown 252c4762a1bSJed Brown # Same tests as above, but with balancing of the shared point partition 253c4762a1bSJed Brown test: 254c4762a1bSJed Brown suffix: 9 255c4762a1bSJed Brown requires: triangle 25630602db0SMatthew G. Knepley args: -dm_coord_space 0 -dm_view ascii:mesh.tex:ascii_latex -partition_balance 257c4762a1bSJed Brown test: 258c4762a1bSJed Brown suffix: 10 259c4762a1bSJed Brown requires: triangle 260c4762a1bSJed Brown nsize: 3 26130602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance 262c4762a1bSJed Brown test: 263c4762a1bSJed Brown suffix: 11 264c4762a1bSJed Brown requires: triangle 265c4762a1bSJed Brown nsize: 8 26630602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -dm_view ascii::ascii_info_detail -partition_balance 267c4762a1bSJed Brown # Parallel, level-1 overlap tests 3-4 268c4762a1bSJed Brown test: 269c4762a1bSJed Brown suffix: 12 270c4762a1bSJed Brown requires: triangle 271c4762a1bSJed Brown nsize: 3 27230602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance 273c4762a1bSJed Brown test: 274c4762a1bSJed Brown suffix: 13 275c4762a1bSJed Brown requires: triangle 276c4762a1bSJed Brown nsize: 8 27730602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance 278c4762a1bSJed Brown # Parallel, level-2 overlap test 5 279c4762a1bSJed Brown test: 280c4762a1bSJed Brown suffix: 14 281c4762a1bSJed Brown requires: triangle 282c4762a1bSJed Brown nsize: 8 28330602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 2 -dm_view ascii::ascii_info_detail -partition_balance 284c4762a1bSJed Brown # Parallel load balancing, test 6-7 285c4762a1bSJed Brown test: 286c4762a1bSJed Brown suffix: 15 287c4762a1bSJed Brown requires: triangle 288c4762a1bSJed Brown nsize: 2 28930602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -dm_view ascii::ascii_info_detail -partition_balance 290c4762a1bSJed Brown test: 291c4762a1bSJed Brown suffix: 16 292c4762a1bSJed Brown requires: triangle 293c4762a1bSJed Brown nsize: 2 29430602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_partition -overlap 1 -load_balance -dm_view ascii::ascii_info_detail -partition_balance 295c4762a1bSJed Brown # Parallel redundant copying, test 8 296c4762a1bSJed Brown test: 297c4762a1bSJed Brown suffix: 17 298c4762a1bSJed Brown requires: triangle 299c4762a1bSJed Brown nsize: 2 30030602db0SMatthew G. Knepley args: -dm_coord_space 0 -test_redundant -dm_view ascii::ascii_info_detail -partition_balance 301c4762a1bSJed Brown test: 302c4762a1bSJed Brown suffix: lb_1 303c4762a1bSJed Brown requires: parmetis 304c4762a1bSJed Brown nsize: 4 30530602db0SMatthew 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 306c4762a1bSJed Brown TEST*/ 307