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