xref: /petsc/src/mat/tests/ex83.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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;
33ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
34*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 4,comm,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 4 processors ");
35ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
36c4762a1bSJed Brown   /*set a small matrix */
37c4762a1bSJed Brown   ierr = PetscMalloc1(5,&ia);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = PetscMalloc1(16,&ja);CHKERRQ(ierr);
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   }
60c4762a1bSJed Brown   ierr = MatCreateMPIAdj(comm,4,16,ia,ja,NULL,&A);CHKERRQ(ierr);
61c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
62c4762a1bSJed Brown   /*
63c4762a1bSJed Brown    Partition the graph of the matrix
64c4762a1bSJed Brown   */
65c4762a1bSJed Brown   ierr = MatPartitioningCreate(comm,&part);CHKERRQ(ierr);
66c4762a1bSJed Brown   ierr = MatPartitioningSetAdjacency(part,A);CHKERRQ(ierr);
67c4762a1bSJed Brown   ierr = MatPartitioningSetType(part,MATPARTITIONINGHIERARCH);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = MatPartitioningHierarchicalSetNcoarseparts(part,2);CHKERRQ(ierr);
69c4762a1bSJed Brown   ierr = MatPartitioningHierarchicalSetNfineparts(part,2);CHKERRQ(ierr);
70c4762a1bSJed Brown   ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr);
71c4762a1bSJed Brown   /* get new processor owner number of each vertex */
72c4762a1bSJed Brown   ierr = MatPartitioningApply(part,&is);CHKERRQ(ierr);
73c4762a1bSJed Brown   /* coarse parts */
74c4762a1bSJed Brown   ierr = MatPartitioningHierarchicalGetCoarseparts(part,&coarseparts);CHKERRQ(ierr);
75c4762a1bSJed Brown   ierr = ISView(coarseparts,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
76c4762a1bSJed Brown   /* fine parts */
77c4762a1bSJed Brown   ierr = MatPartitioningHierarchicalGetFineparts(part,&fineparts);CHKERRQ(ierr);
78c4762a1bSJed Brown   ierr = ISView(fineparts,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
79c4762a1bSJed Brown   /* partitioning */
80c4762a1bSJed Brown   ierr = ISView(is,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
81c4762a1bSJed Brown   /* compute coming rows */
82c4762a1bSJed Brown   ierr = ISBuildTwoSided(is,NULL,&isrows);CHKERRQ(ierr);
83c4762a1bSJed Brown   ierr = ISView(isrows,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
84c4762a1bSJed Brown   /*create a sub-communicator */
85ffc4695bSBarry Smith   ierr = MPI_Comm_split(comm, membershipKey,rank,&scomm);CHKERRMPI(ierr);
86c4762a1bSJed Brown   ierr = ISGetLocalSize(isrows,&isrows_localsize);CHKERRQ(ierr);
87c4762a1bSJed Brown   ierr = PetscMalloc1(isrows_localsize,&indices_sc);CHKERRQ(ierr);
88c4762a1bSJed Brown   ierr = ISGetIndices(isrows,&indices);CHKERRQ(ierr);
89c4762a1bSJed Brown   ierr = PetscArraycpy(indices_sc,indices,isrows_localsize);CHKERRQ(ierr);
90c4762a1bSJed Brown   ierr = ISRestoreIndices(isrows,&indices);CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = ISDestroy(&is);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = ISDestroy(&coarseparts);CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr = ISDestroy(&fineparts);CHKERRQ(ierr);
94c4762a1bSJed Brown   ierr = ISDestroy(&isrows);CHKERRQ(ierr);
95c4762a1bSJed Brown   ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr);
96c4762a1bSJed Brown   /*create a sub-IS on the sub communicator  */
97c4762a1bSJed Brown   ierr = ISCreateGeneral(scomm,isrows_localsize,indices_sc,PETSC_OWN_POINTER,&isrows_sc);CHKERRQ(ierr);
98c4762a1bSJed Brown   ierr = MatConvert(A,MATMPIAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
99c4762a1bSJed Brown #if 1
100c4762a1bSJed Brown   ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
101c4762a1bSJed Brown #endif
102c4762a1bSJed Brown   /*increase overlap */
103c4762a1bSJed Brown   ierr = MatIncreaseOverlapSplit(B,1,&isrows_sc,1);CHKERRQ(ierr);
104c4762a1bSJed Brown   ierr = ISView(isrows_sc,NULL);CHKERRQ(ierr);
105c4762a1bSJed Brown   ierr = ISDestroy(&isrows_sc);CHKERRQ(ierr);
106c4762a1bSJed Brown   /*
107c4762a1bSJed Brown     Free work space.  All PETSc objects should be destroyed when they
108c4762a1bSJed Brown     are no longer needed.
109c4762a1bSJed Brown   */
110c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
111c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
112c4762a1bSJed Brown   ierr = PetscFinalize();
113c4762a1bSJed Brown   return ierr;
114c4762a1bSJed Brown }
115c4762a1bSJed Brown 
116