1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Partition a tiny grid using hierarchical partitioning.\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 11*c4762a1bSJed Brown /* 12*c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices. Note that this file 13*c4762a1bSJed Brown automatically includes: 14*c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 15*c4762a1bSJed Brown petscmat.h - matrices 16*c4762a1bSJed Brown petscis.h - index sets 17*c4762a1bSJed Brown petscviewer.h - viewers 18*c4762a1bSJed Brown */ 19*c4762a1bSJed Brown #include <petscmat.h> 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown int main(int argc,char **args) 22*c4762a1bSJed Brown { 23*c4762a1bSJed Brown Mat A; 24*c4762a1bSJed Brown PetscErrorCode ierr; 25*c4762a1bSJed Brown PetscMPIInt rank,size; 26*c4762a1bSJed Brown PetscInt *ia,*ja; 27*c4762a1bSJed Brown MatPartitioning part; 28*c4762a1bSJed Brown IS is,isn,isrows; 29*c4762a1bSJed Brown IS coarseparts,fineparts; 30*c4762a1bSJed Brown MPI_Comm comm; 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 33*c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 34*c4762a1bSJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 35*c4762a1bSJed Brown if (size != 4) SETERRQ(comm,1,"Must run with 4 processors"); 36*c4762a1bSJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown ierr = PetscMalloc1(5,&ia);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = PetscMalloc1(16,&ja);CHKERRQ(ierr); 40*c4762a1bSJed Brown if (!rank) { 41*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; 42*c4762a1bSJed Brown ja[8] = 2; ja[9] = 7; 43*c4762a1bSJed Brown ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 44*c4762a1bSJed Brown } else if (rank == 1) { 45*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; 46*c4762a1bSJed Brown ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11; 47*c4762a1bSJed Brown ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 48*c4762a1bSJed Brown } else if (rank == 2) { 49*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; 50*c4762a1bSJed Brown ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15; 51*c4762a1bSJed Brown ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 52*c4762a1bSJed Brown } else { 53*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; 54*c4762a1bSJed Brown ja[8] = 11; ja[9] = 14; 55*c4762a1bSJed Brown ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 56*c4762a1bSJed Brown } 57*c4762a1bSJed Brown ierr = MatCreateMPIAdj(comm,4,16,ia,ja,NULL,&A);CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 59*c4762a1bSJed Brown /* 60*c4762a1bSJed Brown Partition the graph of the matrix 61*c4762a1bSJed Brown */ 62*c4762a1bSJed Brown ierr = MatPartitioningCreate(comm,&part);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = MatPartitioningSetAdjacency(part,A);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = MatPartitioningSetType(part,MATPARTITIONINGHIERARCH);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = MatPartitioningHierarchicalSetNcoarseparts(part,2);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = MatPartitioningHierarchicalSetNfineparts(part,2);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr); 68*c4762a1bSJed Brown /* get new processor owner number of each vertex */ 69*c4762a1bSJed Brown ierr = MatPartitioningApply(part,&is);CHKERRQ(ierr); 70*c4762a1bSJed Brown /* coarse parts */ 71*c4762a1bSJed Brown ierr = MatPartitioningHierarchicalGetCoarseparts(part,&coarseparts);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = ISView(coarseparts,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 73*c4762a1bSJed Brown /* fine parts */ 74*c4762a1bSJed Brown ierr = MatPartitioningHierarchicalGetFineparts(part,&fineparts);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = ISView(fineparts,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 76*c4762a1bSJed Brown /* partitioning */ 77*c4762a1bSJed Brown ierr = ISView(is,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 78*c4762a1bSJed Brown /* get new global number of each old global number */ 79*c4762a1bSJed Brown ierr = ISPartitioningToNumbering(is,&isn);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = ISView(isn,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = ISBuildTwoSided(is,NULL,&isrows);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = ISView(isrows,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = ISDestroy(&coarseparts);CHKERRQ(ierr); 85*c4762a1bSJed Brown ierr = ISDestroy(&fineparts);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = ISDestroy(&isrows);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = ISDestroy(&isn);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr); 89*c4762a1bSJed Brown /* 90*c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 91*c4762a1bSJed Brown are no longer needed. 92*c4762a1bSJed Brown */ 93*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 94*c4762a1bSJed Brown ierr = PetscFinalize(); 95*c4762a1bSJed Brown return ierr; 96*c4762a1bSJed Brown } 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown 99*c4762a1bSJed Brown /*TEST 100*c4762a1bSJed Brown 101*c4762a1bSJed Brown test: 102*c4762a1bSJed Brown nsize: 4 103*c4762a1bSJed Brown requires: parmetis 104*c4762a1bSJed Brown TODO: tests cannot use parmetis because it produces different results on different machines 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown TEST*/ 107