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