xref: /petsc/src/mat/graphops/order/spectral.c (revision 8be712e46db5d855f641c6bd97b4543e0efe65bd)
1*8be712e4SBarry Smith #include <petscmat.h> /*I <petscmat.h> I*/
2*8be712e4SBarry Smith #include <petscblaslapack.h>
3*8be712e4SBarry Smith 
4*8be712e4SBarry Smith /*@
5*8be712e4SBarry Smith   MatCreateLaplacian - Create the matrix Laplacian, with all values in the matrix less than the tolerance set to zero
6*8be712e4SBarry Smith 
7*8be712e4SBarry Smith   Input Parameters:
8*8be712e4SBarry Smith + A        - The matrix
9*8be712e4SBarry Smith . tol      - The zero tolerance
10*8be712e4SBarry Smith - weighted - Flag for using edge weights
11*8be712e4SBarry Smith 
12*8be712e4SBarry Smith   Output Parameter:
13*8be712e4SBarry Smith . L - The graph Laplacian matrix
14*8be712e4SBarry Smith 
15*8be712e4SBarry Smith   Level: intermediate
16*8be712e4SBarry Smith 
17*8be712e4SBarry Smith .seealso: `MatFilter()`, `MatGetGraph()`
18*8be712e4SBarry Smith  @*/
19*8be712e4SBarry Smith PetscErrorCode MatCreateLaplacian(Mat A, PetscReal tol, PetscBool weighted, Mat *L)
20*8be712e4SBarry Smith {
21*8be712e4SBarry Smith   PetscScalar *newVals;
22*8be712e4SBarry Smith   PetscInt    *newCols;
23*8be712e4SBarry Smith   PetscInt     rStart, rEnd, r, colMax = 0;
24*8be712e4SBarry Smith   PetscInt    *dnnz, *onnz;
25*8be712e4SBarry Smith   PetscInt     m, n, M, N;
26*8be712e4SBarry Smith 
27*8be712e4SBarry Smith   PetscFunctionBegin;
28*8be712e4SBarry Smith   PetscCheck(!weighted, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Will get to this soon");
29*8be712e4SBarry Smith   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), L));
30*8be712e4SBarry Smith   PetscCall(MatGetSize(A, &M, &N));
31*8be712e4SBarry Smith   PetscCall(MatGetLocalSize(A, &m, &n));
32*8be712e4SBarry Smith   PetscCall(MatSetSizes(*L, m, n, M, N));
33*8be712e4SBarry Smith   PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd));
34*8be712e4SBarry Smith   PetscCall(PetscMalloc2(m, &dnnz, m, &onnz));
35*8be712e4SBarry Smith   for (r = rStart; r < rEnd; ++r) {
36*8be712e4SBarry Smith     const PetscScalar *vals;
37*8be712e4SBarry Smith     const PetscInt    *cols;
38*8be712e4SBarry Smith     PetscInt           ncols, newcols, c;
39*8be712e4SBarry Smith     PetscBool          hasdiag = PETSC_FALSE;
40*8be712e4SBarry Smith 
41*8be712e4SBarry Smith     dnnz[r - rStart] = onnz[r - rStart] = 0;
42*8be712e4SBarry Smith     PetscCall(MatGetRow(A, r, &ncols, &cols, &vals));
43*8be712e4SBarry Smith     for (c = 0, newcols = 0; c < ncols; ++c) {
44*8be712e4SBarry Smith       if (cols[c] == r) {
45*8be712e4SBarry Smith         ++newcols;
46*8be712e4SBarry Smith         hasdiag = PETSC_TRUE;
47*8be712e4SBarry Smith         ++dnnz[r - rStart];
48*8be712e4SBarry Smith       } else if (PetscAbsScalar(vals[c]) >= tol) {
49*8be712e4SBarry Smith         if ((cols[c] >= rStart) && (cols[c] < rEnd)) ++dnnz[r - rStart];
50*8be712e4SBarry Smith         else ++onnz[r - rStart];
51*8be712e4SBarry Smith         ++newcols;
52*8be712e4SBarry Smith       }
53*8be712e4SBarry Smith     }
54*8be712e4SBarry Smith     if (!hasdiag) {
55*8be712e4SBarry Smith       ++newcols;
56*8be712e4SBarry Smith       ++dnnz[r - rStart];
57*8be712e4SBarry Smith     }
58*8be712e4SBarry Smith     colMax = PetscMax(colMax, newcols);
59*8be712e4SBarry Smith     PetscCall(MatRestoreRow(A, r, &ncols, &cols, &vals));
60*8be712e4SBarry Smith   }
61*8be712e4SBarry Smith   PetscCall(MatSetFromOptions(*L));
62*8be712e4SBarry Smith   PetscCall(MatXAIJSetPreallocation(*L, 1, dnnz, onnz, NULL, NULL));
63*8be712e4SBarry Smith   PetscCall(MatSetUp(*L));
64*8be712e4SBarry Smith   PetscCall(PetscMalloc2(colMax, &newCols, colMax, &newVals));
65*8be712e4SBarry Smith   for (r = rStart; r < rEnd; ++r) {
66*8be712e4SBarry Smith     const PetscScalar *vals;
67*8be712e4SBarry Smith     const PetscInt    *cols;
68*8be712e4SBarry Smith     PetscInt           ncols, newcols, c;
69*8be712e4SBarry Smith     PetscBool          hasdiag = PETSC_FALSE;
70*8be712e4SBarry Smith 
71*8be712e4SBarry Smith     PetscCall(MatGetRow(A, r, &ncols, &cols, &vals));
72*8be712e4SBarry Smith     for (c = 0, newcols = 0; c < ncols; ++c) {
73*8be712e4SBarry Smith       if (cols[c] == r) {
74*8be712e4SBarry Smith         newCols[newcols] = cols[c];
75*8be712e4SBarry Smith         newVals[newcols] = dnnz[r - rStart] + onnz[r - rStart] - 1;
76*8be712e4SBarry Smith         ++newcols;
77*8be712e4SBarry Smith         hasdiag = PETSC_TRUE;
78*8be712e4SBarry Smith       } else if (PetscAbsScalar(vals[c]) >= tol) {
79*8be712e4SBarry Smith         newCols[newcols] = cols[c];
80*8be712e4SBarry Smith         newVals[newcols] = -1.0;
81*8be712e4SBarry Smith         ++newcols;
82*8be712e4SBarry Smith       }
83*8be712e4SBarry Smith       PetscCheck(newcols <= colMax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Overran work space");
84*8be712e4SBarry Smith     }
85*8be712e4SBarry Smith     if (!hasdiag) {
86*8be712e4SBarry Smith       newCols[newcols] = r;
87*8be712e4SBarry Smith       newVals[newcols] = dnnz[r - rStart] + onnz[r - rStart] - 1;
88*8be712e4SBarry Smith       ++newcols;
89*8be712e4SBarry Smith     }
90*8be712e4SBarry Smith     PetscCall(MatRestoreRow(A, r, &ncols, &cols, &vals));
91*8be712e4SBarry Smith     PetscCall(MatSetValues(*L, 1, &r, newcols, newCols, newVals, INSERT_VALUES));
92*8be712e4SBarry Smith   }
93*8be712e4SBarry Smith   PetscCall(PetscFree2(dnnz, onnz));
94*8be712e4SBarry Smith   PetscCall(MatAssemblyBegin(*L, MAT_FINAL_ASSEMBLY));
95*8be712e4SBarry Smith   PetscCall(MatAssemblyEnd(*L, MAT_FINAL_ASSEMBLY));
96*8be712e4SBarry Smith   PetscCall(PetscFree2(newCols, newVals));
97*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
98*8be712e4SBarry Smith }
99*8be712e4SBarry Smith 
100*8be712e4SBarry Smith /*
101*8be712e4SBarry Smith   MatGetOrdering_Spectral - Find the symmetric reordering of the graph by .
102*8be712e4SBarry Smith */
103*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_Spectral(Mat A, MatOrderingType type, IS *row, IS *col)
104*8be712e4SBarry Smith {
105*8be712e4SBarry Smith   Mat             L;
106*8be712e4SBarry Smith   const PetscReal eps = 1.0e-12;
107*8be712e4SBarry Smith 
108*8be712e4SBarry Smith   PetscFunctionBegin;
109*8be712e4SBarry Smith   PetscCall(MatCreateLaplacian(A, eps, PETSC_FALSE, &L));
110*8be712e4SBarry Smith   {
111*8be712e4SBarry Smith     /* Check Laplacian */
112*8be712e4SBarry Smith     PetscReal norm;
113*8be712e4SBarry Smith     Vec       x, y;
114*8be712e4SBarry Smith 
115*8be712e4SBarry Smith     PetscCall(MatCreateVecs(L, &x, NULL));
116*8be712e4SBarry Smith     PetscCall(VecDuplicate(x, &y));
117*8be712e4SBarry Smith     PetscCall(VecSet(x, 1.0));
118*8be712e4SBarry Smith     PetscCall(MatMult(L, x, y));
119*8be712e4SBarry Smith     PetscCall(VecNorm(y, NORM_INFINITY, &norm));
120*8be712e4SBarry Smith     PetscCheck(norm <= 1.0e-10, PetscObjectComm((PetscObject)y), PETSC_ERR_PLIB, "Invalid graph Laplacian");
121*8be712e4SBarry Smith     PetscCall(VecDestroy(&x));
122*8be712e4SBarry Smith     PetscCall(VecDestroy(&y));
123*8be712e4SBarry Smith   }
124*8be712e4SBarry Smith   /* Compute Fiedler vector (right now, all eigenvectors) */
125*8be712e4SBarry Smith #if defined(PETSC_USE_COMPLEX)
126*8be712e4SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Spectral partitioning does not support complex numbers");
127*8be712e4SBarry Smith #else
128*8be712e4SBarry Smith   {
129*8be712e4SBarry Smith     Mat          LD;
130*8be712e4SBarry Smith     PetscScalar *a;
131*8be712e4SBarry Smith     PetscReal   *realpart, *imagpart, *eigvec, *work;
132*8be712e4SBarry Smith     PetscReal    sdummy;
133*8be712e4SBarry Smith     PetscBLASInt bn, bN, lwork = 0, lierr, idummy;
134*8be712e4SBarry Smith     PetscInt     n, i, evInd, *perm, tmp;
135*8be712e4SBarry Smith 
136*8be712e4SBarry Smith     PetscCall(MatConvert(L, MATDENSE, MAT_INITIAL_MATRIX, &LD));
137*8be712e4SBarry Smith     PetscCall(MatGetLocalSize(LD, &n, NULL));
138*8be712e4SBarry Smith     PetscCall(MatDenseGetArray(LD, &a));
139*8be712e4SBarry Smith     PetscCall(PetscBLASIntCast(n, &bn));
140*8be712e4SBarry Smith     PetscCall(PetscBLASIntCast(n, &bN));
141*8be712e4SBarry Smith     PetscCall(PetscBLASIntCast(5 * n, &lwork));
142*8be712e4SBarry Smith     PetscCall(PetscBLASIntCast(1, &idummy));
143*8be712e4SBarry Smith     PetscCall(PetscMalloc4(n, &realpart, n, &imagpart, n * n, &eigvec, lwork, &work));
144*8be712e4SBarry Smith     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
145*8be712e4SBarry Smith     PetscCallBLAS("LAPACKgeev", LAPACKgeev_("N", "V", &bn, a, &bN, realpart, imagpart, &sdummy, &idummy, eigvec, &bN, work, &lwork, &lierr));
146*8be712e4SBarry Smith     PetscCheck(!lierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int)lierr);
147*8be712e4SBarry Smith     PetscCall(PetscFPTrapPop());
148*8be712e4SBarry Smith     PetscCall(MatDenseRestoreArray(LD, &a));
149*8be712e4SBarry Smith     PetscCall(MatDestroy(&LD));
150*8be712e4SBarry Smith     /* Check lowest eigenvalue and eigenvector */
151*8be712e4SBarry Smith     PetscCall(PetscMalloc1(n, &perm));
152*8be712e4SBarry Smith     for (i = 0; i < n; ++i) perm[i] = i;
153*8be712e4SBarry Smith     PetscCall(PetscSortRealWithPermutation(n, realpart, perm));
154*8be712e4SBarry Smith     evInd = perm[0];
155*8be712e4SBarry 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");
156*8be712e4SBarry Smith     evInd = perm[1];
157*8be712e4SBarry 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");
158*8be712e4SBarry Smith     evInd = perm[0];
159*8be712e4SBarry Smith     for (i = 0; i < n; ++i) {
160*8be712e4SBarry Smith       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]));
161*8be712e4SBarry Smith     }
162*8be712e4SBarry Smith     /* Construct Fiedler partition */
163*8be712e4SBarry Smith     evInd = perm[1];
164*8be712e4SBarry Smith     for (i = 0; i < n; ++i) perm[i] = i;
165*8be712e4SBarry Smith     PetscCall(PetscSortRealWithPermutation(n, &eigvec[evInd * n], perm));
166*8be712e4SBarry Smith     for (i = 0; i < n / 2; ++i) {
167*8be712e4SBarry Smith       tmp             = perm[n - 1 - i];
168*8be712e4SBarry Smith       perm[n - 1 - i] = perm[i];
169*8be712e4SBarry Smith       perm[i]         = tmp;
170*8be712e4SBarry Smith     }
171*8be712e4SBarry Smith     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, perm, PETSC_OWN_POINTER, row));
172*8be712e4SBarry Smith     PetscCall(PetscObjectReference((PetscObject)*row));
173*8be712e4SBarry Smith     *col = *row;
174*8be712e4SBarry Smith 
175*8be712e4SBarry Smith     PetscCall(PetscFree4(realpart, imagpart, eigvec, work));
176*8be712e4SBarry Smith     PetscCall(MatDestroy(&L));
177*8be712e4SBarry Smith     PetscFunctionReturn(PETSC_SUCCESS);
178*8be712e4SBarry Smith   }
179*8be712e4SBarry Smith #endif
180*8be712e4SBarry Smith }
181