1c4762a1bSJed Brown /* 2c4762a1bSJed Brown * ex193.c 3c4762a1bSJed Brown * 4c4762a1bSJed Brown * Created on: Jul 29, 2015 5c4762a1bSJed Brown * Author: Fande Kong fdkong.jd@gmail.com 6c4762a1bSJed Brown */ 7c4762a1bSJed Brown /* 8c4762a1bSJed Brown * An example demonstrates how to use hierarchical partitioning approach 9c4762a1bSJed Brown */ 10c4762a1bSJed Brown 11c4762a1bSJed Brown #include <petscmat.h> 12c4762a1bSJed Brown 13c4762a1bSJed Brown static char help[] = "Illustrates use of hierarchical partitioning.\n"; 14c4762a1bSJed Brown 15c4762a1bSJed Brown int main(int argc,char **args) 16c4762a1bSJed Brown { 17c4762a1bSJed Brown Mat A; /* matrix */ 18c4762a1bSJed Brown PetscInt m,n; /* mesh dimensions in x- and y- directions */ 19c4762a1bSJed Brown PetscInt i,j,Ii,J,Istart,Iend; 20c4762a1bSJed Brown PetscMPIInt size; 21c4762a1bSJed Brown PetscScalar v; 22c4762a1bSJed Brown MatPartitioning part; 23c4762a1bSJed Brown IS coarseparts,fineparts; 24c4762a1bSJed Brown IS is,isn,isrows; 25c4762a1bSJed Brown MPI_Comm comm; 26c4762a1bSJed Brown 27*327415f7SBarry Smith PetscFunctionBeginUser; 289566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 29c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 31d0609cedSBarry Smith PetscOptionsBegin(comm,NULL,"ex193","hierarchical partitioning"); 32c4762a1bSJed Brown m = 15; 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-M","Number of mesh points in the x-direction","partitioning",m,&m,NULL)); 34c4762a1bSJed Brown n = 15; 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-N","Number of mesh points in the y-direction","partitioning",n,&n,NULL)); 36d0609cedSBarry Smith PetscOptionsEnd(); 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* 39c4762a1bSJed Brown Assemble the matrix for the five point stencil (finite difference), YET AGAIN 40c4762a1bSJed Brown */ 419566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&A)); 429566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 439566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 449566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 459566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); 46c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 47c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 489566063dSJacob Faibussowitsch if (i>0) {J = Ii - n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 499566063dSJacob Faibussowitsch if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 509566063dSJacob Faibussowitsch if (j>0) {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 519566063dSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 529566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 53c4762a1bSJed Brown } 549566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 559566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 569566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 57c4762a1bSJed Brown /* 58c4762a1bSJed Brown Partition the graph of the matrix 59c4762a1bSJed Brown */ 609566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(comm,&part)); 619566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(part,A)); 629566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetType(part,MATPARTITIONINGHIERARCH)); 639566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNcoarseparts(part,2)); 649566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNfineparts(part,4)); 659566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(part)); 66c4762a1bSJed Brown /* get new processor owner number of each vertex */ 679566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(part,&is)); 68c4762a1bSJed Brown /* coarse parts */ 699566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalGetCoarseparts(part,&coarseparts)); 709566063dSJacob Faibussowitsch PetscCall(ISView(coarseparts,PETSC_VIEWER_STDOUT_WORLD)); 71c4762a1bSJed Brown /* fine parts */ 729566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalGetFineparts(part,&fineparts)); 739566063dSJacob Faibussowitsch PetscCall(ISView(fineparts,PETSC_VIEWER_STDOUT_WORLD)); 74c4762a1bSJed Brown /* partitioning */ 759566063dSJacob Faibussowitsch PetscCall(ISView(is,PETSC_VIEWER_STDOUT_WORLD)); 76c4762a1bSJed Brown /* get new global number of each old global number */ 779566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(is,&isn)); 789566063dSJacob Faibussowitsch PetscCall(ISView(isn,PETSC_VIEWER_STDOUT_WORLD)); 799566063dSJacob Faibussowitsch PetscCall(ISBuildTwoSided(is,NULL,&isrows)); 809566063dSJacob Faibussowitsch PetscCall(ISView(isrows,PETSC_VIEWER_STDOUT_WORLD)); 819566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&coarseparts)); 839566063dSJacob Faibussowitsch PetscCall(ISDestroy(&fineparts)); 849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 859566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isn)); 869566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&part)); 879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 889566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 89b122ec5aSJacob Faibussowitsch return 0; 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown /*TEST 93c4762a1bSJed Brown 94c4762a1bSJed Brown test: 95c4762a1bSJed Brown nsize: 4 96c4762a1bSJed Brown args: -mat_partitioning_hierarchical_Nfineparts 2 97c4762a1bSJed Brown requires: parmetis 98c4762a1bSJed Brown TODO: cannot run because parmetis does reproduce across all machines, probably due to nonportable random number generator 99c4762a1bSJed Brown 100c4762a1bSJed Brown TEST*/ 101