1abe9303eSLisandro Dalcin static char help[] = "Tests PetscPartitioner.\n\n"; 2abe9303eSLisandro Dalcin 3abe9303eSLisandro Dalcin #include <petscpartitioner.h> 4abe9303eSLisandro Dalcin 5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 6d71ae5a4SJacob Faibussowitsch { 7abe9303eSLisandro Dalcin PetscPartitioner p; 8d7cc930eSLisandro Dalcin PetscSection partSection, vertexSection = NULL, targetSection = NULL; 9abe9303eSLisandro Dalcin IS partition, is; 10abe9303eSLisandro Dalcin PetscMPIInt size, rank; 11d7cc930eSLisandro Dalcin PetscInt nparts, i; 12abe9303eSLisandro Dalcin PetscInt nv = 4; 13abe9303eSLisandro Dalcin PetscInt vv[5] = {0, 2, 4, 6, 8}; 14abe9303eSLisandro Dalcin PetscInt vadj[8] = {3, 1, 0, 2, 1, 3, 2, 0}; 1556991f02SLisandro Dalcin PetscBool sequential; 16d7cc930eSLisandro Dalcin PetscBool vwgts = PETSC_FALSE; 17d7cc930eSLisandro Dalcin PetscBool pwgts = PETSC_FALSE; 18abe9303eSLisandro Dalcin 19327415f7SBarry Smith PetscFunctionBeginUser; 209566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 23d7cc930eSLisandro Dalcin nparts = size; 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nparts", &nparts, NULL)); 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-vwgts", &vwgts, NULL)); 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-pwgts", &pwgts, NULL)); 27abe9303eSLisandro Dalcin 28abe9303eSLisandro Dalcin /* create PetscPartitioner */ 299566063dSJacob Faibussowitsch PetscCall(PetscPartitionerCreate(PETSC_COMM_WORLD, &p)); 309566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(p, PETSCPARTITIONERSIMPLE)); 319566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(p)); 32abe9303eSLisandro Dalcin 33d7cc930eSLisandro Dalcin /* create partition section */ 349566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PETSC_COMM_WORLD, &partSection)); 35d7cc930eSLisandro Dalcin 36d7cc930eSLisandro Dalcin if (vwgts) { /* create vertex weights section */ 379566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PETSC_COMM_WORLD, &vertexSection)); 389566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(vertexSection, 0, nv)); 399566063dSJacob Faibussowitsch for (i = 0; i < nv; i++) PetscCall(PetscSectionSetDof(vertexSection, i, 1)); 409566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(vertexSection)); 41d7cc930eSLisandro Dalcin } 42d7cc930eSLisandro Dalcin 43d7cc930eSLisandro Dalcin if (pwgts) { /* create partition weights section */ 449566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PETSC_COMM_WORLD, &targetSection)); 459566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(targetSection, 0, nparts)); 469566063dSJacob Faibussowitsch for (i = 0; i < nparts; i++) PetscCall(PetscSectionSetDof(targetSection, i, 1)); 479566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(targetSection)); 48d7cc930eSLisandro Dalcin } 49d7cc930eSLisandro Dalcin 502611ad71SToby Isaac if (PetscDefined(USE_LOG)) { /* Test logging */ 51d7cc930eSLisandro Dalcin PetscLogEvent event; 52d7cc930eSLisandro Dalcin 539566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("MyPartitionerEvent", PETSCPARTITIONER_CLASSID, &event)); 54d7cc930eSLisandro Dalcin { /* PetscLogEventExcludeClass is broken, new events are not deactivated */ 55d7cc930eSLisandro Dalcin char logList[256]; 56d7cc930eSLisandro Dalcin PetscBool opt, pkg; 57d7cc930eSLisandro Dalcin 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-log_exclude", logList, sizeof(logList), &opt)); 59d7cc930eSLisandro Dalcin if (opt) { 609566063dSJacob Faibussowitsch PetscCall(PetscStrInList("partitioner", logList, ',', &pkg)); 619566063dSJacob Faibussowitsch if (pkg) PetscCall(PetscLogEventExcludeClass(PETSCPARTITIONER_CLASSID)); 62d7cc930eSLisandro Dalcin } 63d7cc930eSLisandro Dalcin } 649566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(event, p, NULL, NULL, NULL)); 659566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(event, p, NULL, NULL, NULL)); 66d7cc930eSLisandro Dalcin } 67d7cc930eSLisandro Dalcin 68d7cc930eSLisandro Dalcin /* test setup and reset */ 699566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetUp(p)); 709566063dSJacob Faibussowitsch PetscCall(PetscPartitionerReset(p)); 71d7cc930eSLisandro Dalcin 72abe9303eSLisandro Dalcin /* test partitioning an empty graph */ 7321c92275SMatthew G. Knepley PetscCall(PetscPartitionerPartition(p, nparts, 0, NULL, NULL, vertexSection, NULL, targetSection, partSection, &partition)); 749566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)partSection, "NULL SECTION")); 759566063dSJacob Faibussowitsch PetscCall(PetscSectionView(partSection, NULL)); 769566063dSJacob Faibussowitsch PetscCall(ISOnComm(partition, PETSC_COMM_WORLD, PETSC_USE_POINTER, &is)); 779566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)is, "NULL PARTITION")); 789566063dSJacob Faibussowitsch PetscCall(ISView(is, NULL)); 799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 809566063dSJacob Faibussowitsch PetscCall(ISDestroy(&partition)); 81abe9303eSLisandro Dalcin 82d7cc930eSLisandro Dalcin /* test view from options */ 839566063dSJacob Faibussowitsch PetscCall(PetscPartitionerViewFromOptions(p, NULL, "-part_view")); 84d7cc930eSLisandro Dalcin 859dddd249SSatish Balay /* test partitioning a graph on one process only (not main) */ 86abe9303eSLisandro Dalcin if (rank == size - 1) { 8721c92275SMatthew G. Knepley PetscCall(PetscPartitionerPartition(p, nparts, nv, vv, vadj, vertexSection, NULL, targetSection, partSection, &partition)); 88abe9303eSLisandro Dalcin } else { 8921c92275SMatthew G. Knepley PetscCall(PetscPartitionerPartition(p, nparts, 0, NULL, NULL, vertexSection, NULL, targetSection, partSection, &partition)); 90abe9303eSLisandro Dalcin } 919566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)partSection, "SEQ SECTION")); 929566063dSJacob Faibussowitsch PetscCall(PetscSectionView(partSection, NULL)); 939566063dSJacob Faibussowitsch PetscCall(ISOnComm(partition, PETSC_COMM_WORLD, PETSC_USE_POINTER, &is)); 949566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)is, "SEQ PARTITION")); 959566063dSJacob Faibussowitsch PetscCall(ISView(is, NULL)); 969566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 979566063dSJacob Faibussowitsch PetscCall(ISDestroy(&partition)); 98abe9303eSLisandro Dalcin 999566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)p, &sequential, PETSCPARTITIONERCHACO, NULL)); 10056991f02SLisandro Dalcin if (sequential) goto finally; 10156991f02SLisandro Dalcin 102da81f932SPierre Jolivet /* test partitioning a graph on a subset of the processes only */ 103abe9303eSLisandro Dalcin if (rank % 2) { 10421c92275SMatthew G. Knepley PetscCall(PetscPartitionerPartition(p, nparts, 0, NULL, NULL, NULL, NULL, targetSection, partSection, &partition)); 105abe9303eSLisandro Dalcin } else { 106abe9303eSLisandro Dalcin PetscInt i, totv = nv * ((size + 1) / 2), *pvadj; 107abe9303eSLisandro Dalcin 1089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * nv, &pvadj)); 109abe9303eSLisandro Dalcin for (i = 0; i < nv; i++) { 110abe9303eSLisandro Dalcin pvadj[2 * i] = (nv * (rank / 2) + totv + i - 1) % totv; 111abe9303eSLisandro Dalcin pvadj[2 * i + 1] = (nv * (rank / 2) + totv + i + 1) % totv; 112abe9303eSLisandro Dalcin } 11321c92275SMatthew G. Knepley PetscCall(PetscPartitionerPartition(p, nparts, nv, vv, pvadj, NULL, NULL, targetSection, partSection, &partition)); 1149566063dSJacob Faibussowitsch PetscCall(PetscFree(pvadj)); 115abe9303eSLisandro Dalcin } 1169566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)partSection, "PARVOID SECTION")); 1179566063dSJacob Faibussowitsch PetscCall(PetscSectionView(partSection, NULL)); 1189566063dSJacob Faibussowitsch PetscCall(ISOnComm(partition, PETSC_COMM_WORLD, PETSC_USE_POINTER, &is)); 1199566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)is, "PARVOID PARTITION")); 1209566063dSJacob Faibussowitsch PetscCall(ISView(is, NULL)); 1219566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 1229566063dSJacob Faibussowitsch PetscCall(ISDestroy(&partition)); 123abe9303eSLisandro Dalcin 12456991f02SLisandro Dalcin finally: 1259566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&partSection)); 1269566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&vertexSection)); 1279566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&targetSection)); 1289566063dSJacob Faibussowitsch PetscCall(PetscPartitionerDestroy(&p)); 1299566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 130b122ec5aSJacob Faibussowitsch return 0; 131abe9303eSLisandro Dalcin } 132abe9303eSLisandro Dalcin 133abe9303eSLisandro Dalcin /*TEST 134abe9303eSLisandro Dalcin 135abe9303eSLisandro Dalcin test: 136d7cc930eSLisandro Dalcin suffix: default 137d7cc930eSLisandro Dalcin 138d7cc930eSLisandro Dalcin testset: 139dfd57a17SPierre Jolivet requires: defined(PETSC_USE_LOG) 140d75802c7SJacob Faibussowitsch args: -petscpartitioner_type simple -log_view 14179229e12SLisandro Dalcin filter: grep MyPartitionerEvent | cut -d " " -f 1 142d7cc930eSLisandro Dalcin test: 143d7cc930eSLisandro Dalcin suffix: log_include 144d7cc930eSLisandro Dalcin test: 145d7cc930eSLisandro Dalcin suffix: log_exclude 146d7cc930eSLisandro Dalcin args: -log_exclude partitioner 147*e0008caeSPierre Jolivet output_file: output/empty.out 148d7cc930eSLisandro Dalcin 149d7cc930eSLisandro Dalcin test: 150abe9303eSLisandro Dalcin suffix: simple 151abe9303eSLisandro Dalcin nsize: {{1 2 3}separate output} 152d7cc930eSLisandro Dalcin args: -nparts {{1 2 3}separate output} -pwgts {{false true}separate output} -petscpartitioner_type simple -petscpartitioner_view 153d7cc930eSLisandro Dalcin 154d7cc930eSLisandro Dalcin test: 155d7cc930eSLisandro Dalcin suffix: shell 156d7cc930eSLisandro Dalcin nsize: {{1 2 3}separate output} 157d7cc930eSLisandro Dalcin args: -nparts {{1 2 3}separate output} -petscpartitioner_type shell -petscpartitioner_shell_random -petscpartitioner_view 158abe9303eSLisandro Dalcin 159abe9303eSLisandro Dalcin test: 1605130b6c4SLisandro Dalcin suffix: gather 1615130b6c4SLisandro Dalcin nsize: {{1 2 3}separate output} 1625130b6c4SLisandro Dalcin args: -nparts {{1 2 3}separate output} -petscpartitioner_type gather -petscpartitioner_view -petscpartitioner_view_graph 1635130b6c4SLisandro Dalcin 1645130b6c4SLisandro Dalcin test: 165abe9303eSLisandro Dalcin requires: parmetis 166abe9303eSLisandro Dalcin suffix: parmetis 167abe9303eSLisandro Dalcin nsize: {{1 2 3}separate output} 168d7cc930eSLisandro Dalcin args: -nparts {{1 2 3}separate output} -pwgts {{false true}} -vwgts {{false true}} 169d7cc930eSLisandro Dalcin args: -petscpartitioner_type parmetis -petscpartitioner_view -petscpartitioner_view_graph 170d7cc930eSLisandro Dalcin 171d7cc930eSLisandro Dalcin test: 172d7cc930eSLisandro Dalcin requires: parmetis 173d7cc930eSLisandro Dalcin suffix: parmetis_type 174d7cc930eSLisandro Dalcin nsize: {{1 2}} 175d7cc930eSLisandro Dalcin args: -petscpartitioner_type parmetis -part_view 176d7cc930eSLisandro Dalcin args: -petscpartitioner_parmetis_type {{kway rb}separate output} 177d7cc930eSLisandro Dalcin filter: grep "ParMetis type" 178abe9303eSLisandro Dalcin 179abe9303eSLisandro Dalcin test: 180abe9303eSLisandro Dalcin requires: ptscotch 181abe9303eSLisandro Dalcin suffix: ptscotch 182abe9303eSLisandro Dalcin nsize: {{1 2 3}separate output} 183d7cc930eSLisandro Dalcin args: -nparts {{1 2 3}separate output} -pwgts {{false true}separate output} -vwgts {{false true}} 184d7cc930eSLisandro Dalcin args: -petscpartitioner_type ptscotch -petscpartitioner_view -petscpartitioner_view_graph 185abe9303eSLisandro Dalcin 18656991f02SLisandro Dalcin test: 187d7cc930eSLisandro Dalcin requires: ptscotch 188d7cc930eSLisandro Dalcin suffix: ptscotch_strategy 189d7cc930eSLisandro Dalcin nsize: {{1 2}} 190d7cc930eSLisandro Dalcin args: -petscpartitioner_type ptscotch -part_view 191d7cc930eSLisandro Dalcin args: -petscpartitioner_ptscotch_strategy {{DEFAULT QUALITY SPEED BALANCE SAFETY SCALABILITY RECURSIVE REMAP}separate output} 192d7cc930eSLisandro Dalcin filter: grep "partitioning strategy" 193d7cc930eSLisandro Dalcin 194d7cc930eSLisandro Dalcin test: 19556991f02SLisandro Dalcin requires: chaco 19656991f02SLisandro Dalcin suffix: chaco 19756991f02SLisandro Dalcin nsize: {{1 2 3}separate output} 198d7cc930eSLisandro Dalcin args: -nparts {{1}separate output} -petscpartitioner_type chaco -petscpartitioner_view -petscpartitioner_view_graph 199d7cc930eSLisandro Dalcin 200d7cc930eSLisandro Dalcin test: 201d7cc930eSLisandro Dalcin TODO: non reproducible (uses C stdlib rand()) 202d7cc930eSLisandro Dalcin requires: chaco 203d7cc930eSLisandro Dalcin suffix: chaco 204d7cc930eSLisandro Dalcin nsize: {{1 2 3}separate output} 205d7cc930eSLisandro Dalcin args: -nparts {{2 3}separate output} -petscpartitioner_type chaco -petscpartitioner_view -petscpartitioner_view_graph 20656991f02SLisandro Dalcin 207abe9303eSLisandro Dalcin TEST*/ 208