xref: /petsc/src/mat/tests/ex193.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown  * ex193.c
3c4762a1bSJed Brown  *
4c4762a1bSJed Brown  *  Created on: Jul 29, 2015
5c4762a1bSJed Brown  *      Author: Fande Kong fdkong.jd@gmail.com
6c4762a1bSJed Brown  */
7c4762a1bSJed Brown /*
8c4762a1bSJed Brown  * An example demonstrates how to use hierarchical partitioning approach
9c4762a1bSJed Brown  */
10c4762a1bSJed Brown 
11c4762a1bSJed Brown #include <petscmat.h>
12c4762a1bSJed Brown 
13c4762a1bSJed Brown static char help[] = "Illustrates use of hierarchical partitioning.\n";
14c4762a1bSJed Brown 
15*9371c9d4SSatish Balay int main(int argc, char **args) {
16c4762a1bSJed Brown   Mat             A;    /* matrix */
17c4762a1bSJed Brown   PetscInt        m, n; /* mesh dimensions in x- and y- directions */
18c4762a1bSJed Brown   PetscInt        i, j, Ii, J, Istart, Iend;
19c4762a1bSJed Brown   PetscMPIInt     size;
20c4762a1bSJed Brown   PetscScalar     v;
21c4762a1bSJed Brown   MatPartitioning part;
22c4762a1bSJed Brown   IS              coarseparts, fineparts;
23c4762a1bSJed Brown   IS              is, isn, isrows;
24c4762a1bSJed Brown   MPI_Comm        comm;
25c4762a1bSJed Brown 
26327415f7SBarry Smith   PetscFunctionBeginUser;
279566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
28c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
30d0609cedSBarry Smith   PetscOptionsBegin(comm, NULL, "ex193", "hierarchical partitioning");
31c4762a1bSJed Brown   m = 15;
329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-M", "Number of mesh points in the x-direction", "partitioning", m, &m, NULL));
33c4762a1bSJed Brown   n = 15;
349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-N", "Number of mesh points in the y-direction", "partitioning", n, &n, NULL));
35d0609cedSBarry Smith   PetscOptionsEnd();
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /*
38c4762a1bSJed Brown      Assemble the matrix for the five point stencil (finite difference), YET AGAIN
39c4762a1bSJed Brown   */
409566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &A));
419566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
429566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
439566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
449566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &Istart, &Iend));
45c4762a1bSJed Brown   for (Ii = Istart; Ii < Iend; Ii++) {
46*9371c9d4SSatish Balay     v = -1.0;
47*9371c9d4SSatish Balay     i = Ii / n;
48*9371c9d4SSatish Balay     j = Ii - i * n;
49*9371c9d4SSatish Balay     if (i > 0) {
50*9371c9d4SSatish Balay       J = Ii - n;
51*9371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
52*9371c9d4SSatish Balay     }
53*9371c9d4SSatish Balay     if (i < m - 1) {
54*9371c9d4SSatish Balay       J = Ii + n;
55*9371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
56*9371c9d4SSatish Balay     }
57*9371c9d4SSatish Balay     if (j > 0) {
58*9371c9d4SSatish Balay       J = Ii - 1;
59*9371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
60*9371c9d4SSatish Balay     }
61*9371c9d4SSatish Balay     if (j < n - 1) {
62*9371c9d4SSatish Balay       J = Ii + 1;
63*9371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
64*9371c9d4SSatish Balay     }
65*9371c9d4SSatish Balay     v = 4.0;
66*9371c9d4SSatish Balay     PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
67c4762a1bSJed Brown   }
689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
709566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
71c4762a1bSJed Brown   /*
72c4762a1bSJed Brown    Partition the graph of the matrix
73c4762a1bSJed Brown   */
749566063dSJacob Faibussowitsch   PetscCall(MatPartitioningCreate(comm, &part));
759566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetAdjacency(part, A));
769566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetType(part, MATPARTITIONINGHIERARCH));
779566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalSetNcoarseparts(part, 2));
789566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalSetNfineparts(part, 4));
799566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetFromOptions(part));
80c4762a1bSJed Brown   /* get new processor owner number of each vertex */
819566063dSJacob Faibussowitsch   PetscCall(MatPartitioningApply(part, &is));
82c4762a1bSJed Brown   /* coarse parts */
839566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalGetCoarseparts(part, &coarseparts));
849566063dSJacob Faibussowitsch   PetscCall(ISView(coarseparts, PETSC_VIEWER_STDOUT_WORLD));
85c4762a1bSJed Brown   /* fine parts */
869566063dSJacob Faibussowitsch   PetscCall(MatPartitioningHierarchicalGetFineparts(part, &fineparts));
879566063dSJacob Faibussowitsch   PetscCall(ISView(fineparts, PETSC_VIEWER_STDOUT_WORLD));
88c4762a1bSJed Brown   /* partitioning */
899566063dSJacob Faibussowitsch   PetscCall(ISView(is, PETSC_VIEWER_STDOUT_WORLD));
90c4762a1bSJed Brown   /* get new global number of each old global number */
919566063dSJacob Faibussowitsch   PetscCall(ISPartitioningToNumbering(is, &isn));
929566063dSJacob Faibussowitsch   PetscCall(ISView(isn, PETSC_VIEWER_STDOUT_WORLD));
939566063dSJacob Faibussowitsch   PetscCall(ISBuildTwoSided(is, NULL, &isrows));
949566063dSJacob Faibussowitsch   PetscCall(ISView(isrows, PETSC_VIEWER_STDOUT_WORLD));
959566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
969566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&coarseparts));
979566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&fineparts));
989566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrows));
999566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isn));
1009566063dSJacob Faibussowitsch   PetscCall(MatPartitioningDestroy(&part));
1019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1029566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
103b122ec5aSJacob Faibussowitsch   return 0;
104c4762a1bSJed Brown }
105c4762a1bSJed Brown 
106c4762a1bSJed Brown /*TEST
107c4762a1bSJed Brown 
108c4762a1bSJed Brown    test:
109c4762a1bSJed Brown       nsize: 4
110c4762a1bSJed Brown       args: -mat_partitioning_hierarchical_Nfineparts 2
111c4762a1bSJed Brown       requires: parmetis
112c4762a1bSJed Brown       TODO: cannot run because parmetis does reproduce across all machines, probably due to nonportable random number generator
113c4762a1bSJed Brown 
114c4762a1bSJed Brown TEST*/
115