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