xref: /petsc/src/mat/tests/ex83.c (revision daa037dfd3c3bec8dc8659548d2b20b07c1dc6de)
1 
2 static char help[] = "Partition tiny grid using hierarchical partitioning and increase overlap using MatIncreaseOverlapSplit.\n\n";
3 
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>
13 
14 int main(int argc,char **args)
15 {
16   Mat             A,B;
17   PetscMPIInt     rank,size,membershipKey;
18   PetscInt        *ia,*ja,*indices_sc,isrows_localsize;
19   const PetscInt  *indices;
20   MatPartitioning part;
21   IS              is,isrows,isrows_sc;
22   IS              coarseparts,fineparts;
23   MPI_Comm        comm,scomm;
24 
25   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
26   comm = PETSC_COMM_WORLD;
27   PetscCallMPI(MPI_Comm_size(comm,&size));
28   PetscCheck(size == 4,comm,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 4 processors ");
29   PetscCallMPI(MPI_Comm_rank(comm,&rank));
30   /*set a small matrix */
31   PetscCall(PetscMalloc1(5,&ia));
32   PetscCall(PetscMalloc1(16,&ja));
33   if (rank == 0) {
34     ja[0] = 1; ja[1] = 4; ja[2] = 0; ja[3] = 2; ja[4] = 5; ja[5] = 1; ja[6] = 3; ja[7] = 6;
35     ja[8] = 2; ja[9] = 7;
36     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
37     membershipKey = 0;
38   } else if (rank == 1) {
39     ja[0] = 0; ja[1] = 5; ja[2] = 8; ja[3] = 1; ja[4] = 4; ja[5] = 6; ja[6] = 9; ja[7] = 2;
40     ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11;
41     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
42     membershipKey = 0;
43   } else if (rank == 2) {
44     ja[0] = 4; ja[1] = 9; ja[2] = 12; ja[3] = 5; ja[4] = 8; ja[5] = 10; ja[6] = 13; ja[7] = 6;
45     ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15;
46     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
47     membershipKey = 1;
48   } else {
49     ja[0] = 8; ja[1] = 13; ja[2] = 9; ja[3] = 12; ja[4] = 14; ja[5] = 10; ja[6] = 13; ja[7] = 15;
50     ja[8] = 11; ja[9] = 14;
51     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
52     membershipKey = 1;
53   }
54   PetscCall(MatCreateMPIAdj(comm,4,16,ia,ja,NULL,&A));
55   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
56   /*
57    Partition the graph of the matrix
58   */
59   PetscCall(MatPartitioningCreate(comm,&part));
60   PetscCall(MatPartitioningSetAdjacency(part,A));
61   PetscCall(MatPartitioningSetType(part,MATPARTITIONINGHIERARCH));
62   PetscCall(MatPartitioningHierarchicalSetNcoarseparts(part,2));
63   PetscCall(MatPartitioningHierarchicalSetNfineparts(part,2));
64   PetscCall(MatPartitioningSetFromOptions(part));
65   /* get new processor owner number of each vertex */
66   PetscCall(MatPartitioningApply(part,&is));
67   /* coarse parts */
68   PetscCall(MatPartitioningHierarchicalGetCoarseparts(part,&coarseparts));
69   PetscCall(ISView(coarseparts,PETSC_VIEWER_STDOUT_WORLD));
70   /* fine parts */
71   PetscCall(MatPartitioningHierarchicalGetFineparts(part,&fineparts));
72   PetscCall(ISView(fineparts,PETSC_VIEWER_STDOUT_WORLD));
73   /* partitioning */
74   PetscCall(ISView(is,PETSC_VIEWER_STDOUT_WORLD));
75   /* compute coming rows */
76   PetscCall(ISBuildTwoSided(is,NULL,&isrows));
77   PetscCall(ISView(isrows,PETSC_VIEWER_STDOUT_WORLD));
78   /*create a sub-communicator */
79   PetscCallMPI(MPI_Comm_split(comm, membershipKey,rank,&scomm));
80   PetscCall(ISGetLocalSize(isrows,&isrows_localsize));
81   PetscCall(PetscMalloc1(isrows_localsize,&indices_sc));
82   PetscCall(ISGetIndices(isrows,&indices));
83   PetscCall(PetscArraycpy(indices_sc,indices,isrows_localsize));
84   PetscCall(ISRestoreIndices(isrows,&indices));
85   PetscCall(ISDestroy(&is));
86   PetscCall(ISDestroy(&coarseparts));
87   PetscCall(ISDestroy(&fineparts));
88   PetscCall(ISDestroy(&isrows));
89   PetscCall(MatPartitioningDestroy(&part));
90   /*create a sub-IS on the sub communicator  */
91   PetscCall(ISCreateGeneral(scomm,isrows_localsize,indices_sc,PETSC_OWN_POINTER,&isrows_sc));
92   PetscCall(MatConvert(A,MATMPIAIJ,MAT_INITIAL_MATRIX,&B));
93 #if 1
94   PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
95 #endif
96   /*increase overlap */
97   PetscCall(MatIncreaseOverlapSplit(B,1,&isrows_sc,1));
98   PetscCall(ISView(isrows_sc,NULL));
99   PetscCall(ISDestroy(&isrows_sc));
100   /*
101     Free work space.  All PETSc objects should be destroyed when they
102     are no longer needed.
103   */
104   PetscCall(MatDestroy(&A));
105   PetscCall(MatDestroy(&B));
106   PetscCall(PetscFinalize());
107   return 0;
108 }
109