xref: /petsc/src/mat/tests/ex193.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
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   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
31   ierr = PetscOptionsBegin(comm,NULL,"ex193","hierarchical partitioning");CHKERRQ(ierr);
32   m = 15;
33   ierr = PetscOptionsInt("-M","Number of mesh points in the x-direction","partitioning",m,&m,NULL);CHKERRQ(ierr);
34   n = 15;
35   ierr = PetscOptionsInt("-N","Number of mesh points in the y-direction","partitioning",n,&n,NULL);CHKERRQ(ierr);
36   ierr = PetscOptionsEnd();CHKERRQ(ierr);
37 
38   /*
39      Assemble the matrix for the five point stencil (finite difference), YET AGAIN
40   */
41   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
42   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
43   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
44   ierr = MatSetUp(A);CHKERRQ(ierr);
45   ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);
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; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
49     if (i<m-1) {J = Ii + n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
50     if (j>0)   {J = Ii - 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
51     if (j<n-1) {J = Ii + 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
52     v = 4.0; ierr = MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
53   }
54   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
55   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
56   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
57   /*
58    Partition the graph of the matrix
59   */
60   ierr = MatPartitioningCreate(comm,&part);CHKERRQ(ierr);
61   ierr = MatPartitioningSetAdjacency(part,A);CHKERRQ(ierr);
62   ierr = MatPartitioningSetType(part,MATPARTITIONINGHIERARCH);CHKERRQ(ierr);
63   ierr = MatPartitioningHierarchicalSetNcoarseparts(part,2);CHKERRQ(ierr);
64   ierr = MatPartitioningHierarchicalSetNfineparts(part,4);CHKERRQ(ierr);
65   ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr);
66   /* get new processor owner number of each vertex */
67   ierr = MatPartitioningApply(part,&is);CHKERRQ(ierr);
68   /* coarse parts */
69   ierr = MatPartitioningHierarchicalGetCoarseparts(part,&coarseparts);CHKERRQ(ierr);
70   ierr = ISView(coarseparts,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
71   /* fine parts */
72   ierr = MatPartitioningHierarchicalGetFineparts(part,&fineparts);CHKERRQ(ierr);
73   ierr = ISView(fineparts,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
74   /* partitioning */
75   ierr = ISView(is,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
76   /* get new global number of each old global number */
77   ierr = ISPartitioningToNumbering(is,&isn);CHKERRQ(ierr);
78   ierr = ISView(isn,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
79   ierr = ISBuildTwoSided(is,NULL,&isrows);CHKERRQ(ierr);
80   ierr = ISView(isrows,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
81   ierr = ISDestroy(&is);CHKERRQ(ierr);
82   ierr = ISDestroy(&coarseparts);CHKERRQ(ierr);
83   ierr = ISDestroy(&fineparts);CHKERRQ(ierr);
84   ierr = ISDestroy(&isrows);CHKERRQ(ierr);
85   ierr = ISDestroy(&isn);CHKERRQ(ierr);
86   ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr);
87   ierr = MatDestroy(&A);CHKERRQ(ierr);
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