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*9371c9d4SSatish Balay int main(int argc, char **args) { 15c4762a1bSJed Brown Mat A; 16c4762a1bSJed Brown PetscMPIInt rank, size; 17c4762a1bSJed Brown PetscInt *ia, *ja; 18c4762a1bSJed Brown MatPartitioning part; 19c4762a1bSJed Brown IS is, isn, isrows; 20c4762a1bSJed Brown IS coarseparts, fineparts; 21c4762a1bSJed Brown MPI_Comm comm; 22c4762a1bSJed Brown 23327415f7SBarry Smith PetscFunctionBeginUser; 249566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 25c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 27be096a46SBarry Smith PetscCheck(size == 4, comm, PETSC_ERR_WRONG_MPI_SIZE, "Must run with 4 processors"); 289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 29c4762a1bSJed Brown 309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5, &ia)); 319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(16, &ja)); 32dd400576SPatrick Sanan if (rank == 0) { 33*9371c9d4SSatish Balay ja[0] = 1; 34*9371c9d4SSatish Balay ja[1] = 4; 35*9371c9d4SSatish Balay ja[2] = 0; 36*9371c9d4SSatish Balay ja[3] = 2; 37*9371c9d4SSatish Balay ja[4] = 5; 38*9371c9d4SSatish Balay ja[5] = 1; 39*9371c9d4SSatish Balay ja[6] = 3; 40*9371c9d4SSatish Balay ja[7] = 6; 41*9371c9d4SSatish Balay ja[8] = 2; 42*9371c9d4SSatish Balay ja[9] = 7; 43*9371c9d4SSatish Balay ia[0] = 0; 44*9371c9d4SSatish Balay ia[1] = 2; 45*9371c9d4SSatish Balay ia[2] = 5; 46*9371c9d4SSatish Balay ia[3] = 8; 47*9371c9d4SSatish Balay ia[4] = 10; 48c4762a1bSJed Brown } else if (rank == 1) { 49*9371c9d4SSatish Balay ja[0] = 0; 50*9371c9d4SSatish Balay ja[1] = 5; 51*9371c9d4SSatish Balay ja[2] = 8; 52*9371c9d4SSatish Balay ja[3] = 1; 53*9371c9d4SSatish Balay ja[4] = 4; 54*9371c9d4SSatish Balay ja[5] = 6; 55*9371c9d4SSatish Balay ja[6] = 9; 56*9371c9d4SSatish Balay ja[7] = 2; 57*9371c9d4SSatish Balay ja[8] = 5; 58*9371c9d4SSatish Balay ja[9] = 7; 59*9371c9d4SSatish Balay ja[10] = 10; 60*9371c9d4SSatish Balay ja[11] = 3; 61*9371c9d4SSatish Balay ja[12] = 6; 62*9371c9d4SSatish Balay ja[13] = 11; 63*9371c9d4SSatish Balay ia[0] = 0; 64*9371c9d4SSatish Balay ia[1] = 3; 65*9371c9d4SSatish Balay ia[2] = 7; 66*9371c9d4SSatish Balay ia[3] = 11; 67*9371c9d4SSatish Balay ia[4] = 14; 68c4762a1bSJed Brown } else if (rank == 2) { 69*9371c9d4SSatish Balay ja[0] = 4; 70*9371c9d4SSatish Balay ja[1] = 9; 71*9371c9d4SSatish Balay ja[2] = 12; 72*9371c9d4SSatish Balay ja[3] = 5; 73*9371c9d4SSatish Balay ja[4] = 8; 74*9371c9d4SSatish Balay ja[5] = 10; 75*9371c9d4SSatish Balay ja[6] = 13; 76*9371c9d4SSatish Balay ja[7] = 6; 77*9371c9d4SSatish Balay ja[8] = 9; 78*9371c9d4SSatish Balay ja[9] = 11; 79*9371c9d4SSatish Balay ja[10] = 14; 80*9371c9d4SSatish Balay ja[11] = 7; 81*9371c9d4SSatish Balay ja[12] = 10; 82*9371c9d4SSatish Balay ja[13] = 15; 83*9371c9d4SSatish Balay ia[0] = 0; 84*9371c9d4SSatish Balay ia[1] = 3; 85*9371c9d4SSatish Balay ia[2] = 7; 86*9371c9d4SSatish Balay ia[3] = 11; 87*9371c9d4SSatish Balay ia[4] = 14; 88c4762a1bSJed Brown } else { 89*9371c9d4SSatish Balay ja[0] = 8; 90*9371c9d4SSatish Balay ja[1] = 13; 91*9371c9d4SSatish Balay ja[2] = 9; 92*9371c9d4SSatish Balay ja[3] = 12; 93*9371c9d4SSatish Balay ja[4] = 14; 94*9371c9d4SSatish Balay ja[5] = 10; 95*9371c9d4SSatish Balay ja[6] = 13; 96*9371c9d4SSatish Balay ja[7] = 15; 97*9371c9d4SSatish Balay ja[8] = 11; 98*9371c9d4SSatish Balay ja[9] = 14; 99*9371c9d4SSatish Balay ia[0] = 0; 100*9371c9d4SSatish Balay ia[1] = 2; 101*9371c9d4SSatish Balay ia[2] = 5; 102*9371c9d4SSatish Balay ia[3] = 8; 103*9371c9d4SSatish Balay ia[4] = 10; 104c4762a1bSJed Brown } 1059566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(comm, 4, 16, ia, ja, NULL, &A)); 1069566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 107c4762a1bSJed Brown /* 108c4762a1bSJed Brown Partition the graph of the matrix 109c4762a1bSJed Brown */ 1109566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(comm, &part)); 1119566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(part, A)); 1129566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetType(part, MATPARTITIONINGHIERARCH)); 1139566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNcoarseparts(part, 2)); 1149566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNfineparts(part, 2)); 1159566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(part)); 116c4762a1bSJed Brown /* get new processor owner number of each vertex */ 1179566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(part, &is)); 118c4762a1bSJed Brown /* coarse parts */ 1199566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalGetCoarseparts(part, &coarseparts)); 1209566063dSJacob Faibussowitsch PetscCall(ISView(coarseparts, PETSC_VIEWER_STDOUT_WORLD)); 121c4762a1bSJed Brown /* fine parts */ 1229566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalGetFineparts(part, &fineparts)); 1239566063dSJacob Faibussowitsch PetscCall(ISView(fineparts, PETSC_VIEWER_STDOUT_WORLD)); 124c4762a1bSJed Brown /* partitioning */ 1259566063dSJacob Faibussowitsch PetscCall(ISView(is, PETSC_VIEWER_STDOUT_WORLD)); 126c4762a1bSJed Brown /* get new global number of each old global number */ 1279566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(is, &isn)); 1289566063dSJacob Faibussowitsch PetscCall(ISView(isn, PETSC_VIEWER_STDOUT_WORLD)); 1299566063dSJacob Faibussowitsch PetscCall(ISBuildTwoSided(is, NULL, &isrows)); 1309566063dSJacob Faibussowitsch PetscCall(ISView(isrows, PETSC_VIEWER_STDOUT_WORLD)); 1319566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 1329566063dSJacob Faibussowitsch PetscCall(ISDestroy(&coarseparts)); 1339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&fineparts)); 1349566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 1359566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isn)); 1369566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&part)); 137c4762a1bSJed Brown /* 138c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 139c4762a1bSJed Brown are no longer needed. 140c4762a1bSJed Brown */ 1419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1429566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 143b122ec5aSJacob Faibussowitsch return 0; 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown /*TEST 147c4762a1bSJed Brown 148c4762a1bSJed Brown test: 149c4762a1bSJed Brown nsize: 4 150c4762a1bSJed Brown requires: parmetis 151c4762a1bSJed Brown TODO: tests cannot use parmetis because it produces different results on different machines 152c4762a1bSJed Brown 153c4762a1bSJed Brown TEST*/ 154