xref: /petsc/src/mat/tests/ex170.c (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
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   PetscErrorCode ierr;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   PetscFunctionBegin;
48c4762a1bSJed Brown   ierr = MatCreate(comm, &G);CHKERRQ(ierr);
49c4762a1bSJed Brown   /* The identity matrix */
50c4762a1bSJed Brown   switch (testnum) {
51c4762a1bSJed Brown   case 0:
52c4762a1bSJed Brown   {
53c4762a1bSJed Brown     Vec D;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown     ierr = MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5);CHKERRQ(ierr);
56c4762a1bSJed Brown     ierr = MatSetUp(G);CHKERRQ(ierr);
57c4762a1bSJed Brown     ierr = MatCreateVecs(G, &D, NULL);CHKERRQ(ierr);
58c4762a1bSJed Brown     ierr = VecSet(D, 1.0);CHKERRQ(ierr);
59c4762a1bSJed Brown     ierr = MatDiagonalSet(G, D, INSERT_VALUES);CHKERRQ(ierr);
60c4762a1bSJed Brown     ierr = VecDestroy(&D);CHKERRQ(ierr);
61c4762a1bSJed Brown   }
62c4762a1bSJed Brown   break;
63c4762a1bSJed Brown   case 1:
64c4762a1bSJed Brown   {
65c4762a1bSJed Brown     PetscScalar vals[3] = {1.0, 1.0, 1.0};
66c4762a1bSJed Brown     PetscInt    cols[3];
67c4762a1bSJed Brown     PetscInt    rStart, rEnd, row;
68c4762a1bSJed Brown 
69c4762a1bSJed Brown     ierr = MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5);CHKERRQ(ierr);
70c4762a1bSJed Brown     ierr = MatSetFromOptions(G);CHKERRQ(ierr);
71c4762a1bSJed Brown     ierr = MatSeqAIJSetPreallocation(G, 2, NULL);CHKERRQ(ierr);
72c4762a1bSJed Brown     ierr = MatSetUp(G);CHKERRQ(ierr);
73c4762a1bSJed Brown     ierr = MatGetOwnershipRange(G, &rStart, &rEnd);CHKERRQ(ierr);
74c4762a1bSJed Brown     row  = 0;
75c4762a1bSJed Brown     cols[0] = 0; cols[1] = 1;
76c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
77c4762a1bSJed Brown     row  = 1;
78c4762a1bSJed Brown     cols[0] = 0; cols[1] = 1;
79c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
80c4762a1bSJed Brown     row  = 2;
81c4762a1bSJed Brown     cols[0] = 2; cols[1] = 3;
82c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
83c4762a1bSJed Brown     row  = 3;
84c4762a1bSJed Brown     cols[0] = 3; cols[1] = 4;
85c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
86c4762a1bSJed Brown     row  = 4;
87c4762a1bSJed Brown     cols[0] = 4; cols[1] = 2;
88c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
89c4762a1bSJed Brown     ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90c4762a1bSJed Brown     ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown   break;
93c4762a1bSJed Brown   case 2:
94c4762a1bSJed Brown   {
95c4762a1bSJed Brown     PetscScalar vals[3] = {1.0, 1.0, 1.0};
96c4762a1bSJed Brown     PetscInt    cols[3];
97c4762a1bSJed Brown     PetscInt    rStart, rEnd, row;
98c4762a1bSJed Brown 
99c4762a1bSJed Brown     ierr = MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5);CHKERRQ(ierr);
100c4762a1bSJed Brown     ierr = MatSetFromOptions(G);CHKERRQ(ierr);
101c4762a1bSJed Brown     ierr = MatSeqAIJSetPreallocation(G, 2, NULL);CHKERRQ(ierr);
102c4762a1bSJed Brown     ierr = MatSetUp(G);CHKERRQ(ierr);
103c4762a1bSJed Brown     ierr = MatGetOwnershipRange(G, &rStart, &rEnd);CHKERRQ(ierr);
104c4762a1bSJed Brown     row  = 0;
105c4762a1bSJed Brown     cols[0] = 0; cols[1] = 4;
106c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
107c4762a1bSJed Brown     row  = 1;
108c4762a1bSJed Brown     cols[0] = 1; cols[1] = 2;
109c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
110c4762a1bSJed Brown     row  = 2;
111c4762a1bSJed Brown     cols[0] = 2; cols[1] = 3;
112c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
113c4762a1bSJed Brown     row  = 3;
114c4762a1bSJed Brown     cols[0] = 3; cols[1] = 1;
115c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
116c4762a1bSJed Brown     row  = 4;
117c4762a1bSJed Brown     cols[0] = 0; cols[1] = 4;
118c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
119c4762a1bSJed Brown     ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
120c4762a1bSJed Brown     ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
121c4762a1bSJed Brown   }
122c4762a1bSJed Brown   break;
123c4762a1bSJed Brown   default:
124*98921bdaSJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_PLIB, "Unknown test %d", testnum);
125c4762a1bSJed Brown   }
126c4762a1bSJed Brown   *A = G;
127c4762a1bSJed Brown   PetscFunctionReturn(0);
128c4762a1bSJed Brown }
129c4762a1bSJed Brown 
130c4762a1bSJed Brown int main(int argc, char **argv)
131c4762a1bSJed Brown {
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   PetscErrorCode ierr;
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
143c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
144ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
145c4762a1bSJed Brown   /* Use matrix to encode a graph */
146c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL, "-testnum", &testnum, NULL);CHKERRQ(ierr);
147c4762a1bSJed Brown   ierr = CreateGraph(comm, testnum, &A);CHKERRQ(ierr);
148c4762a1bSJed Brown   ierr = MatGetSize(A, &V, NULL);CHKERRQ(ierr);
149c4762a1bSJed Brown   /* Replace matrix-vector multiplication with one that calculates the minimum rather than the sum */
150c4762a1bSJed Brown   if (size == 1) {
151c4762a1bSJed Brown     ierr = MatShellSetOperation(A, MATOP_MULT, (void (*)) MatMultMax_SeqAIJ);CHKERRQ(ierr);
152c4762a1bSJed Brown   } else {
153c4762a1bSJed Brown     Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
154c4762a1bSJed Brown 
155c4762a1bSJed Brown     ierr = MatShellSetOperation(a->A, MATOP_MULT, (void (*)) MatMultMax_SeqAIJ);CHKERRQ(ierr);
156c4762a1bSJed Brown     ierr = MatShellSetOperation(a->B, MATOP_MULT, (void (*)) MatMultMax_SeqAIJ);CHKERRQ(ierr);
157c4762a1bSJed Brown     ierr = MatShellSetOperation(a->B, MATOP_MULT_ADD, (void (*)) MatMultAddMax_SeqAIJ);CHKERRQ(ierr);
158c4762a1bSJed Brown   }
159c4762a1bSJed Brown   /* Initialize each vertex as a separate component */
160c4762a1bSJed Brown   ierr = MatCreateVecs(A, &c, NULL);CHKERRQ(ierr);
161c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A, &vStart, &vEnd);CHKERRQ(ierr);
162c4762a1bSJed Brown   ierr = VecGetArray(c, &carray);CHKERRQ(ierr);
163c4762a1bSJed Brown   for (v = vStart; v < vEnd; ++v) {
164c4762a1bSJed Brown     carray[v-vStart] = v;
165c4762a1bSJed Brown   }
166c4762a1bSJed Brown   ierr = VecRestoreArray(c, &carray);CHKERRQ(ierr);
167c4762a1bSJed Brown   /* Preprocess in parallel to find local components */
168c4762a1bSJed Brown   /* Multiply until c does not change */
169c4762a1bSJed Brown   ierr = VecDuplicate(c, &cold);CHKERRQ(ierr);
170c4762a1bSJed Brown   for (v = 0; v < V; ++v) {
171c4762a1bSJed Brown     Vec       cnew = cold;
172c4762a1bSJed Brown     PetscBool stop;
173c4762a1bSJed Brown 
174c4762a1bSJed Brown     ierr = MatMult(A, c, cnew);CHKERRQ(ierr);
175c4762a1bSJed Brown     ierr = VecEqual(c, cnew, &stop);CHKERRQ(ierr);
176c4762a1bSJed Brown     if (stop) break;
177c4762a1bSJed Brown     cold = c;
178c4762a1bSJed Brown     c    = cnew;
179c4762a1bSJed Brown   }
180c4762a1bSJed Brown   /* Report */
181c4762a1bSJed Brown   ierr = VecUniqueEntries(c, &n, NULL);CHKERRQ(ierr);
182c4762a1bSJed Brown   ierr = PetscPrintf(comm, "Components: %d Iterations: %d\n", n, v);CHKERRQ(ierr);
183c4762a1bSJed Brown   ierr = VecView(c, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
184c4762a1bSJed Brown   /* Cleanup */
185c4762a1bSJed Brown   ierr = VecDestroy(&c);CHKERRQ(ierr);
186c4762a1bSJed Brown   ierr = VecDestroy(&cold);CHKERRQ(ierr);
187c4762a1bSJed Brown   ierr = PetscFinalize();
188c4762a1bSJed Brown   return ierr;
189c4762a1bSJed Brown }
190