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 14c4762a1bSJed Brown int main(int argc,char **args) 15c4762a1bSJed Brown { 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 24*327415f7SBarry 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) { 34c4762a1bSJed Brown ja[0] = 1; ja[1] = 4; ja[2] = 0; ja[3] = 2; ja[4] = 5; ja[5] = 1; ja[6] = 3; ja[7] = 6; 35c4762a1bSJed Brown ja[8] = 2; ja[9] = 7; 36c4762a1bSJed Brown ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 37c4762a1bSJed Brown } else if (rank == 1) { 38c4762a1bSJed Brown ja[0] = 0; ja[1] = 5; ja[2] = 8; ja[3] = 1; ja[4] = 4; ja[5] = 6; ja[6] = 9; ja[7] = 2; 39c4762a1bSJed Brown ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11; 40c4762a1bSJed Brown ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 41c4762a1bSJed Brown } else if (rank == 2) { 42c4762a1bSJed Brown ja[0] = 4; ja[1] = 9; ja[2] = 12; ja[3] = 5; ja[4] = 8; ja[5] = 10; ja[6] = 13; ja[7] = 6; 43c4762a1bSJed Brown ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15; 44c4762a1bSJed Brown ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 45c4762a1bSJed Brown } else { 46c4762a1bSJed Brown ja[0] = 8; ja[1] = 13; ja[2] = 9; ja[3] = 12; ja[4] = 14; ja[5] = 10; ja[6] = 13; ja[7] = 15; 47c4762a1bSJed Brown ja[8] = 11; ja[9] = 14; 48c4762a1bSJed Brown ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 49c4762a1bSJed Brown } 509566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(comm,4,16,ia,ja,NULL,&A)); 519566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 52c4762a1bSJed Brown /* 53c4762a1bSJed Brown Partition the graph of the matrix 54c4762a1bSJed Brown */ 559566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(comm,&part)); 569566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(part,A)); 579566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetType(part,MATPARTITIONINGHIERARCH)); 589566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNcoarseparts(part,2)); 599566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNfineparts(part,2)); 609566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(part)); 61c4762a1bSJed Brown /* get new processor owner number of each vertex */ 629566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(part,&is)); 63c4762a1bSJed Brown /* coarse parts */ 649566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalGetCoarseparts(part,&coarseparts)); 659566063dSJacob Faibussowitsch PetscCall(ISView(coarseparts,PETSC_VIEWER_STDOUT_WORLD)); 66c4762a1bSJed Brown /* fine parts */ 679566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalGetFineparts(part,&fineparts)); 689566063dSJacob Faibussowitsch PetscCall(ISView(fineparts,PETSC_VIEWER_STDOUT_WORLD)); 69c4762a1bSJed Brown /* partitioning */ 709566063dSJacob Faibussowitsch PetscCall(ISView(is,PETSC_VIEWER_STDOUT_WORLD)); 71c4762a1bSJed Brown /* get new global number of each old global number */ 729566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(is,&isn)); 739566063dSJacob Faibussowitsch PetscCall(ISView(isn,PETSC_VIEWER_STDOUT_WORLD)); 749566063dSJacob Faibussowitsch PetscCall(ISBuildTwoSided(is,NULL,&isrows)); 759566063dSJacob Faibussowitsch PetscCall(ISView(isrows,PETSC_VIEWER_STDOUT_WORLD)); 769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&coarseparts)); 789566063dSJacob Faibussowitsch PetscCall(ISDestroy(&fineparts)); 799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 809566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isn)); 819566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&part)); 82c4762a1bSJed Brown /* 83c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 84c4762a1bSJed Brown are no longer needed. 85c4762a1bSJed Brown */ 869566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 879566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 88b122ec5aSJacob Faibussowitsch return 0; 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91c4762a1bSJed Brown /*TEST 92c4762a1bSJed Brown 93c4762a1bSJed Brown test: 94c4762a1bSJed Brown nsize: 4 95c4762a1bSJed Brown requires: parmetis 96c4762a1bSJed Brown TODO: tests cannot use parmetis because it produces different results on different machines 97c4762a1bSJed Brown 98c4762a1bSJed Brown TEST*/ 99