1c4762a1bSJed Brown static char help[] = "Partition tiny grid using hierarchical partitioning and increase overlap using MatIncreaseOverlapSplit.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices. Note that this file 5c4762a1bSJed Brown automatically includes: 6c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 7c4762a1bSJed Brown petscmat.h - matrices 8c4762a1bSJed Brown petscis.h - index sets 9c4762a1bSJed Brown petscviewer.h - viewers 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown #include <petscmat.h> 12c4762a1bSJed Brown 13d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 14d71ae5a4SJacob Faibussowitsch { 15c4762a1bSJed Brown Mat A, B; 16c4762a1bSJed Brown PetscMPIInt rank, size, membershipKey; 17c4762a1bSJed Brown PetscInt *ia, *ja, *indices_sc, isrows_localsize; 18c4762a1bSJed Brown const PetscInt *indices; 19c4762a1bSJed Brown MatPartitioning part; 20c4762a1bSJed Brown IS is, isrows, isrows_sc; 21c4762a1bSJed Brown IS coarseparts, fineparts; 22c4762a1bSJed Brown MPI_Comm comm, scomm; 23c4762a1bSJed Brown 24327415f7SBarry Smith PetscFunctionBeginUser; 25*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, 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 /*set a small matrix */ 319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5, &ia)); 329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(16, &ja)); 33dd400576SPatrick Sanan if (rank == 0) { 349371c9d4SSatish Balay ja[0] = 1; 359371c9d4SSatish Balay ja[1] = 4; 369371c9d4SSatish Balay ja[2] = 0; 379371c9d4SSatish Balay ja[3] = 2; 389371c9d4SSatish Balay ja[4] = 5; 399371c9d4SSatish Balay ja[5] = 1; 409371c9d4SSatish Balay ja[6] = 3; 419371c9d4SSatish Balay ja[7] = 6; 429371c9d4SSatish Balay ja[8] = 2; 439371c9d4SSatish Balay ja[9] = 7; 449371c9d4SSatish Balay ia[0] = 0; 459371c9d4SSatish Balay ia[1] = 2; 469371c9d4SSatish Balay ia[2] = 5; 479371c9d4SSatish Balay ia[3] = 8; 489371c9d4SSatish Balay ia[4] = 10; 49c4762a1bSJed Brown membershipKey = 0; 50c4762a1bSJed Brown } else if (rank == 1) { 519371c9d4SSatish Balay ja[0] = 0; 529371c9d4SSatish Balay ja[1] = 5; 539371c9d4SSatish Balay ja[2] = 8; 549371c9d4SSatish Balay ja[3] = 1; 559371c9d4SSatish Balay ja[4] = 4; 569371c9d4SSatish Balay ja[5] = 6; 579371c9d4SSatish Balay ja[6] = 9; 589371c9d4SSatish Balay ja[7] = 2; 599371c9d4SSatish Balay ja[8] = 5; 609371c9d4SSatish Balay ja[9] = 7; 619371c9d4SSatish Balay ja[10] = 10; 629371c9d4SSatish Balay ja[11] = 3; 639371c9d4SSatish Balay ja[12] = 6; 649371c9d4SSatish Balay ja[13] = 11; 659371c9d4SSatish Balay ia[0] = 0; 669371c9d4SSatish Balay ia[1] = 3; 679371c9d4SSatish Balay ia[2] = 7; 689371c9d4SSatish Balay ia[3] = 11; 699371c9d4SSatish Balay ia[4] = 14; 70c4762a1bSJed Brown membershipKey = 0; 71c4762a1bSJed Brown } else if (rank == 2) { 729371c9d4SSatish Balay ja[0] = 4; 739371c9d4SSatish Balay ja[1] = 9; 749371c9d4SSatish Balay ja[2] = 12; 759371c9d4SSatish Balay ja[3] = 5; 769371c9d4SSatish Balay ja[4] = 8; 779371c9d4SSatish Balay ja[5] = 10; 789371c9d4SSatish Balay ja[6] = 13; 799371c9d4SSatish Balay ja[7] = 6; 809371c9d4SSatish Balay ja[8] = 9; 819371c9d4SSatish Balay ja[9] = 11; 829371c9d4SSatish Balay ja[10] = 14; 839371c9d4SSatish Balay ja[11] = 7; 849371c9d4SSatish Balay ja[12] = 10; 859371c9d4SSatish Balay ja[13] = 15; 869371c9d4SSatish Balay ia[0] = 0; 879371c9d4SSatish Balay ia[1] = 3; 889371c9d4SSatish Balay ia[2] = 7; 899371c9d4SSatish Balay ia[3] = 11; 909371c9d4SSatish Balay ia[4] = 14; 91c4762a1bSJed Brown membershipKey = 1; 92c4762a1bSJed Brown } else { 939371c9d4SSatish Balay ja[0] = 8; 949371c9d4SSatish Balay ja[1] = 13; 959371c9d4SSatish Balay ja[2] = 9; 969371c9d4SSatish Balay ja[3] = 12; 979371c9d4SSatish Balay ja[4] = 14; 989371c9d4SSatish Balay ja[5] = 10; 999371c9d4SSatish Balay ja[6] = 13; 1009371c9d4SSatish Balay ja[7] = 15; 1019371c9d4SSatish Balay ja[8] = 11; 1029371c9d4SSatish Balay ja[9] = 14; 1039371c9d4SSatish Balay ia[0] = 0; 1049371c9d4SSatish Balay ia[1] = 2; 1059371c9d4SSatish Balay ia[2] = 5; 1069371c9d4SSatish Balay ia[3] = 8; 1079371c9d4SSatish Balay ia[4] = 10; 108c4762a1bSJed Brown membershipKey = 1; 109c4762a1bSJed Brown } 1109566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(comm, 4, 16, ia, ja, NULL, &A)); 1119566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 112c4762a1bSJed Brown /* 113c4762a1bSJed Brown Partition the graph of the matrix 114c4762a1bSJed Brown */ 1159566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(comm, &part)); 1169566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(part, A)); 1179566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetType(part, MATPARTITIONINGHIERARCH)); 1189566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNcoarseparts(part, 2)); 1199566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNfineparts(part, 2)); 1209566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(part)); 121c4762a1bSJed Brown /* get new processor owner number of each vertex */ 1229566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(part, &is)); 123c4762a1bSJed Brown /* coarse parts */ 1249566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalGetCoarseparts(part, &coarseparts)); 1259566063dSJacob Faibussowitsch PetscCall(ISView(coarseparts, PETSC_VIEWER_STDOUT_WORLD)); 126c4762a1bSJed Brown /* fine parts */ 1279566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalGetFineparts(part, &fineparts)); 1289566063dSJacob Faibussowitsch PetscCall(ISView(fineparts, PETSC_VIEWER_STDOUT_WORLD)); 129c4762a1bSJed Brown /* partitioning */ 1309566063dSJacob Faibussowitsch PetscCall(ISView(is, PETSC_VIEWER_STDOUT_WORLD)); 131c4762a1bSJed Brown /* compute coming rows */ 1329566063dSJacob Faibussowitsch PetscCall(ISBuildTwoSided(is, NULL, &isrows)); 1339566063dSJacob Faibussowitsch PetscCall(ISView(isrows, PETSC_VIEWER_STDOUT_WORLD)); 134c4762a1bSJed Brown /*create a sub-communicator */ 1359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_split(comm, membershipKey, rank, &scomm)); 1369566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows, &isrows_localsize)); 1379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(isrows_localsize, &indices_sc)); 1389566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows, &indices)); 1399566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(indices_sc, indices, isrows_localsize)); 1409566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows, &indices)); 1419566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 1429566063dSJacob Faibussowitsch PetscCall(ISDestroy(&coarseparts)); 1439566063dSJacob Faibussowitsch PetscCall(ISDestroy(&fineparts)); 1449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 1459566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&part)); 146c4762a1bSJed Brown /*create a sub-IS on the sub communicator */ 1479566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(scomm, isrows_localsize, indices_sc, PETSC_OWN_POINTER, &isrows_sc)); 1489566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATMPIAIJ, MAT_INITIAL_MATRIX, &B)); 149c4762a1bSJed Brown #if 1 1509566063dSJacob Faibussowitsch PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 151c4762a1bSJed Brown #endif 152c4762a1bSJed Brown /*increase overlap */ 1539566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlapSplit(B, 1, &isrows_sc, 1)); 1549566063dSJacob Faibussowitsch PetscCall(ISView(isrows_sc, NULL)); 1559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows_sc)); 156c4762a1bSJed Brown /* 157c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 158c4762a1bSJed Brown are no longer needed. 159c4762a1bSJed Brown */ 1609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1629566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 163b122ec5aSJacob Faibussowitsch return 0; 164c4762a1bSJed Brown } 165