1c4762a1bSJed Brown static char help[] = "Partition tiny grid.\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 { 156a09307cSBarry Smith Mat A, At; 16c4762a1bSJed Brown PetscMPIInt rank, size; 176a09307cSBarry Smith PetscInt *ia, *ja, row; 18c4762a1bSJed Brown MatPartitioning part; 19c4762a1bSJed Brown IS is, isn; 206a09307cSBarry Smith PetscBool equal; 21c4762a1bSJed Brown 22327415f7SBarry Smith PetscFunctionBeginUser; 23*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 25be096a46SBarry Smith PetscCheck(size == 4, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must run with 4 processors"); 269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 27c4762a1bSJed Brown 289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5, &ia)); 299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(16, &ja)); 30dd400576SPatrick Sanan if (rank == 0) { 319371c9d4SSatish Balay ja[0] = 1; 329371c9d4SSatish Balay ja[1] = 4; 339371c9d4SSatish Balay ja[2] = 0; 349371c9d4SSatish Balay ja[3] = 2; 359371c9d4SSatish Balay ja[4] = 5; 369371c9d4SSatish Balay ja[5] = 1; 379371c9d4SSatish Balay ja[6] = 3; 389371c9d4SSatish Balay ja[7] = 6; 399371c9d4SSatish Balay ja[8] = 2; 409371c9d4SSatish Balay ja[9] = 7; 419371c9d4SSatish Balay ia[0] = 0; 429371c9d4SSatish Balay ia[1] = 2; 439371c9d4SSatish Balay ia[2] = 5; 449371c9d4SSatish Balay ia[3] = 8; 459371c9d4SSatish Balay ia[4] = 10; 46c4762a1bSJed Brown } else if (rank == 1) { 479371c9d4SSatish Balay ja[0] = 0; 489371c9d4SSatish Balay ja[1] = 5; 499371c9d4SSatish Balay ja[2] = 8; 509371c9d4SSatish Balay ja[3] = 1; 519371c9d4SSatish Balay ja[4] = 4; 529371c9d4SSatish Balay ja[5] = 6; 539371c9d4SSatish Balay ja[6] = 9; 549371c9d4SSatish Balay ja[7] = 2; 559371c9d4SSatish Balay ja[8] = 5; 569371c9d4SSatish Balay ja[9] = 7; 579371c9d4SSatish Balay ja[10] = 10; 589371c9d4SSatish Balay ja[11] = 3; 599371c9d4SSatish Balay ja[12] = 6; 609371c9d4SSatish Balay ja[13] = 11; 619371c9d4SSatish Balay ia[0] = 0; 629371c9d4SSatish Balay ia[1] = 3; 639371c9d4SSatish Balay ia[2] = 7; 649371c9d4SSatish Balay ia[3] = 11; 659371c9d4SSatish Balay ia[4] = 14; 66c4762a1bSJed Brown } else if (rank == 2) { 679371c9d4SSatish Balay ja[0] = 4; 689371c9d4SSatish Balay ja[1] = 9; 699371c9d4SSatish Balay ja[2] = 12; 709371c9d4SSatish Balay ja[3] = 5; 719371c9d4SSatish Balay ja[4] = 8; 729371c9d4SSatish Balay ja[5] = 10; 739371c9d4SSatish Balay ja[6] = 13; 749371c9d4SSatish Balay ja[7] = 6; 759371c9d4SSatish Balay ja[8] = 9; 769371c9d4SSatish Balay ja[9] = 11; 779371c9d4SSatish Balay ja[10] = 14; 789371c9d4SSatish Balay ja[11] = 7; 799371c9d4SSatish Balay ja[12] = 10; 809371c9d4SSatish Balay ja[13] = 15; 819371c9d4SSatish Balay ia[0] = 0; 829371c9d4SSatish Balay ia[1] = 3; 839371c9d4SSatish Balay ia[2] = 7; 849371c9d4SSatish Balay ia[3] = 11; 859371c9d4SSatish Balay ia[4] = 14; 86c4762a1bSJed Brown } else { 879371c9d4SSatish Balay ja[0] = 8; 889371c9d4SSatish Balay ja[1] = 13; 899371c9d4SSatish Balay ja[2] = 9; 909371c9d4SSatish Balay ja[3] = 12; 919371c9d4SSatish Balay ja[4] = 14; 929371c9d4SSatish Balay ja[5] = 10; 939371c9d4SSatish Balay ja[6] = 13; 949371c9d4SSatish Balay ja[7] = 15; 959371c9d4SSatish Balay ja[8] = 11; 969371c9d4SSatish Balay ja[9] = 14; 979371c9d4SSatish Balay ia[0] = 0; 989371c9d4SSatish Balay ia[1] = 2; 999371c9d4SSatish Balay ia[2] = 5; 1009371c9d4SSatish Balay ia[3] = 8; 1019371c9d4SSatish Balay ia[4] = 10; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 1049566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_WORLD, 4, 16, ia, ja, NULL, &A)); 1059566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 106c4762a1bSJed Brown 1076a09307cSBarry Smith /* Create the same matrix but using MatSetValues() */ 1086a09307cSBarry Smith PetscCall(MatCreate(PETSC_COMM_WORLD, &At)); 1096a09307cSBarry Smith PetscCall(MatSetSizes(At, 4, 4, 16, 16)); 1106a09307cSBarry Smith PetscCall(MatSetType(At, MATMPIADJ)); 1116a09307cSBarry Smith for (PetscInt i = 0; i < 4; i++) { 1126a09307cSBarry Smith row = i + 4 * rank; 1136a09307cSBarry Smith PetscCall(MatSetValues(At, 1, &row, ia[i + 1] - ia[i], ja + ia[i], NULL, INSERT_VALUES)); 1146a09307cSBarry Smith } 1156a09307cSBarry Smith PetscCall(MatAssemblyBegin(At, MAT_FINAL_ASSEMBLY)); 1166a09307cSBarry Smith PetscCall(MatAssemblyEnd(At, MAT_FINAL_ASSEMBLY)); 1176a09307cSBarry Smith PetscCall(MatEqual(A, At, &equal)); 1186a09307cSBarry Smith PetscCheck(equal, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Matrices are not equal that should be equal"); 1196a09307cSBarry Smith PetscCall(MatDestroy(&At)); 1206a09307cSBarry Smith 121c4762a1bSJed Brown /* 122c4762a1bSJed Brown Partition the graph of the matrix 123c4762a1bSJed Brown */ 1249566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PETSC_COMM_WORLD, &part)); 1259566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(part, A)); 1269566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(part)); 127c4762a1bSJed Brown /* get new processor owner number of each vertex */ 1289566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(part, &is)); 129c4762a1bSJed Brown /* get new global number of each old global number */ 1309566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(is, &isn)); 1319566063dSJacob Faibussowitsch PetscCall(ISView(isn, PETSC_VIEWER_STDOUT_WORLD)); 1329566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 133c4762a1bSJed Brown 1349566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isn)); 1359566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&part)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* 138c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 139c4762a1bSJed Brown are no longer needed. 140c4762a1bSJed Brown */ 1419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 142c4762a1bSJed Brown 1439566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 144b122ec5aSJacob Faibussowitsch return 0; 145c4762a1bSJed Brown } 1466a09307cSBarry Smith /* 1476a09307cSBarry Smith test: 1486a09307cSBarry Smith requires: parmetis 1496a09307cSBarry Smith args: -mat_view 1506a09307cSBarry Smith nsize: 4 1516a09307cSBarry Smith */ 152