xref: /petsc/src/mat/tests/ex193.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 /*
2  * ex193.c
3  *
4  *  Created on: Jul 29, 2015
5  *      Author: Fande Kong fdkong.jd@gmail.com
6  */
7 /*
8  * An example demonstrates how to use hierarchical partitioning approach
9  */
10 
11 #include <petscmat.h>
12 
13 static char help[] = "Illustrates use of hierarchical partitioning.\n";
14 
15 int main(int argc,char **args)
16 {
17   Mat             A;                      /* matrix */
18   PetscInt        m,n;                    /* mesh dimensions in x- and y- directions */
19   PetscInt        i,j,Ii,J,Istart,Iend;
20   PetscErrorCode  ierr;
21   PetscMPIInt     size;
22   PetscScalar     v;
23   MatPartitioning part;
24   IS              coarseparts,fineparts;
25   IS              is,isn,isrows;
26   MPI_Comm        comm;
27 
28   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
29   comm = PETSC_COMM_WORLD;
30   CHKERRMPI(MPI_Comm_size(comm,&size));
31   ierr = PetscOptionsBegin(comm,NULL,"ex193","hierarchical partitioning");CHKERRQ(ierr);
32   m = 15;
33   CHKERRQ(PetscOptionsInt("-M","Number of mesh points in the x-direction","partitioning",m,&m,NULL));
34   n = 15;
35   CHKERRQ(PetscOptionsInt("-N","Number of mesh points in the y-direction","partitioning",n,&n,NULL));
36   ierr = PetscOptionsEnd();CHKERRQ(ierr);
37 
38   /*
39      Assemble the matrix for the five point stencil (finite difference), YET AGAIN
40   */
41   CHKERRQ(MatCreate(comm,&A));
42   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
43   CHKERRQ(MatSetFromOptions(A));
44   CHKERRQ(MatSetUp(A));
45   CHKERRQ(MatGetOwnershipRange(A,&Istart,&Iend));
46   for (Ii=Istart; Ii<Iend; Ii++) {
47     v = -1.0; i = Ii/n; j = Ii - i*n;
48     if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
49     if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
50     if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
51     if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
52     v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
53   }
54   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
55   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
56   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
57   /*
58    Partition the graph of the matrix
59   */
60   CHKERRQ(MatPartitioningCreate(comm,&part));
61   CHKERRQ(MatPartitioningSetAdjacency(part,A));
62   CHKERRQ(MatPartitioningSetType(part,MATPARTITIONINGHIERARCH));
63   CHKERRQ(MatPartitioningHierarchicalSetNcoarseparts(part,2));
64   CHKERRQ(MatPartitioningHierarchicalSetNfineparts(part,4));
65   CHKERRQ(MatPartitioningSetFromOptions(part));
66   /* get new processor owner number of each vertex */
67   CHKERRQ(MatPartitioningApply(part,&is));
68   /* coarse parts */
69   CHKERRQ(MatPartitioningHierarchicalGetCoarseparts(part,&coarseparts));
70   CHKERRQ(ISView(coarseparts,PETSC_VIEWER_STDOUT_WORLD));
71   /* fine parts */
72   CHKERRQ(MatPartitioningHierarchicalGetFineparts(part,&fineparts));
73   CHKERRQ(ISView(fineparts,PETSC_VIEWER_STDOUT_WORLD));
74   /* partitioning */
75   CHKERRQ(ISView(is,PETSC_VIEWER_STDOUT_WORLD));
76   /* get new global number of each old global number */
77   CHKERRQ(ISPartitioningToNumbering(is,&isn));
78   CHKERRQ(ISView(isn,PETSC_VIEWER_STDOUT_WORLD));
79   CHKERRQ(ISBuildTwoSided(is,NULL,&isrows));
80   CHKERRQ(ISView(isrows,PETSC_VIEWER_STDOUT_WORLD));
81   CHKERRQ(ISDestroy(&is));
82   CHKERRQ(ISDestroy(&coarseparts));
83   CHKERRQ(ISDestroy(&fineparts));
84   CHKERRQ(ISDestroy(&isrows));
85   CHKERRQ(ISDestroy(&isn));
86   CHKERRQ(MatPartitioningDestroy(&part));
87   CHKERRQ(MatDestroy(&A));
88   ierr = PetscFinalize();
89   return ierr;
90 }
91 
92 /*TEST
93 
94    test:
95       nsize: 4
96       args: -mat_partitioning_hierarchical_Nfineparts 2
97       requires: parmetis
98       TODO: cannot run because parmetis does reproduce across all machines, probably due to nonportable random number generator
99 
100 TEST*/
101