xref: /petsc/src/mat/tests/ex170.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Scalable algorithm for Connected Components problem.\n\
2*c4762a1bSJed Brown Entails changing the MatMult() for this matrix.\n\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode MatMultMax_SeqAIJ(Mat,Vec,Vec);
7*c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode MatMultAddMax_SeqAIJ(Mat,Vec,Vec,Vec);
8*c4762a1bSJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown /*
11*c4762a1bSJed Brown   Paper with Ananth: Frbenius norm of band was good proxy, but really want to know the rank outside
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown   LU for diagonal blocks must do shifting instead of pivoting, preferably shifting individual rows (like Pardiso)
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown   Draw picture of flow of reordering
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown   Measure Forbenius norm of the blocks being dropped by Truncated SPIKE (might be contaminated by pivoting in LU)
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown   Report on using Florida matrices (Maxim, Murat)
20*c4762a1bSJed Brown */
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown /*
23*c4762a1bSJed Brown I have thought about how to do this. Here is a prototype algorithm. Let A be
24*c4762a1bSJed Brown the adjacency matrix (0 or 1), and let each component be identified by the
25*c4762a1bSJed Brown lowest numbered vertex in it. We initialize a vector c so that each vertex is
26*c4762a1bSJed Brown a component, c_i = i. Now we act on c with A, using a special product
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown   c = A * c
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown where we replace addition with min. The fixed point of this operation is a vector
31*c4762a1bSJed Brown c which is the component for each vertex. The number of iterates is
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown   max_{components} depth of BFS tree for component
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown We can accelerate this algorithm by preprocessing all locals domains using the
36*c4762a1bSJed Brown same algorithm. Then the number of iterations is bounded the depth of the BFS
37*c4762a1bSJed Brown tree for the graph on supervertices defined over local components, which is
38*c4762a1bSJed Brown bounded by p. In practice, this should be very fast.
39*c4762a1bSJed Brown */
40*c4762a1bSJed Brown 
41*c4762a1bSJed Brown /* Only isolated vertices get a 1 on the diagonal */
42*c4762a1bSJed Brown PetscErrorCode CreateGraph(MPI_Comm comm, PetscInt testnum, Mat *A)
43*c4762a1bSJed Brown {
44*c4762a1bSJed Brown   Mat            G;
45*c4762a1bSJed Brown   PetscErrorCode ierr;
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown   PetscFunctionBegin;
48*c4762a1bSJed Brown   ierr = MatCreate(comm, &G);CHKERRQ(ierr);
49*c4762a1bSJed Brown   /* The identity matrix */
50*c4762a1bSJed Brown   switch (testnum) {
51*c4762a1bSJed Brown   case 0:
52*c4762a1bSJed Brown   {
53*c4762a1bSJed Brown     Vec D;
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown     ierr = MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5);CHKERRQ(ierr);
56*c4762a1bSJed Brown     ierr = MatSetUp(G);CHKERRQ(ierr);
57*c4762a1bSJed Brown     ierr = MatCreateVecs(G, &D, NULL);CHKERRQ(ierr);
58*c4762a1bSJed Brown     ierr = VecSet(D, 1.0);CHKERRQ(ierr);
59*c4762a1bSJed Brown     ierr = MatDiagonalSet(G, D, INSERT_VALUES);CHKERRQ(ierr);
60*c4762a1bSJed Brown     ierr = VecDestroy(&D);CHKERRQ(ierr);
61*c4762a1bSJed Brown   }
62*c4762a1bSJed Brown   break;
63*c4762a1bSJed Brown   case 1:
64*c4762a1bSJed Brown   {
65*c4762a1bSJed Brown     PetscScalar vals[3] = {1.0, 1.0, 1.0};
66*c4762a1bSJed Brown     PetscInt    cols[3];
67*c4762a1bSJed Brown     PetscInt    rStart, rEnd, row;
68*c4762a1bSJed Brown 
69*c4762a1bSJed Brown     ierr = MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5);CHKERRQ(ierr);
70*c4762a1bSJed Brown     ierr = MatSetFromOptions(G);CHKERRQ(ierr);
71*c4762a1bSJed Brown     ierr = MatSeqAIJSetPreallocation(G, 2, NULL);CHKERRQ(ierr);
72*c4762a1bSJed Brown     ierr = MatSetUp(G);CHKERRQ(ierr);
73*c4762a1bSJed Brown     ierr = MatGetOwnershipRange(G, &rStart, &rEnd);CHKERRQ(ierr);
74*c4762a1bSJed Brown     row  = 0;
75*c4762a1bSJed Brown     cols[0] = 0; cols[1] = 1;
76*c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
77*c4762a1bSJed Brown     row  = 1;
78*c4762a1bSJed Brown     cols[0] = 0; cols[1] = 1;
79*c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
80*c4762a1bSJed Brown     row  = 2;
81*c4762a1bSJed Brown     cols[0] = 2; cols[1] = 3;
82*c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
83*c4762a1bSJed Brown     row  = 3;
84*c4762a1bSJed Brown     cols[0] = 3; cols[1] = 4;
85*c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
86*c4762a1bSJed Brown     row  = 4;
87*c4762a1bSJed Brown     cols[0] = 4; cols[1] = 2;
88*c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
89*c4762a1bSJed Brown     ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90*c4762a1bSJed Brown     ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91*c4762a1bSJed Brown   }
92*c4762a1bSJed Brown   break;
93*c4762a1bSJed Brown   case 2:
94*c4762a1bSJed Brown   {
95*c4762a1bSJed Brown     PetscScalar vals[3] = {1.0, 1.0, 1.0};
96*c4762a1bSJed Brown     PetscInt    cols[3];
97*c4762a1bSJed Brown     PetscInt    rStart, rEnd, row;
98*c4762a1bSJed Brown 
99*c4762a1bSJed Brown     ierr = MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5);CHKERRQ(ierr);
100*c4762a1bSJed Brown     ierr = MatSetFromOptions(G);CHKERRQ(ierr);
101*c4762a1bSJed Brown     ierr = MatSeqAIJSetPreallocation(G, 2, NULL);CHKERRQ(ierr);
102*c4762a1bSJed Brown     ierr = MatSetUp(G);CHKERRQ(ierr);
103*c4762a1bSJed Brown     ierr = MatGetOwnershipRange(G, &rStart, &rEnd);CHKERRQ(ierr);
104*c4762a1bSJed Brown     row  = 0;
105*c4762a1bSJed Brown     cols[0] = 0; cols[1] = 4;
106*c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
107*c4762a1bSJed Brown     row  = 1;
108*c4762a1bSJed Brown     cols[0] = 1; cols[1] = 2;
109*c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
110*c4762a1bSJed Brown     row  = 2;
111*c4762a1bSJed Brown     cols[0] = 2; cols[1] = 3;
112*c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
113*c4762a1bSJed Brown     row  = 3;
114*c4762a1bSJed Brown     cols[0] = 3; cols[1] = 1;
115*c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
116*c4762a1bSJed Brown     row  = 4;
117*c4762a1bSJed Brown     cols[0] = 0; cols[1] = 4;
118*c4762a1bSJed Brown     if ((row >= rStart) && (row < rEnd)) {ierr = MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);}
119*c4762a1bSJed Brown     ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
120*c4762a1bSJed Brown     ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
121*c4762a1bSJed Brown   }
122*c4762a1bSJed Brown   break;
123*c4762a1bSJed Brown   default:
124*c4762a1bSJed Brown     SETERRQ1(comm, PETSC_ERR_PLIB, "Unknown test %d", testnum);
125*c4762a1bSJed Brown   }
126*c4762a1bSJed Brown   *A = G;
127*c4762a1bSJed Brown   PetscFunctionReturn(0);
128*c4762a1bSJed Brown }
129*c4762a1bSJed Brown 
130*c4762a1bSJed Brown int main(int argc, char **argv)
131*c4762a1bSJed Brown {
132*c4762a1bSJed Brown   MPI_Comm       comm;
133*c4762a1bSJed Brown   Mat            A;    /* A graph */
134*c4762a1bSJed Brown   Vec            c;    /* A vector giving the component of each vertex */
135*c4762a1bSJed Brown   Vec            cold; /* The vector c from the last iteration */
136*c4762a1bSJed Brown   PetscScalar   *carray;
137*c4762a1bSJed Brown   PetscInt       testnum = 0;
138*c4762a1bSJed Brown   PetscInt       V, vStart, vEnd, v, n;
139*c4762a1bSJed Brown   PetscMPIInt    size;
140*c4762a1bSJed Brown   PetscErrorCode ierr;
141*c4762a1bSJed Brown 
142*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
143*c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
144*c4762a1bSJed Brown   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
145*c4762a1bSJed Brown   /* Use matrix to encode a graph */
146*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL, "-testnum", &testnum, NULL);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = CreateGraph(comm, testnum, &A);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = MatGetSize(A, &V, NULL);CHKERRQ(ierr);
149*c4762a1bSJed Brown   /* Replace matrix-vector multiplication with one that calculates the minimum rather than the sum */
150*c4762a1bSJed Brown   if (size == 1) {
151*c4762a1bSJed Brown     ierr = MatShellSetOperation(A, MATOP_MULT, (void (*)) MatMultMax_SeqAIJ);CHKERRQ(ierr);
152*c4762a1bSJed Brown   } else {
153*c4762a1bSJed Brown     Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data;
154*c4762a1bSJed Brown 
155*c4762a1bSJed Brown     ierr = MatShellSetOperation(a->A, MATOP_MULT, (void (*)) MatMultMax_SeqAIJ);CHKERRQ(ierr);
156*c4762a1bSJed Brown     ierr = MatShellSetOperation(a->B, MATOP_MULT, (void (*)) MatMultMax_SeqAIJ);CHKERRQ(ierr);
157*c4762a1bSJed Brown     ierr = MatShellSetOperation(a->B, MATOP_MULT_ADD, (void (*)) MatMultAddMax_SeqAIJ);CHKERRQ(ierr);
158*c4762a1bSJed Brown   }
159*c4762a1bSJed Brown   /* Initialize each vertex as a separate component */
160*c4762a1bSJed Brown   ierr = MatCreateVecs(A, &c, NULL);CHKERRQ(ierr);
161*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A, &vStart, &vEnd);CHKERRQ(ierr);
162*c4762a1bSJed Brown   ierr = VecGetArray(c, &carray);CHKERRQ(ierr);
163*c4762a1bSJed Brown   for (v = vStart; v < vEnd; ++v) {
164*c4762a1bSJed Brown     carray[v-vStart] = v;
165*c4762a1bSJed Brown   }
166*c4762a1bSJed Brown   ierr = VecRestoreArray(c, &carray);CHKERRQ(ierr);
167*c4762a1bSJed Brown   /* Preprocess in parallel to find local components */
168*c4762a1bSJed Brown   /* Multiply until c does not change */
169*c4762a1bSJed Brown   ierr = VecDuplicate(c, &cold);CHKERRQ(ierr);
170*c4762a1bSJed Brown   for (v = 0; v < V; ++v) {
171*c4762a1bSJed Brown     Vec       cnew = cold;
172*c4762a1bSJed Brown     PetscBool stop;
173*c4762a1bSJed Brown 
174*c4762a1bSJed Brown     ierr = MatMult(A, c, cnew);CHKERRQ(ierr);
175*c4762a1bSJed Brown     ierr = VecEqual(c, cnew, &stop);CHKERRQ(ierr);
176*c4762a1bSJed Brown     if (stop) break;
177*c4762a1bSJed Brown     cold = c;
178*c4762a1bSJed Brown     c    = cnew;
179*c4762a1bSJed Brown   }
180*c4762a1bSJed Brown   /* Report */
181*c4762a1bSJed Brown   ierr = VecUniqueEntries(c, &n, NULL);CHKERRQ(ierr);
182*c4762a1bSJed Brown   ierr = PetscPrintf(comm, "Components: %d Iterations: %d\n", n, v);CHKERRQ(ierr);
183*c4762a1bSJed Brown   ierr = VecView(c, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
184*c4762a1bSJed Brown   /* Cleanup */
185*c4762a1bSJed Brown   ierr = VecDestroy(&c);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = VecDestroy(&cold);CHKERRQ(ierr);
187*c4762a1bSJed Brown   ierr = PetscFinalize();
188*c4762a1bSJed Brown   return ierr;
189*c4762a1bSJed Brown }
190