Actual source code: ex82.c


  2: static char help[] = "Partition a tiny grid using hierarchical partitioning.\n\n";

  4: /*
  5:   Include "petscmat.h" so that we can use matrices.  Note that this file
  6:   automatically includes:
  7:      petscsys.h       - base PETSc routines   petscvec.h - vectors
  8:      petscmat.h - matrices
  9:      petscis.h     - index sets
 10:      petscviewer.h - viewers
 11: */
 12: #include <petscmat.h>

 14: int main(int argc,char **args)
 15: {
 16:   Mat             A;
 17:   PetscMPIInt     rank,size;
 18:   PetscInt        *ia,*ja;
 19:   MatPartitioning part;
 20:   IS              is,isn,isrows;
 21:   IS              coarseparts,fineparts;
 22:   MPI_Comm        comm;

 24:   PetscInitialize(&argc,&args,(char*)0,help);
 25:   comm = PETSC_COMM_WORLD;
 26:   MPI_Comm_size(comm,&size);
 28:   MPI_Comm_rank(comm,&rank);

 30:   PetscMalloc1(5,&ia);
 31:   PetscMalloc1(16,&ja);
 32:   if (rank == 0) {
 33:     ja[0] = 1; ja[1] = 4; ja[2] = 0; ja[3] = 2; ja[4] = 5; ja[5] = 1; ja[6] = 3; ja[7] = 6;
 34:     ja[8] = 2; ja[9] = 7;
 35:     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
 36:   } else if (rank == 1) {
 37:     ja[0] = 0; ja[1] = 5; ja[2] = 8; ja[3] = 1; ja[4] = 4; ja[5] = 6; ja[6] = 9; ja[7] = 2;
 38:     ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11;
 39:     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
 40:   } else if (rank == 2) {
 41:     ja[0] = 4; ja[1] = 9; ja[2] = 12; ja[3] = 5; ja[4] = 8; ja[5] = 10; ja[6] = 13; ja[7] = 6;
 42:     ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15;
 43:     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
 44:   } else {
 45:     ja[0] = 8; ja[1] = 13; ja[2] = 9; ja[3] = 12; ja[4] = 14; ja[5] = 10; ja[6] = 13; ja[7] = 15;
 46:     ja[8] = 11; ja[9] = 14;
 47:     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
 48:   }
 49:   MatCreateMPIAdj(comm,4,16,ia,ja,NULL,&A);
 50:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);
 51:   /*
 52:    Partition the graph of the matrix
 53:   */
 54:   MatPartitioningCreate(comm,&part);
 55:   MatPartitioningSetAdjacency(part,A);
 56:   MatPartitioningSetType(part,MATPARTITIONINGHIERARCH);
 57:   MatPartitioningHierarchicalSetNcoarseparts(part,2);
 58:   MatPartitioningHierarchicalSetNfineparts(part,2);
 59:   MatPartitioningSetFromOptions(part);
 60:   /* get new processor owner number of each vertex */
 61:   MatPartitioningApply(part,&is);
 62:   /* coarse parts */
 63:   MatPartitioningHierarchicalGetCoarseparts(part,&coarseparts);
 64:   ISView(coarseparts,PETSC_VIEWER_STDOUT_WORLD);
 65:   /* fine parts */
 66:   MatPartitioningHierarchicalGetFineparts(part,&fineparts);
 67:   ISView(fineparts,PETSC_VIEWER_STDOUT_WORLD);
 68:   /* partitioning */
 69:   ISView(is,PETSC_VIEWER_STDOUT_WORLD);
 70:   /* get new global number of each old global number */
 71:   ISPartitioningToNumbering(is,&isn);
 72:   ISView(isn,PETSC_VIEWER_STDOUT_WORLD);
 73:   ISBuildTwoSided(is,NULL,&isrows);
 74:   ISView(isrows,PETSC_VIEWER_STDOUT_WORLD);
 75:   ISDestroy(&is);
 76:   ISDestroy(&coarseparts);
 77:   ISDestroy(&fineparts);
 78:   ISDestroy(&isrows);
 79:   ISDestroy(&isn);
 80:   MatPartitioningDestroy(&part);
 81:   /*
 82:     Free work space.  All PETSc objects should be destroyed when they
 83:     are no longer needed.
 84:   */
 85:   MatDestroy(&A);
 86:   PetscFinalize();
 87:   return 0;
 88: }

 90: /*TEST

 92:    test:
 93:       nsize: 4
 94:       requires: parmetis
 95:       TODO: tests cannot use parmetis because it produces different results on different machines

 97: TEST*/