1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Partition a tiny grid using hierarchical partitioning.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices. Note that this file 6c4762a1bSJed Brown automatically includes: 7c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 8c4762a1bSJed Brown petscmat.h - matrices 9c4762a1bSJed Brown petscis.h - index sets 10c4762a1bSJed Brown petscviewer.h - viewers 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown #include <petscmat.h> 13c4762a1bSJed Brown 14*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 15*d71ae5a4SJacob Faibussowitsch { 16c4762a1bSJed Brown Mat A; 17c4762a1bSJed Brown PetscMPIInt rank, size; 18c4762a1bSJed Brown PetscInt *ia, *ja; 19c4762a1bSJed Brown MatPartitioning part; 20c4762a1bSJed Brown IS is, isn, isrows; 21c4762a1bSJed Brown IS coarseparts, fineparts; 22c4762a1bSJed Brown MPI_Comm comm; 23c4762a1bSJed Brown 24327415f7SBarry Smith PetscFunctionBeginUser; 259566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 26c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 28be096a46SBarry Smith PetscCheck(size == 4, comm, PETSC_ERR_WRONG_MPI_SIZE, "Must run with 4 processors"); 299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 30c4762a1bSJed Brown 319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5, &ia)); 329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(16, &ja)); 33dd400576SPatrick Sanan if (rank == 0) { 349371c9d4SSatish Balay ja[0] = 1; 359371c9d4SSatish Balay ja[1] = 4; 369371c9d4SSatish Balay ja[2] = 0; 379371c9d4SSatish Balay ja[3] = 2; 389371c9d4SSatish Balay ja[4] = 5; 399371c9d4SSatish Balay ja[5] = 1; 409371c9d4SSatish Balay ja[6] = 3; 419371c9d4SSatish Balay ja[7] = 6; 429371c9d4SSatish Balay ja[8] = 2; 439371c9d4SSatish Balay ja[9] = 7; 449371c9d4SSatish Balay ia[0] = 0; 459371c9d4SSatish Balay ia[1] = 2; 469371c9d4SSatish Balay ia[2] = 5; 479371c9d4SSatish Balay ia[3] = 8; 489371c9d4SSatish Balay ia[4] = 10; 49c4762a1bSJed Brown } else if (rank == 1) { 509371c9d4SSatish Balay ja[0] = 0; 519371c9d4SSatish Balay ja[1] = 5; 529371c9d4SSatish Balay ja[2] = 8; 539371c9d4SSatish Balay ja[3] = 1; 549371c9d4SSatish Balay ja[4] = 4; 559371c9d4SSatish Balay ja[5] = 6; 569371c9d4SSatish Balay ja[6] = 9; 579371c9d4SSatish Balay ja[7] = 2; 589371c9d4SSatish Balay ja[8] = 5; 599371c9d4SSatish Balay ja[9] = 7; 609371c9d4SSatish Balay ja[10] = 10; 619371c9d4SSatish Balay ja[11] = 3; 629371c9d4SSatish Balay ja[12] = 6; 639371c9d4SSatish Balay ja[13] = 11; 649371c9d4SSatish Balay ia[0] = 0; 659371c9d4SSatish Balay ia[1] = 3; 669371c9d4SSatish Balay ia[2] = 7; 679371c9d4SSatish Balay ia[3] = 11; 689371c9d4SSatish Balay ia[4] = 14; 69c4762a1bSJed Brown } else if (rank == 2) { 709371c9d4SSatish Balay ja[0] = 4; 719371c9d4SSatish Balay ja[1] = 9; 729371c9d4SSatish Balay ja[2] = 12; 739371c9d4SSatish Balay ja[3] = 5; 749371c9d4SSatish Balay ja[4] = 8; 759371c9d4SSatish Balay ja[5] = 10; 769371c9d4SSatish Balay ja[6] = 13; 779371c9d4SSatish Balay ja[7] = 6; 789371c9d4SSatish Balay ja[8] = 9; 799371c9d4SSatish Balay ja[9] = 11; 809371c9d4SSatish Balay ja[10] = 14; 819371c9d4SSatish Balay ja[11] = 7; 829371c9d4SSatish Balay ja[12] = 10; 839371c9d4SSatish Balay ja[13] = 15; 849371c9d4SSatish Balay ia[0] = 0; 859371c9d4SSatish Balay ia[1] = 3; 869371c9d4SSatish Balay ia[2] = 7; 879371c9d4SSatish Balay ia[3] = 11; 889371c9d4SSatish Balay ia[4] = 14; 89c4762a1bSJed Brown } else { 909371c9d4SSatish Balay ja[0] = 8; 919371c9d4SSatish Balay ja[1] = 13; 929371c9d4SSatish Balay ja[2] = 9; 939371c9d4SSatish Balay ja[3] = 12; 949371c9d4SSatish Balay ja[4] = 14; 959371c9d4SSatish Balay ja[5] = 10; 969371c9d4SSatish Balay ja[6] = 13; 979371c9d4SSatish Balay ja[7] = 15; 989371c9d4SSatish Balay ja[8] = 11; 999371c9d4SSatish Balay ja[9] = 14; 1009371c9d4SSatish Balay ia[0] = 0; 1019371c9d4SSatish Balay ia[1] = 2; 1029371c9d4SSatish Balay ia[2] = 5; 1039371c9d4SSatish Balay ia[3] = 8; 1049371c9d4SSatish Balay ia[4] = 10; 105c4762a1bSJed Brown } 1069566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(comm, 4, 16, ia, ja, NULL, &A)); 1079566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 108c4762a1bSJed Brown /* 109c4762a1bSJed Brown Partition the graph of the matrix 110c4762a1bSJed Brown */ 1119566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(comm, &part)); 1129566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(part, A)); 1139566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetType(part, MATPARTITIONINGHIERARCH)); 1149566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNcoarseparts(part, 2)); 1159566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNfineparts(part, 2)); 1169566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(part)); 117c4762a1bSJed Brown /* get new processor owner number of each vertex */ 1189566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(part, &is)); 119c4762a1bSJed Brown /* coarse parts */ 1209566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalGetCoarseparts(part, &coarseparts)); 1219566063dSJacob Faibussowitsch PetscCall(ISView(coarseparts, PETSC_VIEWER_STDOUT_WORLD)); 122c4762a1bSJed Brown /* fine parts */ 1239566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalGetFineparts(part, &fineparts)); 1249566063dSJacob Faibussowitsch PetscCall(ISView(fineparts, PETSC_VIEWER_STDOUT_WORLD)); 125c4762a1bSJed Brown /* partitioning */ 1269566063dSJacob Faibussowitsch PetscCall(ISView(is, PETSC_VIEWER_STDOUT_WORLD)); 127c4762a1bSJed Brown /* get new global number of each old global number */ 1289566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(is, &isn)); 1299566063dSJacob Faibussowitsch PetscCall(ISView(isn, PETSC_VIEWER_STDOUT_WORLD)); 1309566063dSJacob Faibussowitsch PetscCall(ISBuildTwoSided(is, NULL, &isrows)); 1319566063dSJacob Faibussowitsch PetscCall(ISView(isrows, PETSC_VIEWER_STDOUT_WORLD)); 1329566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 1339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&coarseparts)); 1349566063dSJacob Faibussowitsch PetscCall(ISDestroy(&fineparts)); 1359566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 1369566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isn)); 1379566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&part)); 138c4762a1bSJed Brown /* 139c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 140c4762a1bSJed Brown are no longer needed. 141c4762a1bSJed Brown */ 1429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1439566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 144b122ec5aSJacob Faibussowitsch return 0; 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /*TEST 148c4762a1bSJed Brown 149c4762a1bSJed Brown test: 150c4762a1bSJed Brown nsize: 4 151c4762a1bSJed Brown requires: parmetis 152c4762a1bSJed Brown TODO: tests cannot use parmetis because it produces different results on different machines 153c4762a1bSJed Brown 154c4762a1bSJed Brown TEST*/ 155