1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Partition tiny grid using hierarchical partitioning and increase overlap using MatIncreaseOverlapSplit.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*T 5c4762a1bSJed Brown Concepts: partitioning 6c4762a1bSJed Brown Processors: 4 7c4762a1bSJed Brown T*/ 8c4762a1bSJed Brown 9c4762a1bSJed Brown /* 10c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices. Note that this file 11c4762a1bSJed Brown automatically includes: 12c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 13c4762a1bSJed Brown petscmat.h - matrices 14c4762a1bSJed Brown petscis.h - index sets 15c4762a1bSJed Brown petscviewer.h - viewers 16c4762a1bSJed Brown */ 17c4762a1bSJed Brown #include <petscmat.h> 18c4762a1bSJed Brown 19c4762a1bSJed Brown int main(int argc,char **args) 20c4762a1bSJed Brown { 21c4762a1bSJed Brown Mat A,B; 22c4762a1bSJed Brown PetscMPIInt rank,size,membershipKey; 23c4762a1bSJed Brown PetscInt *ia,*ja,*indices_sc,isrows_localsize; 24c4762a1bSJed Brown const PetscInt *indices; 25c4762a1bSJed Brown MatPartitioning part; 26c4762a1bSJed Brown IS is,isrows,isrows_sc; 27c4762a1bSJed Brown IS coarseparts,fineparts; 28c4762a1bSJed Brown MPI_Comm comm,scomm; 29c4762a1bSJed Brown 30*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 31c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 325f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm,&size)); 332c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 4,comm,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 4 processors "); 345f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm,&rank)); 35c4762a1bSJed Brown /*set a small matrix */ 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(5,&ia)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(16,&ja)); 38dd400576SPatrick Sanan if (rank == 0) { 39c4762a1bSJed 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; 40c4762a1bSJed Brown ja[8] = 2; ja[9] = 7; 41c4762a1bSJed Brown ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 42c4762a1bSJed Brown membershipKey = 0; 43c4762a1bSJed Brown } else if (rank == 1) { 44c4762a1bSJed 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; 45c4762a1bSJed Brown ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11; 46c4762a1bSJed Brown ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 47c4762a1bSJed Brown membershipKey = 0; 48c4762a1bSJed Brown } else if (rank == 2) { 49c4762a1bSJed 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; 50c4762a1bSJed Brown ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15; 51c4762a1bSJed Brown ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 52c4762a1bSJed Brown membershipKey = 1; 53c4762a1bSJed Brown } else { 54c4762a1bSJed 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; 55c4762a1bSJed Brown ja[8] = 11; ja[9] = 14; 56c4762a1bSJed Brown ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 57c4762a1bSJed Brown membershipKey = 1; 58c4762a1bSJed Brown } 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMPIAdj(comm,4,16,ia,ja,NULL,&A)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 61c4762a1bSJed Brown /* 62c4762a1bSJed Brown Partition the graph of the matrix 63c4762a1bSJed Brown */ 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningCreate(comm,&part)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetAdjacency(part,A)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetType(part,MATPARTITIONINGHIERARCH)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningHierarchicalSetNcoarseparts(part,2)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningHierarchicalSetNfineparts(part,2)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetFromOptions(part)); 70c4762a1bSJed Brown /* get new processor owner number of each vertex */ 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningApply(part,&is)); 72c4762a1bSJed Brown /* coarse parts */ 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningHierarchicalGetCoarseparts(part,&coarseparts)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(coarseparts,PETSC_VIEWER_STDOUT_WORLD)); 75c4762a1bSJed Brown /* fine parts */ 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningHierarchicalGetFineparts(part,&fineparts)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(fineparts,PETSC_VIEWER_STDOUT_WORLD)); 78c4762a1bSJed Brown /* partitioning */ 795f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(is,PETSC_VIEWER_STDOUT_WORLD)); 80c4762a1bSJed Brown /* compute coming rows */ 815f80ce2aSJacob Faibussowitsch CHKERRQ(ISBuildTwoSided(is,NULL,&isrows)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(isrows,PETSC_VIEWER_STDOUT_WORLD)); 83c4762a1bSJed Brown /*create a sub-communicator */ 845f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_split(comm, membershipKey,rank,&scomm)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrows,&isrows_localsize)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(isrows_localsize,&indices_sc)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrows,&indices)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(indices_sc,indices,isrows_localsize)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrows,&indices)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&coarseparts)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&fineparts)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningDestroy(&part)); 95c4762a1bSJed Brown /*create a sub-IS on the sub communicator */ 965f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(scomm,isrows_localsize,indices_sc,PETSC_OWN_POINTER,&isrows_sc)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATMPIAIJ,MAT_INITIAL_MATRIX,&B)); 98c4762a1bSJed Brown #if 1 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 100c4762a1bSJed Brown #endif 101c4762a1bSJed Brown /*increase overlap */ 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatIncreaseOverlapSplit(B,1,&isrows_sc,1)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(isrows_sc,NULL)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows_sc)); 105c4762a1bSJed Brown /* 106c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 107c4762a1bSJed Brown are no longer needed. 108c4762a1bSJed Brown */ 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 111*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 112*b122ec5aSJacob Faibussowitsch return 0; 113c4762a1bSJed Brown } 114