1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Partition tiny grid.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*T 5c4762a1bSJed Brown Concepts: partitioning 6c4762a1bSJed Brown Processors: 4 7c4762a1bSJed Brown T*/ 8c4762a1bSJed Brown 9c4762a1bSJed Brown /* 10c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices. Note that this file 11c4762a1bSJed Brown automatically includes: 12c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 13c4762a1bSJed Brown petscmat.h - matrices 14c4762a1bSJed Brown petscis.h - index sets 15c4762a1bSJed Brown petscviewer.h - viewers 16c4762a1bSJed Brown */ 17c4762a1bSJed Brown #include <petscmat.h> 18c4762a1bSJed Brown 19c4762a1bSJed Brown int main(int argc,char **args) 20c4762a1bSJed Brown { 21c4762a1bSJed Brown Mat A; 22c4762a1bSJed Brown PetscMPIInt rank,size; 23c4762a1bSJed Brown PetscInt *ia,*ja; 24c4762a1bSJed Brown MatPartitioning part; 25c4762a1bSJed Brown IS is,isn; 26c4762a1bSJed Brown 27*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 285f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 292c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 4,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 4 processors"); 305f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 31c4762a1bSJed Brown 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(5,&ia)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(16,&ja)); 34dd400576SPatrick Sanan if (rank == 0) { 35c4762a1bSJed 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; 36c4762a1bSJed Brown ja[8] = 2; ja[9] = 7; 37c4762a1bSJed Brown ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 38c4762a1bSJed Brown } else if (rank == 1) { 39c4762a1bSJed 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; 40c4762a1bSJed Brown ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11; 41c4762a1bSJed Brown ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 42c4762a1bSJed Brown } else if (rank == 2) { 43c4762a1bSJed 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; 44c4762a1bSJed Brown ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15; 45c4762a1bSJed Brown ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 46c4762a1bSJed Brown } else { 47c4762a1bSJed 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; 48c4762a1bSJed Brown ja[8] = 11; ja[9] = 14; 49c4762a1bSJed Brown ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 50c4762a1bSJed Brown } 51c4762a1bSJed Brown 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMPIAdj(PETSC_COMM_WORLD,4,16,ia,ja,NULL,&A)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* 56c4762a1bSJed Brown Partition the graph of the matrix 57c4762a1bSJed Brown */ 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningCreate(PETSC_COMM_WORLD,&part)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetAdjacency(part,A)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetFromOptions(part)); 61c4762a1bSJed Brown /* get new processor owner number of each vertex */ 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningApply(part,&is)); 63c4762a1bSJed Brown /* get new global number of each old global number */ 645f80ce2aSJacob Faibussowitsch CHKERRQ(ISPartitioningToNumbering(is,&isn)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(isn,PETSC_VIEWER_STDOUT_WORLD)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 67c4762a1bSJed Brown 685f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isn)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningDestroy(&part)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* 72c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 73c4762a1bSJed Brown are no longer needed. 74c4762a1bSJed Brown */ 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 76c4762a1bSJed Brown 77*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 78*b122ec5aSJacob Faibussowitsch return 0; 79c4762a1bSJed Brown } 80