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