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