1c4762a1bSJed Brown static char help[] = "Test that MatPartitioning and PetscPartitioner interfaces are equivalent when using PETSCPARTITIONERMATPARTITIONING\n\n"; 2c4762a1bSJed Brown static char FILENAME[] = "ex24.c"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdmplex.h> 5c4762a1bSJed Brown #include <petscviewerhdf5.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown #if defined(PETSC_HAVE_PTSCOTCH) 8c4762a1bSJed Brown EXTERN_C_BEGIN 9c4762a1bSJed Brown #include <ptscotch.h> 10c4762a1bSJed Brown EXTERN_C_END 11c4762a1bSJed Brown #endif 12c4762a1bSJed Brown 13c4762a1bSJed Brown typedef struct { 14c4762a1bSJed Brown PetscBool compare_is; /* Compare ISs and PetscSections */ 15c4762a1bSJed Brown PetscBool compare_dm; /* Compare DM */ 16c4762a1bSJed Brown PetscBool tpw; /* Use target partition weights */ 17c4762a1bSJed Brown char partitioning[64]; 18c4762a1bSJed Brown char repartitioning[64]; 19c4762a1bSJed Brown } AppCtx; 20c4762a1bSJed Brown 21c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 22c4762a1bSJed Brown { 23c4762a1bSJed Brown PetscBool repartition = PETSC_TRUE; 24c4762a1bSJed Brown PetscErrorCode ierr; 25c4762a1bSJed Brown 26c4762a1bSJed Brown PetscFunctionBegin; 27c4762a1bSJed Brown options->compare_is = PETSC_FALSE; 28c4762a1bSJed Brown options->compare_dm = PETSC_FALSE; 2930602db0SMatthew G. Knepley 30c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscOptionsBool("-compare_is", "Compare ISs and PetscSections?", FILENAME, options->compare_is, &options->compare_is, NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = PetscOptionsBool("-compare_dm", "Compare DMs?", FILENAME, options->compare_dm, &options->compare_dm, NULL);CHKERRQ(ierr); 33589a23caSBarry Smith ierr = PetscStrncpy(options->partitioning,MATPARTITIONINGPARMETIS,sizeof(options->partitioning));CHKERRQ(ierr); 34589a23caSBarry Smith ierr = PetscOptionsString("-partitioning","The mat partitioning type to test","None",options->partitioning, options->partitioning,sizeof(options->partitioning),NULL);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = PetscOptionsBool("-repartition", "Partition again after the first partition?", FILENAME, repartition, &repartition, NULL);CHKERRQ(ierr); 36c4762a1bSJed Brown if (repartition) { 37c4762a1bSJed Brown ierr = PetscStrncpy(options->repartitioning,MATPARTITIONINGPARMETIS,64);CHKERRQ(ierr); 38589a23caSBarry Smith ierr = PetscOptionsString("-repartitioning","The mat partitioning type to test (second partitioning)","None", options->repartitioning, options->repartitioning,sizeof(options->repartitioning),NULL);CHKERRQ(ierr); 39c4762a1bSJed Brown } else { 40c4762a1bSJed Brown options->repartitioning[0] = '\0'; 41c4762a1bSJed Brown } 42c4762a1bSJed Brown ierr = PetscOptionsBool("-tpweight", "Use target partition weights", FILENAME, options->tpw, &options->tpw, NULL);CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = PetscOptionsEnd(); 44c4762a1bSJed Brown PetscFunctionReturn(0); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown static PetscErrorCode ScotchResetRandomSeed() 48c4762a1bSJed Brown { 49*362febeeSStefano Zampini PetscFunctionBegin; 50c4762a1bSJed Brown #if defined(PETSC_HAVE_PTSCOTCH) 51c4762a1bSJed Brown SCOTCH_randomReset(); 52c4762a1bSJed Brown #endif 53c4762a1bSJed Brown PetscFunctionReturn(0); 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 57c4762a1bSJed Brown { 58c4762a1bSJed Brown PetscErrorCode ierr; 59c4762a1bSJed Brown 60c4762a1bSJed Brown PetscFunctionBegin; 6130602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 6230602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 6330602db0SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 6430602db0SMatthew G. Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 65c4762a1bSJed Brown PetscFunctionReturn(0); 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68c4762a1bSJed Brown int main(int argc, char **argv) 69c4762a1bSJed Brown { 70c4762a1bSJed Brown MPI_Comm comm; 71c4762a1bSJed Brown DM dm1, dm2, dmdist1, dmdist2; 7230602db0SMatthew G. Knepley DMPlexInterpolatedFlag interp; 73c4762a1bSJed Brown MatPartitioning mp; 74c4762a1bSJed Brown PetscPartitioner part1, part2; 75c4762a1bSJed Brown AppCtx user; 76c4762a1bSJed Brown IS is1=NULL, is2=NULL; 77c4762a1bSJed Brown IS is1g, is2g; 78c4762a1bSJed Brown PetscSection s1=NULL, s2=NULL, tpws = NULL; 79c4762a1bSJed Brown PetscInt i; 80c4762a1bSJed Brown PetscBool flg; 81c4762a1bSJed Brown PetscErrorCode ierr; 82c4762a1bSJed Brown PetscMPIInt size; 83c4762a1bSJed Brown 84c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 85c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 86ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 87c4762a1bSJed Brown ierr = ProcessOptions(comm, &user);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = CreateMesh(comm, &user, &dm1);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = CreateMesh(comm, &user, &dm2);CHKERRQ(ierr); 90c4762a1bSJed Brown 91c4762a1bSJed Brown if (user.tpw) { 92c4762a1bSJed Brown ierr = PetscSectionCreate(comm, &tpws);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = PetscSectionSetChart(tpws, 0, size);CHKERRQ(ierr); 94c4762a1bSJed Brown for (i=0;i<size;i++) { 95c4762a1bSJed Brown PetscInt tdof = i%2 ? 2*i -1 : i+2; 96c4762a1bSJed Brown ierr = PetscSectionSetDof(tpws, i, tdof);CHKERRQ(ierr); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown if (size > 1) { /* test zero tpw entry */ 99c4762a1bSJed Brown ierr = PetscSectionSetDof(tpws, 0, 0);CHKERRQ(ierr); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown ierr = PetscSectionSetUp(tpws);CHKERRQ(ierr); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* partition dm1 using PETSCPARTITIONERPARMETIS */ 105c4762a1bSJed Brown ierr = ScotchResetRandomSeed();CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = DMPlexGetPartitioner(dm1, &part1);CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)part1,"p1_");CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = PetscPartitionerSetType(part1, user.partitioning);CHKERRQ(ierr); 109c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part1);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = PetscSectionCreate(comm, &s1);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = PetscPartitionerDMPlexPartition(part1, dm1, tpws, s1, &is1);CHKERRQ(ierr); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* partition dm2 using PETSCPARTITIONERMATPARTITIONING with MATPARTITIONINGPARMETIS */ 114c4762a1bSJed Brown ierr = ScotchResetRandomSeed();CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = DMPlexGetPartitioner(dm2, &part2);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)part2,"p2_");CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = PetscPartitionerSetType(part2, PETSCPARTITIONERMATPARTITIONING);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = PetscPartitionerMatPartitioningGetMatPartitioning(part2, &mp);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = MatPartitioningSetType(mp, user.partitioning);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part2);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = PetscSectionCreate(comm, &s2);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = PetscPartitionerDMPlexPartition(part2, dm2, tpws, s2, &is2);CHKERRQ(ierr); 123c4762a1bSJed Brown 124c4762a1bSJed Brown ierr = ISOnComm(is1, comm, PETSC_USE_POINTER, &is1g);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = ISOnComm(is2, comm, PETSC_USE_POINTER, &is2g);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = ISViewFromOptions(is1g, NULL, "-seq_is1_view");CHKERRQ(ierr); 127c4762a1bSJed Brown ierr = ISViewFromOptions(is2g, NULL, "-seq_is2_view");CHKERRQ(ierr); 128c4762a1bSJed Brown /* compare the two ISs */ 129c4762a1bSJed Brown if (user.compare_is) { 130c4762a1bSJed Brown ierr = ISEqualUnsorted(is1g, is2g, &flg);CHKERRQ(ierr); 131b5675b0fSBarry Smith if (!flg) {ierr = PetscPrintf(comm, "ISs are not equal with type %s with size %d.\n",user.partitioning,size);CHKERRQ(ierr);} 132c4762a1bSJed Brown } 133c4762a1bSJed Brown ierr = ISDestroy(&is1g);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = ISDestroy(&is2g);CHKERRQ(ierr); 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* compare the two PetscSections */ 137c4762a1bSJed Brown ierr = PetscSectionViewFromOptions(s1, NULL, "-seq_s1_view");CHKERRQ(ierr); 138c4762a1bSJed Brown ierr = PetscSectionViewFromOptions(s2, NULL, "-seq_s2_view");CHKERRQ(ierr); 139c4762a1bSJed Brown if (user.compare_is) { 140c4762a1bSJed Brown ierr = PetscSectionCompare(s1, s2, &flg);CHKERRQ(ierr); 141b5675b0fSBarry Smith if (!flg) {ierr = PetscPrintf(comm, "PetscSections are not equal with %s with size %d.\n",user.partitioning,size);CHKERRQ(ierr);} 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* distribute both DMs */ 145c4762a1bSJed Brown ierr = ScotchResetRandomSeed();CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = DMPlexDistribute(dm1, 0, NULL, &dmdist1);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = ScotchResetRandomSeed();CHKERRQ(ierr); 148c4762a1bSJed Brown ierr = DMPlexDistribute(dm2, 0, NULL, &dmdist2);CHKERRQ(ierr); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* cleanup */ 151c4762a1bSJed Brown ierr = PetscSectionDestroy(&tpws);CHKERRQ(ierr); 152c4762a1bSJed Brown ierr = PetscSectionDestroy(&s1);CHKERRQ(ierr); 153c4762a1bSJed Brown ierr = PetscSectionDestroy(&s2);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = ISDestroy(&is1);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = ISDestroy(&is2);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = DMDestroy(&dm1);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = DMDestroy(&dm2);CHKERRQ(ierr); 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* if distributed DMs are NULL (sequential case), then quit */ 160c4762a1bSJed Brown if (!dmdist1 && !dmdist2) return ierr; 161c4762a1bSJed Brown 162c4762a1bSJed Brown ierr = DMViewFromOptions(dmdist1, NULL, "-dm_dist1_view");CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = DMViewFromOptions(dmdist2, NULL, "-dm_dist2_view");CHKERRQ(ierr); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* compare the two distributed DMs */ 166c4762a1bSJed Brown if (user.compare_dm) { 167c4762a1bSJed Brown ierr = DMPlexEqual(dmdist1, dmdist2, &flg);CHKERRQ(ierr); 168b5675b0fSBarry Smith if (!flg) {ierr = PetscPrintf(comm, "Distributed DMs are not equal %s with size %d.\n",user.partitioning,size);CHKERRQ(ierr);} 169c4762a1bSJed Brown } 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* if repartitioning is disabled, then quit */ 172c4762a1bSJed Brown if (user.repartitioning[0] == '\0') return ierr; 173c4762a1bSJed Brown 174c4762a1bSJed Brown if (user.tpw) { 175c4762a1bSJed Brown ierr = PetscSectionCreate(comm, &tpws);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = PetscSectionSetChart(tpws, 0, size);CHKERRQ(ierr); 177c4762a1bSJed Brown for (i=0;i<size;i++) { 178c4762a1bSJed Brown PetscInt tdof = i%2 ? i+1 : size - i; 179c4762a1bSJed Brown ierr = PetscSectionSetDof(tpws, i, tdof);CHKERRQ(ierr); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown ierr = PetscSectionSetUp(tpws);CHKERRQ(ierr); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* repartition distributed DM dmdist1 */ 185c4762a1bSJed Brown ierr = ScotchResetRandomSeed();CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = DMPlexGetPartitioner(dmdist1, &part1);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)part1,"dp1_");CHKERRQ(ierr); 188c4762a1bSJed Brown ierr = PetscPartitionerSetType(part1, user.repartitioning);CHKERRQ(ierr); 189c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part1);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = PetscSectionCreate(comm, &s1);CHKERRQ(ierr); 191c4762a1bSJed Brown ierr = PetscPartitionerDMPlexPartition(part1, dmdist1, tpws, s1, &is1);CHKERRQ(ierr); 192c4762a1bSJed Brown 193c4762a1bSJed Brown /* repartition distributed DM dmdist2 */ 194c4762a1bSJed Brown ierr = ScotchResetRandomSeed();CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = DMPlexGetPartitioner(dmdist2, &part2);CHKERRQ(ierr); 196c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)part2,"dp2_");CHKERRQ(ierr); 197c4762a1bSJed Brown ierr = PetscPartitionerSetType(part2, PETSCPARTITIONERMATPARTITIONING);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = PetscPartitionerMatPartitioningGetMatPartitioning(part2, &mp);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = MatPartitioningSetType(mp, user.repartitioning);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part2);CHKERRQ(ierr); 201c4762a1bSJed Brown ierr = PetscSectionCreate(comm, &s2);CHKERRQ(ierr); 202c4762a1bSJed Brown ierr = PetscPartitionerDMPlexPartition(part2, dmdist2, tpws, s2, &is2);CHKERRQ(ierr); 203c4762a1bSJed Brown 204c4762a1bSJed Brown /* compare the two ISs */ 205c4762a1bSJed Brown ierr = ISOnComm(is1, comm, PETSC_USE_POINTER, &is1g);CHKERRQ(ierr); 206c4762a1bSJed Brown ierr = ISOnComm(is2, comm, PETSC_USE_POINTER, &is2g);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = ISViewFromOptions(is1g, NULL, "-dist_is1_view");CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = ISViewFromOptions(is2g, NULL, "-dist_is2_view");CHKERRQ(ierr); 209c4762a1bSJed Brown if (user.compare_is) { 210c4762a1bSJed Brown ierr = ISEqualUnsorted(is1g, is2g, &flg);CHKERRQ(ierr); 211b5675b0fSBarry Smith if (!flg) {ierr = PetscPrintf(comm, "Distributed ISs are not equal, with %s with size %d.\n",user.repartitioning,size);CHKERRQ(ierr);} 212c4762a1bSJed Brown } 213c4762a1bSJed Brown ierr = ISDestroy(&is1g);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = ISDestroy(&is2g);CHKERRQ(ierr); 215c4762a1bSJed Brown 216c4762a1bSJed Brown /* compare the two PetscSections */ 217c4762a1bSJed Brown ierr = PetscSectionViewFromOptions(s1, NULL, "-dist_s1_view");CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = PetscSectionViewFromOptions(s2, NULL, "-dist_s2_view");CHKERRQ(ierr); 219c4762a1bSJed Brown if (user.compare_is) { 220c4762a1bSJed Brown ierr = PetscSectionCompare(s1, s2, &flg);CHKERRQ(ierr); 221b5675b0fSBarry Smith if (!flg) {ierr = PetscPrintf(comm, "Distributed PetscSections are not equal, with %s with size %d.\n",user.repartitioning,size);CHKERRQ(ierr);} 222c4762a1bSJed Brown } 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* redistribute both distributed DMs */ 225c4762a1bSJed Brown ierr = ScotchResetRandomSeed();CHKERRQ(ierr); 226c4762a1bSJed Brown ierr = DMPlexDistribute(dmdist1, 0, NULL, &dm1);CHKERRQ(ierr); 227c4762a1bSJed Brown ierr = ScotchResetRandomSeed();CHKERRQ(ierr); 228c4762a1bSJed Brown ierr = DMPlexDistribute(dmdist2, 0, NULL, &dm2);CHKERRQ(ierr); 229c4762a1bSJed Brown 230c4762a1bSJed Brown /* compare the two distributed DMs */ 23130602db0SMatthew G. Knepley ierr = DMPlexIsInterpolated(dm1, &interp);CHKERRQ(ierr); 23230602db0SMatthew G. Knepley if (interp == DMPLEX_INTERPOLATED_NONE) { 233c4762a1bSJed Brown ierr = DMPlexEqual(dm1, dm2, &flg);CHKERRQ(ierr); 234b5675b0fSBarry Smith if (!flg) {ierr = PetscPrintf(comm, "Redistributed DMs are not equal, with %s with size %d.\n",user.repartitioning,size);CHKERRQ(ierr);} 235c4762a1bSJed Brown } 236c4762a1bSJed Brown 237c4762a1bSJed Brown /* cleanup */ 238c4762a1bSJed Brown ierr = PetscSectionDestroy(&tpws);CHKERRQ(ierr); 239c4762a1bSJed Brown ierr = PetscSectionDestroy(&s1);CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = PetscSectionDestroy(&s2);CHKERRQ(ierr); 241c4762a1bSJed Brown ierr = ISDestroy(&is1);CHKERRQ(ierr); 242c4762a1bSJed Brown ierr = ISDestroy(&is2);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = DMDestroy(&dm1);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = DMDestroy(&dm2);CHKERRQ(ierr); 245c4762a1bSJed Brown ierr = DMDestroy(&dmdist1);CHKERRQ(ierr); 246c4762a1bSJed Brown ierr = DMDestroy(&dmdist2);CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = PetscFinalize(); 248c4762a1bSJed Brown return ierr; 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251c4762a1bSJed Brown /*TEST 252c4762a1bSJed Brown 253c4762a1bSJed Brown test: 254c4762a1bSJed Brown # partition sequential mesh loaded from Exodus file 255c4762a1bSJed Brown suffix: 0 256c4762a1bSJed Brown nsize: {{1 2 3 4 8}} 257c4762a1bSJed Brown requires: chaco parmetis ptscotch exodusii 25830602db0SMatthew G. Knepley args: -dm_plex_filename ${PETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo 259c4762a1bSJed Brown args: -partitioning {{chaco parmetis ptscotch}} -repartitioning {{parmetis ptscotch}} -tpweight {{0 1}} 260c4762a1bSJed Brown test: 261c4762a1bSJed Brown # repartition mesh already partitioned naively by MED loader 262c4762a1bSJed Brown suffix: 1 263c4762a1bSJed Brown nsize: {{1 2 3 4 8}} 264c4762a1bSJed Brown TODO: MED 265c4762a1bSJed Brown requires: parmetis ptscotch med 26630602db0SMatthew G. Knepley args: -dm_plex_filename ${PETSC_DIR}/share/petsc/datafiles/meshes/cylinder.med 267c4762a1bSJed Brown args: -repartition 0 -partitioning {{parmetis ptscotch}} 268c4762a1bSJed Brown test: 269c4762a1bSJed Brown # partition mesh generated by ctetgen using scotch, then repartition with scotch, diff view 270c4762a1bSJed Brown suffix: 3 271c4762a1bSJed Brown nsize: 4 272c4762a1bSJed Brown requires: ptscotch ctetgen 27330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 2,3,2 -partitioning ptscotch -repartitioning ptscotch 274c4762a1bSJed Brown args: -p1_petscpartitioner_view -p2_petscpartitioner_view -dp1_petscpartitioner_view -dp2_petscpartitioner_view -tpweight {{0 1}} 275c4762a1bSJed Brown test: 276c4762a1bSJed Brown # partition mesh generated by ctetgen using partitioners supported both by MatPartitioning and PetscPartitioner 277c4762a1bSJed Brown suffix: 4 278c4762a1bSJed Brown nsize: {{1 2 3 4 8}} 279c4762a1bSJed Brown requires: chaco parmetis ptscotch ctetgen 28030602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces {{2,3,4 5,4,3 7,11,5}} -partitioning {{chaco parmetis ptscotch}} -repartitioning {{parmetis ptscotch}} -tpweight {{0 1}} 281c4762a1bSJed Brown 282c4762a1bSJed Brown TEST*/ 283