1 2 static char help[] = "Partition tiny grid using hierarchical partitioning and increase overlap using MatIncreaseOverlapSplit.\n\n"; 3 4 /*T 5 Concepts: partitioning 6 Processors: 4 7 T*/ 8 9 /* 10 Include "petscmat.h" so that we can use matrices. Note that this file 11 automatically includes: 12 petscsys.h - base PETSc routines petscvec.h - vectors 13 petscmat.h - matrices 14 petscis.h - index sets 15 petscviewer.h - viewers 16 */ 17 #include <petscmat.h> 18 19 int main(int argc,char **args) 20 { 21 Mat A,B; 22 PetscErrorCode ierr; 23 PetscMPIInt rank,size,membershipKey; 24 PetscInt *ia,*ja,*indices_sc,isrows_localsize; 25 const PetscInt *indices; 26 MatPartitioning part; 27 IS is,isrows,isrows_sc; 28 IS coarseparts,fineparts; 29 MPI_Comm comm,scomm; 30 31 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 32 comm = PETSC_COMM_WORLD; 33 CHKERRMPI(MPI_Comm_size(comm,&size)); 34 PetscCheckFalse(size != 4,comm,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 4 processors "); 35 CHKERRMPI(MPI_Comm_rank(comm,&rank)); 36 /*set a small matrix */ 37 CHKERRQ(PetscMalloc1(5,&ia)); 38 CHKERRQ(PetscMalloc1(16,&ja)); 39 if (rank == 0) { 40 ja[0] = 1; ja[1] = 4; ja[2] = 0; ja[3] = 2; ja[4] = 5; ja[5] = 1; ja[6] = 3; ja[7] = 6; 41 ja[8] = 2; ja[9] = 7; 42 ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 43 membershipKey = 0; 44 } else if (rank == 1) { 45 ja[0] = 0; ja[1] = 5; ja[2] = 8; ja[3] = 1; ja[4] = 4; ja[5] = 6; ja[6] = 9; ja[7] = 2; 46 ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11; 47 ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 48 membershipKey = 0; 49 } else if (rank == 2) { 50 ja[0] = 4; ja[1] = 9; ja[2] = 12; ja[3] = 5; ja[4] = 8; ja[5] = 10; ja[6] = 13; ja[7] = 6; 51 ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15; 52 ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 53 membershipKey = 1; 54 } else { 55 ja[0] = 8; ja[1] = 13; ja[2] = 9; ja[3] = 12; ja[4] = 14; ja[5] = 10; ja[6] = 13; ja[7] = 15; 56 ja[8] = 11; ja[9] = 14; 57 ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 58 membershipKey = 1; 59 } 60 CHKERRQ(MatCreateMPIAdj(comm,4,16,ia,ja,NULL,&A)); 61 CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 62 /* 63 Partition the graph of the matrix 64 */ 65 CHKERRQ(MatPartitioningCreate(comm,&part)); 66 CHKERRQ(MatPartitioningSetAdjacency(part,A)); 67 CHKERRQ(MatPartitioningSetType(part,MATPARTITIONINGHIERARCH)); 68 CHKERRQ(MatPartitioningHierarchicalSetNcoarseparts(part,2)); 69 CHKERRQ(MatPartitioningHierarchicalSetNfineparts(part,2)); 70 CHKERRQ(MatPartitioningSetFromOptions(part)); 71 /* get new processor owner number of each vertex */ 72 CHKERRQ(MatPartitioningApply(part,&is)); 73 /* coarse parts */ 74 CHKERRQ(MatPartitioningHierarchicalGetCoarseparts(part,&coarseparts)); 75 CHKERRQ(ISView(coarseparts,PETSC_VIEWER_STDOUT_WORLD)); 76 /* fine parts */ 77 CHKERRQ(MatPartitioningHierarchicalGetFineparts(part,&fineparts)); 78 CHKERRQ(ISView(fineparts,PETSC_VIEWER_STDOUT_WORLD)); 79 /* partitioning */ 80 CHKERRQ(ISView(is,PETSC_VIEWER_STDOUT_WORLD)); 81 /* compute coming rows */ 82 CHKERRQ(ISBuildTwoSided(is,NULL,&isrows)); 83 CHKERRQ(ISView(isrows,PETSC_VIEWER_STDOUT_WORLD)); 84 /*create a sub-communicator */ 85 CHKERRMPI(MPI_Comm_split(comm, membershipKey,rank,&scomm)); 86 CHKERRQ(ISGetLocalSize(isrows,&isrows_localsize)); 87 CHKERRQ(PetscMalloc1(isrows_localsize,&indices_sc)); 88 CHKERRQ(ISGetIndices(isrows,&indices)); 89 CHKERRQ(PetscArraycpy(indices_sc,indices,isrows_localsize)); 90 CHKERRQ(ISRestoreIndices(isrows,&indices)); 91 CHKERRQ(ISDestroy(&is)); 92 CHKERRQ(ISDestroy(&coarseparts)); 93 CHKERRQ(ISDestroy(&fineparts)); 94 CHKERRQ(ISDestroy(&isrows)); 95 CHKERRQ(MatPartitioningDestroy(&part)); 96 /*create a sub-IS on the sub communicator */ 97 CHKERRQ(ISCreateGeneral(scomm,isrows_localsize,indices_sc,PETSC_OWN_POINTER,&isrows_sc)); 98 CHKERRQ(MatConvert(A,MATMPIAIJ,MAT_INITIAL_MATRIX,&B)); 99 #if 1 100 CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 101 #endif 102 /*increase overlap */ 103 CHKERRQ(MatIncreaseOverlapSplit(B,1,&isrows_sc,1)); 104 CHKERRQ(ISView(isrows_sc,NULL)); 105 CHKERRQ(ISDestroy(&isrows_sc)); 106 /* 107 Free work space. All PETSc objects should be destroyed when they 108 are no longer needed. 109 */ 110 CHKERRQ(MatDestroy(&A)); 111 CHKERRQ(MatDestroy(&B)); 112 ierr = PetscFinalize(); 113 return ierr; 114 } 115