xref: /petsc/src/mat/tests/ex170.c (revision 613ce9fe8f5e2bcdf7c72d914e9769b5d960fb4c)
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 
13*1d27aa22SBarry Smith   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 */
CreateGraph(MPI_Comm comm,PetscInt testnum,Mat * A)42d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateGraph(MPI_Comm comm, PetscInt testnum, Mat *A)
43d71ae5a4SJacob Faibussowitsch {
44c4762a1bSJed Brown   Mat G;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   PetscFunctionBegin;
479566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &G));
48c4762a1bSJed Brown   /* The identity matrix */
49c4762a1bSJed Brown   switch (testnum) {
509371c9d4SSatish Balay   case 0: {
51c4762a1bSJed Brown     Vec D;
52c4762a1bSJed Brown 
539566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5));
549566063dSJacob Faibussowitsch     PetscCall(MatSetUp(G));
559566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(G, &D, NULL));
569566063dSJacob Faibussowitsch     PetscCall(VecSet(D, 1.0));
579566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(G, D, INSERT_VALUES));
589566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&D));
599371c9d4SSatish Balay   } break;
609371c9d4SSatish Balay   case 1: {
61c4762a1bSJed Brown     PetscScalar vals[3] = {1.0, 1.0, 1.0};
62c4762a1bSJed Brown     PetscInt    cols[3];
63c4762a1bSJed Brown     PetscInt    rStart, rEnd, row;
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5));
669566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(G));
679566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(G, 2, NULL));
689566063dSJacob Faibussowitsch     PetscCall(MatSetUp(G));
699566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(G, &rStart, &rEnd));
70c4762a1bSJed Brown     row     = 0;
719371c9d4SSatish Balay     cols[0] = 0;
729371c9d4SSatish Balay     cols[1] = 1;
739566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
74c4762a1bSJed Brown     row     = 1;
759371c9d4SSatish Balay     cols[0] = 0;
769371c9d4SSatish Balay     cols[1] = 1;
779566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
78c4762a1bSJed Brown     row     = 2;
799371c9d4SSatish Balay     cols[0] = 2;
809371c9d4SSatish Balay     cols[1] = 3;
819566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
82c4762a1bSJed Brown     row     = 3;
839371c9d4SSatish Balay     cols[0] = 3;
849371c9d4SSatish Balay     cols[1] = 4;
859566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
86c4762a1bSJed Brown     row     = 4;
879371c9d4SSatish Balay     cols[0] = 4;
889371c9d4SSatish Balay     cols[1] = 2;
899566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
909566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
919566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));
929371c9d4SSatish Balay   } break;
939371c9d4SSatish Balay   case 2: {
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;
1049371c9d4SSatish Balay     cols[0] = 0;
1059371c9d4SSatish Balay     cols[1] = 4;
1069566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
107c4762a1bSJed Brown     row     = 1;
1089371c9d4SSatish Balay     cols[0] = 1;
1099371c9d4SSatish Balay     cols[1] = 2;
1109566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
111c4762a1bSJed Brown     row     = 2;
1129371c9d4SSatish Balay     cols[0] = 2;
1139371c9d4SSatish Balay     cols[1] = 3;
1149566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
115c4762a1bSJed Brown     row     = 3;
1169371c9d4SSatish Balay     cols[0] = 3;
1179371c9d4SSatish Balay     cols[1] = 1;
1189566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
119c4762a1bSJed Brown     row     = 4;
1209371c9d4SSatish Balay     cols[0] = 0;
1219371c9d4SSatish Balay     cols[1] = 4;
1229566063dSJacob Faibussowitsch     if ((row >= rStart) && (row < rEnd)) PetscCall(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES));
1239566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY));
1249566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY));
1259371c9d4SSatish Balay   } break;
126d71ae5a4SJacob Faibussowitsch   default:
127d71ae5a4SJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_PLIB, "Unknown test %d", testnum);
128c4762a1bSJed Brown   }
129c4762a1bSJed Brown   *A = G;
1303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
131c4762a1bSJed Brown }
132c4762a1bSJed Brown 
main(int argc,char ** argv)133d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
134d71ae5a4SJacob Faibussowitsch {
135c4762a1bSJed Brown   MPI_Comm     comm;
136c4762a1bSJed Brown   Mat          A;    /* A graph */
137c4762a1bSJed Brown   Vec          c;    /* A vector giving the component of each vertex */
138c4762a1bSJed Brown   Vec          cold; /* The vector c from the last iteration */
139c4762a1bSJed Brown   PetscScalar *carray;
140c4762a1bSJed Brown   PetscInt     testnum = 0;
141c4762a1bSJed Brown   PetscInt     V, vStart, vEnd, v, n;
142c4762a1bSJed Brown   PetscMPIInt  size;
143c4762a1bSJed Brown 
144327415f7SBarry Smith   PetscFunctionBeginUser;
1459566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
146c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
1479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
148c4762a1bSJed Brown   /* Use matrix to encode a graph */
1499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-testnum", &testnum, NULL));
1509566063dSJacob Faibussowitsch   PetscCall(CreateGraph(comm, testnum, &A));
1519566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &V, NULL));
152c4762a1bSJed Brown   /* Replace matrix-vector multiplication with one that calculates the minimum rather than the sum */
153c4762a1bSJed Brown   if (size == 1) {
1549566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(A, MATOP_MULT, (void(*))MatMultMax_SeqAIJ));
155c4762a1bSJed Brown   } else {
156c4762a1bSJed Brown     Mat_MPIAIJ *a = (Mat_MPIAIJ *)A->data;
157c4762a1bSJed Brown 
1589566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(a->A, MATOP_MULT, (void(*))MatMultMax_SeqAIJ));
1599566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(a->B, MATOP_MULT, (void(*))MatMultMax_SeqAIJ));
1609566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(a->B, MATOP_MULT_ADD, (void(*))MatMultAddMax_SeqAIJ));
161c4762a1bSJed Brown   }
162c4762a1bSJed Brown   /* Initialize each vertex as a separate component */
1639566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &c, NULL));
1649566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &vStart, &vEnd));
1659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(c, &carray));
166ad540459SPierre Jolivet   for (v = vStart; v < vEnd; ++v) carray[v - vStart] = v;
1679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(c, &carray));
168c4762a1bSJed Brown   /* Preprocess in parallel to find local components */
169c4762a1bSJed Brown   /* Multiply until c does not change */
1709566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(c, &cold));
171c4762a1bSJed Brown   for (v = 0; v < V; ++v) {
172c4762a1bSJed Brown     Vec       cnew = cold;
173c4762a1bSJed Brown     PetscBool stop;
174c4762a1bSJed Brown 
1759566063dSJacob Faibussowitsch     PetscCall(MatMult(A, c, cnew));
1769566063dSJacob Faibussowitsch     PetscCall(VecEqual(c, cnew, &stop));
177c4762a1bSJed Brown     if (stop) break;
178c4762a1bSJed Brown     cold = c;
179c4762a1bSJed Brown     c    = cnew;
180c4762a1bSJed Brown   }
181c4762a1bSJed Brown   /* Report */
1829566063dSJacob Faibussowitsch   PetscCall(VecUniqueEntries(c, &n, NULL));
1839566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Components: %d Iterations: %d\n", n, v));
1849566063dSJacob Faibussowitsch   PetscCall(VecView(c, PETSC_VIEWER_STDOUT_WORLD));
185c4762a1bSJed Brown   /* Cleanup */
1869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&c));
1879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cold));
1889566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
189b122ec5aSJacob Faibussowitsch   return 0;
190c4762a1bSJed Brown }
191