xref: /petsc/src/mat/tests/ex170.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Scalable algorithm for Connected Components problem.\n\
2c4762a1bSJed Brown Entails changing the MatMult() for this matrix.\n\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode MatMultMax_SeqAIJ(Mat, Vec, Vec);
7c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode MatMultAddMax_SeqAIJ(Mat, Vec, Vec, Vec);
8c4762a1bSJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
9c4762a1bSJed Brown 
10c4762a1bSJed Brown /*
11c4762a1bSJed Brown   Paper with Ananth: Frbenius norm of band was good proxy, but really want to know the rank outside
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   LU for diagonal blocks must do shifting instead of pivoting, preferably shifting individual rows (like Pardiso)
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   Draw picture of flow of reordering
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   Measure Forbenius norm of the blocks being dropped by Truncated SPIKE (might be contaminated by pivoting in LU)
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   Report on using Florida matrices (Maxim, Murat)
20c4762a1bSJed Brown */
21c4762a1bSJed Brown 
22c4762a1bSJed Brown /*
23c4762a1bSJed Brown I have thought about how to do this. Here is a prototype algorithm. Let A be
24c4762a1bSJed Brown the adjacency matrix (0 or 1), and let each component be identified by the
25c4762a1bSJed Brown lowest numbered vertex in it. We initialize a vector c so that each vertex is
26c4762a1bSJed Brown a component, c_i = i. Now we act on c with A, using a special product
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   c = A * c
29c4762a1bSJed Brown 
30c4762a1bSJed Brown where we replace addition with min. The fixed point of this operation is a vector
31c4762a1bSJed Brown c which is the component for each vertex. The number of iterates is
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   max_{components} depth of BFS tree for component
34c4762a1bSJed Brown 
35c4762a1bSJed Brown We can accelerate this algorithm by preprocessing all locals domains using the
36c4762a1bSJed Brown same algorithm. Then the number of iterations is bounded the depth of the BFS
37c4762a1bSJed Brown tree for the graph on supervertices defined over local components, which is
38c4762a1bSJed Brown bounded by p. In practice, this should be very fast.
39c4762a1bSJed Brown */
40c4762a1bSJed Brown 
41c4762a1bSJed Brown /* Only isolated vertices get a 1 on the diagonal */
42*9371c9d4SSatish Balay PetscErrorCode CreateGraph(MPI_Comm comm, PetscInt testnum, Mat *A) {
43c4762a1bSJed Brown   Mat G;
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   PetscFunctionBegin;
469566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &G));
47c4762a1bSJed Brown   /* The identity matrix */
48c4762a1bSJed Brown   switch (testnum) {
49*9371c9d4SSatish Balay   case 0: {
50c4762a1bSJed Brown     Vec D;
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5));
539566063dSJacob Faibussowitsch     PetscCall(MatSetUp(G));
549566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(G, &D, NULL));
559566063dSJacob Faibussowitsch     PetscCall(VecSet(D, 1.0));
569566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(G, D, INSERT_VALUES));
579566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&D));
58*9371c9d4SSatish Balay   } break;
59*9371c9d4SSatish Balay   case 1: {
60c4762a1bSJed Brown     PetscScalar vals[3] = {1.0, 1.0, 1.0};
61c4762a1bSJed Brown     PetscInt    cols[3];
62c4762a1bSJed Brown     PetscInt    rStart, rEnd, row;
63c4762a1bSJed Brown 
649566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5));
659566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(G));
669566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(G, 2, NULL));
679566063dSJacob Faibussowitsch     PetscCall(MatSetUp(G));
689566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(G, &rStart, &rEnd));
69c4762a1bSJed Brown     row     = 0;
70*9371c9d4SSatish Balay     cols[0] = 0;
71*9371c9d4SSatish Balay     cols[1] = 1;
729566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
73c4762a1bSJed Brown     row     = 1;
74*9371c9d4SSatish Balay     cols[0] = 0;
75*9371c9d4SSatish Balay     cols[1] = 1;
769566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
77c4762a1bSJed Brown     row     = 2;
78*9371c9d4SSatish Balay     cols[0] = 2;
79*9371c9d4SSatish Balay     cols[1] = 3;
809566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
81c4762a1bSJed Brown     row     = 3;
82*9371c9d4SSatish Balay     cols[0] = 3;
83*9371c9d4SSatish Balay     cols[1] = 4;
849566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
85c4762a1bSJed Brown     row     = 4;
86*9371c9d4SSatish Balay     cols[0] = 4;
87*9371c9d4SSatish Balay     cols[1] = 2;
889566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
899566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
909566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));
91*9371c9d4SSatish Balay   } break;
92*9371c9d4SSatish Balay   case 2: {
93c4762a1bSJed Brown     PetscScalar vals[3] = {1.0, 1.0, 1.0};
94c4762a1bSJed Brown     PetscInt    cols[3];
95c4762a1bSJed Brown     PetscInt    rStart, rEnd, row;
96c4762a1bSJed Brown 
979566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5));
989566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(G));
999566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(G, 2, NULL));
1009566063dSJacob Faibussowitsch     PetscCall(MatSetUp(G));
1019566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(G, &rStart, &rEnd));
102c4762a1bSJed Brown     row     = 0;
103*9371c9d4SSatish Balay     cols[0] = 0;
104*9371c9d4SSatish Balay     cols[1] = 4;
1059566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
106c4762a1bSJed Brown     row     = 1;
107*9371c9d4SSatish Balay     cols[0] = 1;
108*9371c9d4SSatish Balay     cols[1] = 2;
1099566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
110c4762a1bSJed Brown     row     = 2;
111*9371c9d4SSatish Balay     cols[0] = 2;
112*9371c9d4SSatish Balay     cols[1] = 3;
1139566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
114c4762a1bSJed Brown     row     = 3;
115*9371c9d4SSatish Balay     cols[0] = 3;
116*9371c9d4SSatish Balay     cols[1] = 1;
1179566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
118c4762a1bSJed Brown     row     = 4;
119*9371c9d4SSatish Balay     cols[0] = 0;
120*9371c9d4SSatish Balay     cols[1] = 4;
1219566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
1229566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
1239566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));
124*9371c9d4SSatish Balay   } break;
125*9371c9d4SSatish Balay   default: SETERRQ(comm, PETSC_ERR_PLIB, "Unknown test %d", testnum);
126c4762a1bSJed Brown   }
127c4762a1bSJed Brown   *A = G;
128c4762a1bSJed Brown   PetscFunctionReturn(0);
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131*9371c9d4SSatish Balay int main(int argc, char **argv) {
132c4762a1bSJed Brown   MPI_Comm     comm;
133c4762a1bSJed Brown   Mat          A;    /* A graph */
134c4762a1bSJed Brown   Vec          c;    /* A vector giving the component of each vertex */
135c4762a1bSJed Brown   Vec          cold; /* The vector c from the last iteration */
136c4762a1bSJed Brown   PetscScalar *carray;
137c4762a1bSJed Brown   PetscInt     testnum = 0;
138c4762a1bSJed Brown   PetscInt     V, vStart, vEnd, v, n;
139c4762a1bSJed Brown   PetscMPIInt  size;
140c4762a1bSJed Brown 
141327415f7SBarry Smith   PetscFunctionBeginUser;
1429566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
143c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
1449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
145c4762a1bSJed Brown   /* Use matrix to encode a graph */
1469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-testnum", &testnum, NULL));
1479566063dSJacob Faibussowitsch   PetscCall(CreateGraph(comm, testnum, &A));
1489566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &V, NULL));
149c4762a1bSJed Brown   /* Replace matrix-vector multiplication with one that calculates the minimum rather than the sum */
150c4762a1bSJed Brown   if (size == 1) {
1519566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(A, MATOP_MULT, (void(*))MatMultMax_SeqAIJ));
152c4762a1bSJed Brown   } else {
153c4762a1bSJed Brown     Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data;
154c4762a1bSJed Brown 
1559566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(a->A, MATOP_MULT, (void(*))MatMultMax_SeqAIJ));
1569566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(a->B, MATOP_MULT, (void(*))MatMultMax_SeqAIJ));
1579566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(a->B, MATOP_MULT_ADD, (void(*))MatMultAddMax_SeqAIJ));
158c4762a1bSJed Brown   }
159c4762a1bSJed Brown   /* Initialize each vertex as a separate component */
1609566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &c, NULL));
1619566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &vStart, &vEnd));
1629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(c, &carray));
163*9371c9d4SSatish Balay   for (v = vStart; v < vEnd; ++v) { carray[v - vStart] = v; }
1649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(c, &carray));
165c4762a1bSJed Brown   /* Preprocess in parallel to find local components */
166c4762a1bSJed Brown   /* Multiply until c does not change */
1679566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(c, &cold));
168c4762a1bSJed Brown   for (v = 0; v < V; ++v) {
169c4762a1bSJed Brown     Vec       cnew = cold;
170c4762a1bSJed Brown     PetscBool stop;
171c4762a1bSJed Brown 
1729566063dSJacob Faibussowitsch     PetscCall(MatMult(A, c, cnew));
1739566063dSJacob Faibussowitsch     PetscCall(VecEqual(c, cnew, &stop));
174c4762a1bSJed Brown     if (stop) break;
175c4762a1bSJed Brown     cold = c;
176c4762a1bSJed Brown     c    = cnew;
177c4762a1bSJed Brown   }
178c4762a1bSJed Brown   /* Report */
1799566063dSJacob Faibussowitsch   PetscCall(VecUniqueEntries(c, &n, NULL));
1809566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Components: %d Iterations: %d\n", n, v));
1819566063dSJacob Faibussowitsch   PetscCall(VecView(c, PETSC_VIEWER_STDOUT_WORLD));
182c4762a1bSJed Brown   /* Cleanup */
1839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&c));
1849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cold));
1859566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
186b122ec5aSJacob Faibussowitsch   return 0;
187c4762a1bSJed Brown }
188