1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Partition tiny grid.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices. Note that this file 6c4762a1bSJed Brown automatically includes: 7c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 8c4762a1bSJed Brown petscmat.h - matrices 9c4762a1bSJed Brown petscis.h - index sets 10c4762a1bSJed Brown petscviewer.h - viewers 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown #include <petscmat.h> 13c4762a1bSJed Brown 14c4762a1bSJed Brown int main(int argc,char **args) 15c4762a1bSJed Brown { 16*6a09307cSBarry Smith Mat A,At; 17c4762a1bSJed Brown PetscMPIInt rank,size; 18*6a09307cSBarry Smith PetscInt *ia,*ja,row; 19c4762a1bSJed Brown MatPartitioning part; 20c4762a1bSJed Brown IS is,isn; 21*6a09307cSBarry Smith PetscBool equal; 22c4762a1bSJed Brown 239566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,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) { 31c4762a1bSJed 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; 32c4762a1bSJed Brown ja[8] = 2; ja[9] = 7; 33c4762a1bSJed Brown ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 34c4762a1bSJed Brown } else if (rank == 1) { 35c4762a1bSJed 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; 36c4762a1bSJed Brown ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11; 37c4762a1bSJed Brown ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 38c4762a1bSJed Brown } else if (rank == 2) { 39c4762a1bSJed 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; 40c4762a1bSJed Brown ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15; 41c4762a1bSJed Brown ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 42c4762a1bSJed Brown } else { 43c4762a1bSJed 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; 44c4762a1bSJed Brown ja[8] = 11; ja[9] = 14; 45c4762a1bSJed Brown ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 489566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_WORLD,4,16,ia,ja,NULL,&A)); 499566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 50c4762a1bSJed Brown 51*6a09307cSBarry Smith /* Create the same matrix but using MatSetValues() */ 52*6a09307cSBarry Smith PetscCall(MatCreate(PETSC_COMM_WORLD,&At)); 53*6a09307cSBarry Smith PetscCall(MatSetSizes(At,4,4,16,16)); 54*6a09307cSBarry Smith PetscCall(MatSetType(At,MATMPIADJ)); 55*6a09307cSBarry Smith for (PetscInt i=0; i<4; i++) { 56*6a09307cSBarry Smith row = i + 4*rank; 57*6a09307cSBarry Smith PetscCall(MatSetValues(At,1,&row,ia[i+1]-ia[i],ja + ia[i],NULL,INSERT_VALUES)); 58*6a09307cSBarry Smith } 59*6a09307cSBarry Smith PetscCall(MatAssemblyBegin(At,MAT_FINAL_ASSEMBLY)); 60*6a09307cSBarry Smith PetscCall(MatAssemblyEnd(At,MAT_FINAL_ASSEMBLY)); 61*6a09307cSBarry Smith PetscCall(MatEqual(A,At,&equal)); 62*6a09307cSBarry Smith PetscCheck(equal,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Matrices are not equal that should be equal"); 63*6a09307cSBarry Smith PetscCall(MatDestroy(&At)); 64*6a09307cSBarry Smith 65c4762a1bSJed Brown /* 66c4762a1bSJed Brown Partition the graph of the matrix 67c4762a1bSJed Brown */ 689566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PETSC_COMM_WORLD,&part)); 699566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(part,A)); 709566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(part)); 71c4762a1bSJed Brown /* get new processor owner number of each vertex */ 729566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(part,&is)); 73c4762a1bSJed Brown /* get new global number of each old global number */ 749566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(is,&isn)); 759566063dSJacob Faibussowitsch PetscCall(ISView(isn,PETSC_VIEWER_STDOUT_WORLD)); 769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 77c4762a1bSJed Brown 789566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isn)); 799566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&part)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* 82c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 83c4762a1bSJed Brown are no longer needed. 84c4762a1bSJed Brown */ 859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 86c4762a1bSJed Brown 879566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 88b122ec5aSJacob Faibussowitsch return 0; 89c4762a1bSJed Brown } 90*6a09307cSBarry Smith /* 91*6a09307cSBarry Smith test: 92*6a09307cSBarry Smith requires: parmetis 93*6a09307cSBarry Smith args: -mat_view 94*6a09307cSBarry Smith nsize: 4 95*6a09307cSBarry Smith */ 96