18be712e4SBarry Smith #include <petscmat.h> /*I <petscmat.h> I*/ 28be712e4SBarry Smith #include <petscblaslapack.h> 38be712e4SBarry Smith 48be712e4SBarry Smith /*@ 58be712e4SBarry Smith MatCreateLaplacian - Create the matrix Laplacian, with all values in the matrix less than the tolerance set to zero 68be712e4SBarry Smith 78be712e4SBarry Smith Input Parameters: 88be712e4SBarry Smith + A - The matrix 98be712e4SBarry Smith . tol - The zero tolerance 108be712e4SBarry Smith - weighted - Flag for using edge weights 118be712e4SBarry Smith 128be712e4SBarry Smith Output Parameter: 138be712e4SBarry Smith . L - The graph Laplacian matrix 148be712e4SBarry Smith 158be712e4SBarry Smith Level: intermediate 168be712e4SBarry Smith 178be712e4SBarry Smith .seealso: `MatFilter()`, `MatGetGraph()` 188be712e4SBarry Smith @*/ 198be712e4SBarry Smith PetscErrorCode MatCreateLaplacian(Mat A, PetscReal tol, PetscBool weighted, Mat *L) 208be712e4SBarry Smith { 218be712e4SBarry Smith PetscScalar *newVals; 228be712e4SBarry Smith PetscInt *newCols; 238be712e4SBarry Smith PetscInt rStart, rEnd, r, colMax = 0; 248be712e4SBarry Smith PetscInt *dnnz, *onnz; 258be712e4SBarry Smith PetscInt m, n, M, N; 268be712e4SBarry Smith 278be712e4SBarry Smith PetscFunctionBegin; 288be712e4SBarry Smith PetscCheck(!weighted, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Will get to this soon"); 298be712e4SBarry Smith PetscCall(MatCreate(PetscObjectComm((PetscObject)A), L)); 308be712e4SBarry Smith PetscCall(MatGetSize(A, &M, &N)); 318be712e4SBarry Smith PetscCall(MatGetLocalSize(A, &m, &n)); 328be712e4SBarry Smith PetscCall(MatSetSizes(*L, m, n, M, N)); 338be712e4SBarry Smith PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd)); 348be712e4SBarry Smith PetscCall(PetscMalloc2(m, &dnnz, m, &onnz)); 358be712e4SBarry Smith for (r = rStart; r < rEnd; ++r) { 368be712e4SBarry Smith const PetscScalar *vals; 378be712e4SBarry Smith const PetscInt *cols; 388be712e4SBarry Smith PetscInt ncols, newcols, c; 398be712e4SBarry Smith PetscBool hasdiag = PETSC_FALSE; 408be712e4SBarry Smith 418be712e4SBarry Smith dnnz[r - rStart] = onnz[r - rStart] = 0; 428be712e4SBarry Smith PetscCall(MatGetRow(A, r, &ncols, &cols, &vals)); 438be712e4SBarry Smith for (c = 0, newcols = 0; c < ncols; ++c) { 448be712e4SBarry Smith if (cols[c] == r) { 458be712e4SBarry Smith ++newcols; 468be712e4SBarry Smith hasdiag = PETSC_TRUE; 478be712e4SBarry Smith ++dnnz[r - rStart]; 488be712e4SBarry Smith } else if (PetscAbsScalar(vals[c]) >= tol) { 498be712e4SBarry Smith if ((cols[c] >= rStart) && (cols[c] < rEnd)) ++dnnz[r - rStart]; 508be712e4SBarry Smith else ++onnz[r - rStart]; 518be712e4SBarry Smith ++newcols; 528be712e4SBarry Smith } 538be712e4SBarry Smith } 548be712e4SBarry Smith if (!hasdiag) { 558be712e4SBarry Smith ++newcols; 568be712e4SBarry Smith ++dnnz[r - rStart]; 578be712e4SBarry Smith } 588be712e4SBarry Smith colMax = PetscMax(colMax, newcols); 598be712e4SBarry Smith PetscCall(MatRestoreRow(A, r, &ncols, &cols, &vals)); 608be712e4SBarry Smith } 618be712e4SBarry Smith PetscCall(MatSetFromOptions(*L)); 628be712e4SBarry Smith PetscCall(MatXAIJSetPreallocation(*L, 1, dnnz, onnz, NULL, NULL)); 638be712e4SBarry Smith PetscCall(MatSetUp(*L)); 648be712e4SBarry Smith PetscCall(PetscMalloc2(colMax, &newCols, colMax, &newVals)); 658be712e4SBarry Smith for (r = rStart; r < rEnd; ++r) { 668be712e4SBarry Smith const PetscScalar *vals; 678be712e4SBarry Smith const PetscInt *cols; 688be712e4SBarry Smith PetscInt ncols, newcols, c; 698be712e4SBarry Smith PetscBool hasdiag = PETSC_FALSE; 708be712e4SBarry Smith 718be712e4SBarry Smith PetscCall(MatGetRow(A, r, &ncols, &cols, &vals)); 728be712e4SBarry Smith for (c = 0, newcols = 0; c < ncols; ++c) { 738be712e4SBarry Smith if (cols[c] == r) { 748be712e4SBarry Smith newCols[newcols] = cols[c]; 758be712e4SBarry Smith newVals[newcols] = dnnz[r - rStart] + onnz[r - rStart] - 1; 768be712e4SBarry Smith ++newcols; 778be712e4SBarry Smith hasdiag = PETSC_TRUE; 788be712e4SBarry Smith } else if (PetscAbsScalar(vals[c]) >= tol) { 798be712e4SBarry Smith newCols[newcols] = cols[c]; 808be712e4SBarry Smith newVals[newcols] = -1.0; 818be712e4SBarry Smith ++newcols; 828be712e4SBarry Smith } 838be712e4SBarry Smith PetscCheck(newcols <= colMax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Overran work space"); 848be712e4SBarry Smith } 858be712e4SBarry Smith if (!hasdiag) { 868be712e4SBarry Smith newCols[newcols] = r; 878be712e4SBarry Smith newVals[newcols] = dnnz[r - rStart] + onnz[r - rStart] - 1; 888be712e4SBarry Smith ++newcols; 898be712e4SBarry Smith } 908be712e4SBarry Smith PetscCall(MatRestoreRow(A, r, &ncols, &cols, &vals)); 918be712e4SBarry Smith PetscCall(MatSetValues(*L, 1, &r, newcols, newCols, newVals, INSERT_VALUES)); 928be712e4SBarry Smith } 938be712e4SBarry Smith PetscCall(PetscFree2(dnnz, onnz)); 948be712e4SBarry Smith PetscCall(MatAssemblyBegin(*L, MAT_FINAL_ASSEMBLY)); 958be712e4SBarry Smith PetscCall(MatAssemblyEnd(*L, MAT_FINAL_ASSEMBLY)); 968be712e4SBarry Smith PetscCall(PetscFree2(newCols, newVals)); 978be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 988be712e4SBarry Smith } 998be712e4SBarry Smith 1008be712e4SBarry Smith /* 1018be712e4SBarry Smith MatGetOrdering_Spectral - Find the symmetric reordering of the graph by . 1028be712e4SBarry Smith */ 1038be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_Spectral(Mat A, MatOrderingType type, IS *row, IS *col) 1048be712e4SBarry Smith { 1058be712e4SBarry Smith Mat L; 1068be712e4SBarry Smith const PetscReal eps = 1.0e-12; 1078be712e4SBarry Smith 1088be712e4SBarry Smith PetscFunctionBegin; 1098be712e4SBarry Smith PetscCall(MatCreateLaplacian(A, eps, PETSC_FALSE, &L)); 1108be712e4SBarry Smith { 1118be712e4SBarry Smith /* Check Laplacian */ 1128be712e4SBarry Smith PetscReal norm; 1138be712e4SBarry Smith Vec x, y; 1148be712e4SBarry Smith 1158be712e4SBarry Smith PetscCall(MatCreateVecs(L, &x, NULL)); 1168be712e4SBarry Smith PetscCall(VecDuplicate(x, &y)); 1178be712e4SBarry Smith PetscCall(VecSet(x, 1.0)); 1188be712e4SBarry Smith PetscCall(MatMult(L, x, y)); 1198be712e4SBarry Smith PetscCall(VecNorm(y, NORM_INFINITY, &norm)); 1208be712e4SBarry Smith PetscCheck(norm <= 1.0e-10, PetscObjectComm((PetscObject)y), PETSC_ERR_PLIB, "Invalid graph Laplacian"); 1218be712e4SBarry Smith PetscCall(VecDestroy(&x)); 1228be712e4SBarry Smith PetscCall(VecDestroy(&y)); 1238be712e4SBarry Smith } 1248be712e4SBarry Smith /* Compute Fiedler vector (right now, all eigenvectors) */ 1258be712e4SBarry Smith #if defined(PETSC_USE_COMPLEX) 1268be712e4SBarry Smith SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Spectral partitioning does not support complex numbers"); 1278be712e4SBarry Smith #else 1288be712e4SBarry Smith { 1298be712e4SBarry Smith Mat LD; 1308be712e4SBarry Smith PetscScalar *a; 1318be712e4SBarry Smith PetscReal *realpart, *imagpart, *eigvec, *work; 1328be712e4SBarry Smith PetscReal sdummy; 1338be712e4SBarry Smith PetscBLASInt bn, bN, lwork = 0, lierr, idummy; 1348be712e4SBarry Smith PetscInt n, i, evInd, *perm, tmp; 1358be712e4SBarry Smith 1368be712e4SBarry Smith PetscCall(MatConvert(L, MATDENSE, MAT_INITIAL_MATRIX, &LD)); 1378be712e4SBarry Smith PetscCall(MatGetLocalSize(LD, &n, NULL)); 1388be712e4SBarry Smith PetscCall(MatDenseGetArray(LD, &a)); 1398be712e4SBarry Smith PetscCall(PetscBLASIntCast(n, &bn)); 1408be712e4SBarry Smith PetscCall(PetscBLASIntCast(n, &bN)); 1418be712e4SBarry Smith PetscCall(PetscBLASIntCast(5 * n, &lwork)); 1428be712e4SBarry Smith PetscCall(PetscBLASIntCast(1, &idummy)); 1438be712e4SBarry Smith PetscCall(PetscMalloc4(n, &realpart, n, &imagpart, n * n, &eigvec, lwork, &work)); 1448be712e4SBarry Smith PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 1458be712e4SBarry Smith PetscCallBLAS("LAPACKgeev", LAPACKgeev_("N", "V", &bn, a, &bN, realpart, imagpart, &sdummy, &idummy, eigvec, &bN, work, &lwork, &lierr)); 1468be712e4SBarry Smith PetscCheck(!lierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr); 1478be712e4SBarry Smith PetscCall(PetscFPTrapPop()); 1488be712e4SBarry Smith PetscCall(MatDenseRestoreArray(LD, &a)); 1498be712e4SBarry Smith PetscCall(MatDestroy(&LD)); 1508be712e4SBarry Smith /* Check lowest eigenvalue and eigenvector */ 1518be712e4SBarry Smith PetscCall(PetscMalloc1(n, &perm)); 1528be712e4SBarry Smith for (i = 0; i < n; ++i) perm[i] = i; 1538be712e4SBarry Smith PetscCall(PetscSortRealWithPermutation(n, realpart, perm)); 1548be712e4SBarry Smith evInd = perm[0]; 1558be712e4SBarry Smith PetscCheck(realpart[evInd] <= 1.0e-12 && imagpart[evInd] <= 1.0e-12, PetscObjectComm((PetscObject)L), PETSC_ERR_PLIB, "Graph Laplacian must have lowest eigenvalue 0"); 1568be712e4SBarry Smith evInd = perm[1]; 1578be712e4SBarry Smith PetscCheck(realpart[evInd] >= 1.0e-12 || imagpart[evInd] >= 1.0e-12, PetscObjectComm((PetscObject)L), PETSC_ERR_PLIB, "Graph Laplacian must have only one zero eigenvalue"); 1588be712e4SBarry Smith evInd = perm[0]; 1598be712e4SBarry Smith for (i = 0; i < n; ++i) { 160*f4f49eeaSPierre Jolivet PetscCheck(PetscAbsReal(eigvec[evInd * n + i] - eigvec[evInd * n + 0]) <= 1.0e-10, PetscObjectComm((PetscObject)L), PETSC_ERR_PLIB, "Graph Laplacian must have constant lowest eigenvector ev_%" PetscInt_FMT " %g != ev_0 %g", i, (double)eigvec[evInd * n + i], (double)eigvec[evInd * n + 0]); 1618be712e4SBarry Smith } 1628be712e4SBarry Smith /* Construct Fiedler partition */ 1638be712e4SBarry Smith evInd = perm[1]; 1648be712e4SBarry Smith for (i = 0; i < n; ++i) perm[i] = i; 1658be712e4SBarry Smith PetscCall(PetscSortRealWithPermutation(n, &eigvec[evInd * n], perm)); 1668be712e4SBarry Smith for (i = 0; i < n / 2; ++i) { 1678be712e4SBarry Smith tmp = perm[n - 1 - i]; 1688be712e4SBarry Smith perm[n - 1 - i] = perm[i]; 1698be712e4SBarry Smith perm[i] = tmp; 1708be712e4SBarry Smith } 1718be712e4SBarry Smith PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, perm, PETSC_OWN_POINTER, row)); 1728be712e4SBarry Smith PetscCall(PetscObjectReference((PetscObject)*row)); 1738be712e4SBarry Smith *col = *row; 1748be712e4SBarry Smith 1758be712e4SBarry Smith PetscCall(PetscFree4(realpart, imagpart, eigvec, work)); 1768be712e4SBarry Smith PetscCall(MatDestroy(&L)); 1778be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 1788be712e4SBarry Smith } 1798be712e4SBarry Smith #endif 1808be712e4SBarry Smith } 181