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; 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(comm, &G)); 49c4762a1bSJed Brown /* The identity matrix */ 50c4762a1bSJed Brown switch (testnum) { 51c4762a1bSJed Brown case 0: 52c4762a1bSJed Brown { 53c4762a1bSJed Brown Vec D; 54c4762a1bSJed Brown 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(G)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(G, &D, NULL)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(D, 1.0)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalSet(G, D, INSERT_VALUES)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&D)); 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 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(G)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(G, 2, NULL)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(G)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(G, &rStart, &rEnd)); 74c4762a1bSJed Brown row = 0; 75c4762a1bSJed Brown cols[0] = 0; cols[1] = 1; 765f80ce2aSJacob Faibussowitsch if ((row >= rStart) && (row < rEnd)) CHKERRQ(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES)); 77c4762a1bSJed Brown row = 1; 78c4762a1bSJed Brown cols[0] = 0; cols[1] = 1; 795f80ce2aSJacob Faibussowitsch if ((row >= rStart) && (row < rEnd)) CHKERRQ(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES)); 80c4762a1bSJed Brown row = 2; 81c4762a1bSJed Brown cols[0] = 2; cols[1] = 3; 825f80ce2aSJacob Faibussowitsch if ((row >= rStart) && (row < rEnd)) CHKERRQ(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES)); 83c4762a1bSJed Brown row = 3; 84c4762a1bSJed Brown cols[0] = 3; cols[1] = 4; 855f80ce2aSJacob Faibussowitsch if ((row >= rStart) && (row < rEnd)) CHKERRQ(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES)); 86c4762a1bSJed Brown row = 4; 87c4762a1bSJed Brown cols[0] = 4; cols[1] = 2; 885f80ce2aSJacob Faibussowitsch if ((row >= rStart) && (row < rEnd)) CHKERRQ(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY)); 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 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(G, PETSC_DETERMINE, PETSC_DETERMINE, 5, 5)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(G)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(G, 2, NULL)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(G)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(G, &rStart, &rEnd)); 104c4762a1bSJed Brown row = 0; 105c4762a1bSJed Brown cols[0] = 0; cols[1] = 4; 1065f80ce2aSJacob Faibussowitsch if ((row >= rStart) && (row < rEnd)) CHKERRQ(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES)); 107c4762a1bSJed Brown row = 1; 108c4762a1bSJed Brown cols[0] = 1; cols[1] = 2; 1095f80ce2aSJacob Faibussowitsch if ((row >= rStart) && (row < rEnd)) CHKERRQ(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES)); 110c4762a1bSJed Brown row = 2; 111c4762a1bSJed Brown cols[0] = 2; cols[1] = 3; 1125f80ce2aSJacob Faibussowitsch if ((row >= rStart) && (row < rEnd)) CHKERRQ(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES)); 113c4762a1bSJed Brown row = 3; 114c4762a1bSJed Brown cols[0] = 3; cols[1] = 1; 1155f80ce2aSJacob Faibussowitsch if ((row >= rStart) && (row < rEnd)) CHKERRQ(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES)); 116c4762a1bSJed Brown row = 4; 117c4762a1bSJed Brown cols[0] = 0; cols[1] = 4; 1185f80ce2aSJacob Faibussowitsch if ((row >= rStart) && (row < rEnd)) CHKERRQ(MatSetValues(G, 1, &row, 2, cols, vals, INSERT_VALUES)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY)); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown break; 123c4762a1bSJed Brown default: 12498921bdaSJacob 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 142*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL,help)); 143c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1445f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm, &size)); 145c4762a1bSJed Brown /* Use matrix to encode a graph */ 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-testnum", &testnum, NULL)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(CreateGraph(comm, testnum, &A)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 1515f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(A, MATOP_MULT, (void (*)) MatMultMax_SeqAIJ)); 152c4762a1bSJed Brown } else { 153c4762a1bSJed Brown Mat_MPIAIJ *a = (Mat_MPIAIJ *) A->data; 154c4762a1bSJed Brown 1555f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(a->A, MATOP_MULT, (void (*)) MatMultMax_SeqAIJ)); 1565f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(a->B, MATOP_MULT, (void (*)) MatMultMax_SeqAIJ)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(a->B, MATOP_MULT_ADD, (void (*)) MatMultAddMax_SeqAIJ)); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown /* Initialize each vertex as a separate component */ 1605f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A, &c, NULL)); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A, &vStart, &vEnd)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(c, &carray)); 163c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 164c4762a1bSJed Brown carray[v-vStart] = v; 165c4762a1bSJed Brown } 1665f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(c, &carray)); 167c4762a1bSJed Brown /* Preprocess in parallel to find local components */ 168c4762a1bSJed Brown /* Multiply until c does not change */ 1695f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(c, &cold)); 170c4762a1bSJed Brown for (v = 0; v < V; ++v) { 171c4762a1bSJed Brown Vec cnew = cold; 172c4762a1bSJed Brown PetscBool stop; 173c4762a1bSJed Brown 1745f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A, c, cnew)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(VecEqual(c, cnew, &stop)); 176c4762a1bSJed Brown if (stop) break; 177c4762a1bSJed Brown cold = c; 178c4762a1bSJed Brown c = cnew; 179c4762a1bSJed Brown } 180c4762a1bSJed Brown /* Report */ 1815f80ce2aSJacob Faibussowitsch CHKERRQ(VecUniqueEntries(c, &n, NULL)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(comm, "Components: %d Iterations: %d\n", n, v)); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(c, PETSC_VIEWER_STDOUT_WORLD)); 184c4762a1bSJed Brown /* Cleanup */ 1855f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&c)); 1865f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&cold)); 187*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 188*b122ec5aSJacob Faibussowitsch return 0; 189c4762a1bSJed Brown } 190