xref: /petsc/src/mat/tests/ex170.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 */
42c4762a1bSJed Brown PetscErrorCode CreateGraph(MPI_Comm comm, PetscInt testnum, Mat *A)
43c4762a1bSJed Brown {
44c4762a1bSJed Brown   Mat            G;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   PetscFunctionBegin;
479566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &G));
48c4762a1bSJed Brown   /* The identity matrix */
49c4762a1bSJed Brown   switch (testnum) {
50c4762a1bSJed Brown   case 0:
51c4762a1bSJed Brown   {
52c4762a1bSJed Brown     Vec D;
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5));
559566063dSJacob Faibussowitsch     PetscCall(MatSetUp(G));
569566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(G, &D, NULL));
579566063dSJacob Faibussowitsch     PetscCall(VecSet(D, 1.0));
589566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(G, D, INSERT_VALUES));
599566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&D));
60c4762a1bSJed Brown   }
61c4762a1bSJed Brown   break;
62c4762a1bSJed Brown   case 1:
63c4762a1bSJed Brown   {
64c4762a1bSJed Brown     PetscScalar vals[3] = {1.0, 1.0, 1.0};
65c4762a1bSJed Brown     PetscInt    cols[3];
66c4762a1bSJed Brown     PetscInt    rStart, rEnd, row;
67c4762a1bSJed Brown 
689566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5));
699566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(G));
709566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(G, 2, NULL));
719566063dSJacob Faibussowitsch     PetscCall(MatSetUp(G));
729566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(G, &rStart, &rEnd));
73c4762a1bSJed Brown     row  = 0;
74c4762a1bSJed Brown     cols[0] = 0; cols[1] = 1;
759566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
76c4762a1bSJed Brown     row  = 1;
77c4762a1bSJed Brown     cols[0] = 0; cols[1] = 1;
789566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
79c4762a1bSJed Brown     row  = 2;
80c4762a1bSJed Brown     cols[0] = 2; cols[1] = 3;
819566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
82c4762a1bSJed Brown     row  = 3;
83c4762a1bSJed Brown     cols[0] = 3; cols[1] = 4;
849566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
85c4762a1bSJed Brown     row  = 4;
86c4762a1bSJed Brown     cols[0] = 4; cols[1] = 2;
879566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
889566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
899566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));
90c4762a1bSJed Brown   }
91c4762a1bSJed Brown   break;
92c4762a1bSJed Brown   case 2:
93c4762a1bSJed Brown   {
94c4762a1bSJed Brown     PetscScalar vals[3] = {1.0, 1.0, 1.0};
95c4762a1bSJed Brown     PetscInt    cols[3];
96c4762a1bSJed Brown     PetscInt    rStart, rEnd, row;
97c4762a1bSJed Brown 
989566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5));
999566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(G));
1009566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(G, 2, NULL));
1019566063dSJacob Faibussowitsch     PetscCall(MatSetUp(G));
1029566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(G, &rStart, &rEnd));
103c4762a1bSJed Brown     row  = 0;
104c4762a1bSJed Brown     cols[0] = 0; cols[1] = 4;
1059566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
106c4762a1bSJed Brown     row  = 1;
107c4762a1bSJed Brown     cols[0] = 1; cols[1] = 2;
1089566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
109c4762a1bSJed Brown     row  = 2;
110c4762a1bSJed Brown     cols[0] = 2; cols[1] = 3;
1119566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
112c4762a1bSJed Brown     row  = 3;
113c4762a1bSJed Brown     cols[0] = 3; cols[1] = 1;
1149566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
115c4762a1bSJed Brown     row  = 4;
116c4762a1bSJed Brown     cols[0] = 0; cols[1] = 4;
1179566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
1189566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
1199566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));
120c4762a1bSJed Brown   }
121c4762a1bSJed Brown   break;
122c4762a1bSJed Brown   default:
12398921bdaSJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_PLIB, "Unknown test %d", testnum);
124c4762a1bSJed Brown   }
125c4762a1bSJed Brown   *A = G;
126c4762a1bSJed Brown   PetscFunctionReturn(0);
127c4762a1bSJed Brown }
128c4762a1bSJed Brown 
129c4762a1bSJed Brown int main(int argc, char **argv)
130c4762a1bSJed Brown {
131c4762a1bSJed Brown   MPI_Comm       comm;
132c4762a1bSJed Brown   Mat            A;    /* A graph */
133c4762a1bSJed Brown   Vec            c;    /* A vector giving the component of each vertex */
134c4762a1bSJed Brown   Vec            cold; /* The vector c from the last iteration */
135c4762a1bSJed Brown   PetscScalar   *carray;
136c4762a1bSJed Brown   PetscInt       testnum = 0;
137c4762a1bSJed Brown   PetscInt       V, vStart, vEnd, v, n;
138c4762a1bSJed Brown   PetscMPIInt    size;
139c4762a1bSJed Brown 
140*327415f7SBarry Smith   PetscFunctionBeginUser;
1419566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
142c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
1439566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
144c4762a1bSJed Brown   /* Use matrix to encode a graph */
1459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL, "-testnum", &testnum, NULL));
1469566063dSJacob Faibussowitsch   PetscCall(CreateGraph(comm, testnum, &A));
1479566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &V, NULL));
148c4762a1bSJed Brown   /* Replace matrix-vector multiplication with one that calculates the minimum rather than the sum */
149c4762a1bSJed Brown   if (size == 1) {
1509566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)) MatMultMax_SeqAIJ));
151c4762a1bSJed Brown   } else {
152c4762a1bSJed Brown     Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
153c4762a1bSJed Brown 
1549566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(a->A, MATOP_MULT, (void (*)) MatMultMax_SeqAIJ));
1559566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(a->B, MATOP_MULT, (void (*)) MatMultMax_SeqAIJ));
1569566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(a->B, MATOP_MULT_ADD, (void (*)) MatMultAddMax_SeqAIJ));
157c4762a1bSJed Brown   }
158c4762a1bSJed Brown   /* Initialize each vertex as a separate component */
1599566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &c, NULL));
1609566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &vStart, &vEnd));
1619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(c, &carray));
162c4762a1bSJed Brown   for (v = vStart; v < vEnd; ++v) {
163c4762a1bSJed Brown     carray[v-vStart] = v;
164c4762a1bSJed Brown   }
1659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(c, &carray));
166c4762a1bSJed Brown   /* Preprocess in parallel to find local components */
167c4762a1bSJed Brown   /* Multiply until c does not change */
1689566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(c, &cold));
169c4762a1bSJed Brown   for (v = 0; v < V; ++v) {
170c4762a1bSJed Brown     Vec       cnew = cold;
171c4762a1bSJed Brown     PetscBool stop;
172c4762a1bSJed Brown 
1739566063dSJacob Faibussowitsch     PetscCall(MatMult(A, c, cnew));
1749566063dSJacob Faibussowitsch     PetscCall(VecEqual(c, cnew, &stop));
175c4762a1bSJed Brown     if (stop) break;
176c4762a1bSJed Brown     cold = c;
177c4762a1bSJed Brown     c    = cnew;
178c4762a1bSJed Brown   }
179c4762a1bSJed Brown   /* Report */
1809566063dSJacob Faibussowitsch   PetscCall(VecUniqueEntries(c, &n, NULL));
1819566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Components: %d Iterations: %d\n", n, v));
1829566063dSJacob Faibussowitsch   PetscCall(VecView(c, PETSC_VIEWER_STDOUT_WORLD));
183c4762a1bSJed Brown   /* Cleanup */
1849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&c));
1859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cold));
1869566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
187b122ec5aSJacob Faibussowitsch   return 0;
188c4762a1bSJed Brown }
189