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