xref: /petsc/src/mat/tests/ex83.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Partition tiny grid using hierarchical partitioning and increase overlap using MatIncreaseOverlapSplit.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown   Include "petscmat.h" so that we can use matrices.  Note that this file
6c4762a1bSJed Brown   automatically includes:
7c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
8c4762a1bSJed Brown      petscmat.h - matrices
9c4762a1bSJed Brown      petscis.h     - index sets
10c4762a1bSJed Brown      petscviewer.h - viewers
11c4762a1bSJed Brown */
12c4762a1bSJed Brown #include <petscmat.h>
13c4762a1bSJed Brown 
14c4762a1bSJed Brown int main(int argc,char **args)
15c4762a1bSJed Brown {
16c4762a1bSJed Brown   Mat             A,B;
17c4762a1bSJed Brown   PetscMPIInt     rank,size,membershipKey;
18c4762a1bSJed Brown   PetscInt        *ia,*ja,*indices_sc,isrows_localsize;
19c4762a1bSJed Brown   const PetscInt  *indices;
20c4762a1bSJed Brown   MatPartitioning part;
21c4762a1bSJed Brown   IS              is,isrows,isrows_sc;
22c4762a1bSJed Brown   IS              coarseparts,fineparts;
23c4762a1bSJed Brown   MPI_Comm        comm,scomm;
24c4762a1bSJed Brown 
25*327415f7SBarry Smith   PetscFunctionBeginUser;
269566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
27c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
29be096a46SBarry Smith   PetscCheck(size == 4,comm,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 4 processors ");
309566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
31c4762a1bSJed Brown   /*set a small matrix */
329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(5,&ia));
339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(16,&ja));
34dd400576SPatrick Sanan   if (rank == 0) {
35c4762a1bSJed 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;
36c4762a1bSJed Brown     ja[8] = 2; ja[9] = 7;
37c4762a1bSJed Brown     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
38c4762a1bSJed Brown     membershipKey = 0;
39c4762a1bSJed Brown   } else if (rank == 1) {
40c4762a1bSJed 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;
41c4762a1bSJed Brown     ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11;
42c4762a1bSJed Brown     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
43c4762a1bSJed Brown     membershipKey = 0;
44c4762a1bSJed Brown   } else if (rank == 2) {
45c4762a1bSJed 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;
46c4762a1bSJed Brown     ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15;
47c4762a1bSJed Brown     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
48c4762a1bSJed Brown     membershipKey = 1;
49c4762a1bSJed Brown   } else {
50c4762a1bSJed 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;
51c4762a1bSJed Brown     ja[8] = 11; ja[9] = 14;
52c4762a1bSJed Brown     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
53c4762a1bSJed Brown     membershipKey = 1;
54c4762a1bSJed Brown   }
559566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIAdj(comm,4,16,ia,ja,NULL,&A));
569566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
57c4762a1bSJed Brown   /*
58c4762a1bSJed Brown    Partition the graph of the matrix
59c4762a1bSJed Brown   */
609566063dSJacob Faibussowitsch   PetscCall(MatPartitioningCreate(comm,&part));
619566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetAdjacency(part,A));
629566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetType(part,MATPARTITIONINGHIERARCH));
639566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalSetNcoarseparts(part,2));
649566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalSetNfineparts(part,2));
659566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetFromOptions(part));
66c4762a1bSJed Brown   /* get new processor owner number of each vertex */
679566063dSJacob Faibussowitsch   PetscCall(MatPartitioningApply(part,&is));
68c4762a1bSJed Brown   /* coarse parts */
699566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalGetCoarseparts(part,&coarseparts));
709566063dSJacob Faibussowitsch   PetscCall(ISView(coarseparts,PETSC_VIEWER_STDOUT_WORLD));
71c4762a1bSJed Brown   /* fine parts */
729566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalGetFineparts(part,&fineparts));
739566063dSJacob Faibussowitsch   PetscCall(ISView(fineparts,PETSC_VIEWER_STDOUT_WORLD));
74c4762a1bSJed Brown   /* partitioning */
759566063dSJacob Faibussowitsch   PetscCall(ISView(is,PETSC_VIEWER_STDOUT_WORLD));
76c4762a1bSJed Brown   /* compute coming rows */
779566063dSJacob Faibussowitsch   PetscCall(ISBuildTwoSided(is,NULL,&isrows));
789566063dSJacob Faibussowitsch   PetscCall(ISView(isrows,PETSC_VIEWER_STDOUT_WORLD));
79c4762a1bSJed Brown   /*create a sub-communicator */
809566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_split(comm, membershipKey,rank,&scomm));
819566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrows,&isrows_localsize));
829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(isrows_localsize,&indices_sc));
839566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrows,&indices));
849566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(indices_sc,indices,isrows_localsize));
859566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrows,&indices));
869566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
879566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&coarseparts));
889566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&fineparts));
899566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrows));
909566063dSJacob Faibussowitsch   PetscCall(MatPartitioningDestroy(&part));
91c4762a1bSJed Brown   /*create a sub-IS on the sub communicator  */
929566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(scomm,isrows_localsize,indices_sc,PETSC_OWN_POINTER,&isrows_sc));
939566063dSJacob Faibussowitsch   PetscCall(MatConvert(A,MATMPIAIJ,MAT_INITIAL_MATRIX,&B));
94c4762a1bSJed Brown #if 1
959566063dSJacob Faibussowitsch   PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
96c4762a1bSJed Brown #endif
97c4762a1bSJed Brown   /*increase overlap */
989566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlapSplit(B,1,&isrows_sc,1));
999566063dSJacob Faibussowitsch   PetscCall(ISView(isrows_sc,NULL));
1009566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrows_sc));
101c4762a1bSJed Brown   /*
102c4762a1bSJed Brown     Free work space.  All PETSc objects should be destroyed when they
103c4762a1bSJed Brown     are no longer needed.
104c4762a1bSJed Brown   */
1059566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1079566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
108b122ec5aSJacob Faibussowitsch   return 0;
109c4762a1bSJed Brown }
110