xref: /petsc/src/mat/tests/ex83.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
14*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
15*d71ae5a4SJacob Faibussowitsch {
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 
25327415f7SBarry 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) {
359371c9d4SSatish Balay     ja[0]         = 1;
369371c9d4SSatish Balay     ja[1]         = 4;
379371c9d4SSatish Balay     ja[2]         = 0;
389371c9d4SSatish Balay     ja[3]         = 2;
399371c9d4SSatish Balay     ja[4]         = 5;
409371c9d4SSatish Balay     ja[5]         = 1;
419371c9d4SSatish Balay     ja[6]         = 3;
429371c9d4SSatish Balay     ja[7]         = 6;
439371c9d4SSatish Balay     ja[8]         = 2;
449371c9d4SSatish Balay     ja[9]         = 7;
459371c9d4SSatish Balay     ia[0]         = 0;
469371c9d4SSatish Balay     ia[1]         = 2;
479371c9d4SSatish Balay     ia[2]         = 5;
489371c9d4SSatish Balay     ia[3]         = 8;
499371c9d4SSatish Balay     ia[4]         = 10;
50c4762a1bSJed Brown     membershipKey = 0;
51c4762a1bSJed Brown   } else if (rank == 1) {
529371c9d4SSatish Balay     ja[0]         = 0;
539371c9d4SSatish Balay     ja[1]         = 5;
549371c9d4SSatish Balay     ja[2]         = 8;
559371c9d4SSatish Balay     ja[3]         = 1;
569371c9d4SSatish Balay     ja[4]         = 4;
579371c9d4SSatish Balay     ja[5]         = 6;
589371c9d4SSatish Balay     ja[6]         = 9;
599371c9d4SSatish Balay     ja[7]         = 2;
609371c9d4SSatish Balay     ja[8]         = 5;
619371c9d4SSatish Balay     ja[9]         = 7;
629371c9d4SSatish Balay     ja[10]        = 10;
639371c9d4SSatish Balay     ja[11]        = 3;
649371c9d4SSatish Balay     ja[12]        = 6;
659371c9d4SSatish Balay     ja[13]        = 11;
669371c9d4SSatish Balay     ia[0]         = 0;
679371c9d4SSatish Balay     ia[1]         = 3;
689371c9d4SSatish Balay     ia[2]         = 7;
699371c9d4SSatish Balay     ia[3]         = 11;
709371c9d4SSatish Balay     ia[4]         = 14;
71c4762a1bSJed Brown     membershipKey = 0;
72c4762a1bSJed Brown   } else if (rank == 2) {
739371c9d4SSatish Balay     ja[0]         = 4;
749371c9d4SSatish Balay     ja[1]         = 9;
759371c9d4SSatish Balay     ja[2]         = 12;
769371c9d4SSatish Balay     ja[3]         = 5;
779371c9d4SSatish Balay     ja[4]         = 8;
789371c9d4SSatish Balay     ja[5]         = 10;
799371c9d4SSatish Balay     ja[6]         = 13;
809371c9d4SSatish Balay     ja[7]         = 6;
819371c9d4SSatish Balay     ja[8]         = 9;
829371c9d4SSatish Balay     ja[9]         = 11;
839371c9d4SSatish Balay     ja[10]        = 14;
849371c9d4SSatish Balay     ja[11]        = 7;
859371c9d4SSatish Balay     ja[12]        = 10;
869371c9d4SSatish Balay     ja[13]        = 15;
879371c9d4SSatish Balay     ia[0]         = 0;
889371c9d4SSatish Balay     ia[1]         = 3;
899371c9d4SSatish Balay     ia[2]         = 7;
909371c9d4SSatish Balay     ia[3]         = 11;
919371c9d4SSatish Balay     ia[4]         = 14;
92c4762a1bSJed Brown     membershipKey = 1;
93c4762a1bSJed Brown   } else {
949371c9d4SSatish Balay     ja[0]         = 8;
959371c9d4SSatish Balay     ja[1]         = 13;
969371c9d4SSatish Balay     ja[2]         = 9;
979371c9d4SSatish Balay     ja[3]         = 12;
989371c9d4SSatish Balay     ja[4]         = 14;
999371c9d4SSatish Balay     ja[5]         = 10;
1009371c9d4SSatish Balay     ja[6]         = 13;
1019371c9d4SSatish Balay     ja[7]         = 15;
1029371c9d4SSatish Balay     ja[8]         = 11;
1039371c9d4SSatish Balay     ja[9]         = 14;
1049371c9d4SSatish Balay     ia[0]         = 0;
1059371c9d4SSatish Balay     ia[1]         = 2;
1069371c9d4SSatish Balay     ia[2]         = 5;
1079371c9d4SSatish Balay     ia[3]         = 8;
1089371c9d4SSatish Balay     ia[4]         = 10;
109c4762a1bSJed Brown     membershipKey = 1;
110c4762a1bSJed Brown   }
1119566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIAdj(comm, 4, 16, ia, ja, NULL, &A));
1129566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
113c4762a1bSJed Brown   /*
114c4762a1bSJed Brown    Partition the graph of the matrix
115c4762a1bSJed Brown   */
1169566063dSJacob Faibussowitsch   PetscCall(MatPartitioningCreate(comm, &part));
1179566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetAdjacency(part, A));
1189566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetType(part, MATPARTITIONINGHIERARCH));
1199566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalSetNcoarseparts(part, 2));
1209566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalSetNfineparts(part, 2));
1219566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetFromOptions(part));
122c4762a1bSJed Brown   /* get new processor owner number of each vertex */
1239566063dSJacob Faibussowitsch   PetscCall(MatPartitioningApply(part, &is));
124c4762a1bSJed Brown   /* coarse parts */
1259566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalGetCoarseparts(part, &coarseparts));
1269566063dSJacob Faibussowitsch   PetscCall(ISView(coarseparts, PETSC_VIEWER_STDOUT_WORLD));
127c4762a1bSJed Brown   /* fine parts */
1289566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalGetFineparts(part, &fineparts));
1299566063dSJacob Faibussowitsch   PetscCall(ISView(fineparts, PETSC_VIEWER_STDOUT_WORLD));
130c4762a1bSJed Brown   /* partitioning */
1319566063dSJacob Faibussowitsch   PetscCall(ISView(is, PETSC_VIEWER_STDOUT_WORLD));
132c4762a1bSJed Brown   /* compute coming rows */
1339566063dSJacob Faibussowitsch   PetscCall(ISBuildTwoSided(is, NULL, &isrows));
1349566063dSJacob Faibussowitsch   PetscCall(ISView(isrows, PETSC_VIEWER_STDOUT_WORLD));
135c4762a1bSJed Brown   /*create a sub-communicator */
1369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_split(comm, membershipKey, rank, &scomm));
1379566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrows, &isrows_localsize));
1389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(isrows_localsize, &indices_sc));
1399566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrows, &indices));
1409566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(indices_sc, indices, isrows_localsize));
1419566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrows, &indices));
1429566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
1439566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&coarseparts));
1449566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&fineparts));
1459566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrows));
1469566063dSJacob Faibussowitsch   PetscCall(MatPartitioningDestroy(&part));
147c4762a1bSJed Brown   /*create a sub-IS on the sub communicator  */
1489566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(scomm, isrows_localsize, indices_sc, PETSC_OWN_POINTER, &isrows_sc));
1499566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATMPIAIJ, MAT_INITIAL_MATRIX, &B));
150c4762a1bSJed Brown #if 1
1519566063dSJacob Faibussowitsch   PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD));
152c4762a1bSJed Brown #endif
153c4762a1bSJed Brown   /*increase overlap */
1549566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlapSplit(B, 1, &isrows_sc, 1));
1559566063dSJacob Faibussowitsch   PetscCall(ISView(isrows_sc, NULL));
1569566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrows_sc));
157c4762a1bSJed Brown   /*
158c4762a1bSJed Brown     Free work space.  All PETSc objects should be destroyed when they
159c4762a1bSJed Brown     are no longer needed.
160c4762a1bSJed Brown   */
1619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1639566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
164b122ec5aSJacob Faibussowitsch   return 0;
165c4762a1bSJed Brown }
166