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