xref: /petsc/src/mat/graphops/order/spqmd.c (revision 53673ba54f5aaba04b9d49ab22cf56c7a7461fe9)
1*8be712e4SBarry Smith #include <petscmat.h>
2*8be712e4SBarry Smith #include <petsc/private/matorderimpl.h>
3*8be712e4SBarry Smith 
4*8be712e4SBarry Smith /*
5*8be712e4SBarry Smith     MatGetOrdering_QMD - Find the Quotient Minimum Degree ordering of a given matrix.
6*8be712e4SBarry Smith */
MatGetOrdering_QMD(Mat mat,MatOrderingType type,IS * row,IS * col)7*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_QMD(Mat mat, MatOrderingType type, IS *row, IS *col)
8*8be712e4SBarry Smith {
9*8be712e4SBarry Smith   PetscInt        i, *deg, *marker, *rchset, *nbrhd, *qsize, *qlink, nofsub, *iperm, nrow, *perm;
10*8be712e4SBarry Smith   const PetscInt *ia, *ja;
11*8be712e4SBarry Smith   PetscBool       done;
12*8be712e4SBarry Smith 
13*8be712e4SBarry Smith   PetscFunctionBegin;
14*8be712e4SBarry Smith   PetscCall(MatGetRowIJ(mat, 1, PETSC_TRUE, PETSC_TRUE, &nrow, &ia, &ja, &done));
15*8be712e4SBarry Smith   PetscCheck(done, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot get rows for matrix");
16*8be712e4SBarry Smith 
17*8be712e4SBarry Smith   PetscCall(PetscMalloc1(nrow, &perm));
18*8be712e4SBarry Smith   PetscCall(PetscMalloc5(nrow, &iperm, nrow, &deg, nrow, &marker, nrow, &rchset, nrow, &nbrhd));
19*8be712e4SBarry Smith   PetscCall(PetscMalloc2(nrow, &qsize, nrow, &qlink));
20*8be712e4SBarry Smith   /* WARNING - genqmd trashes ja */
21*8be712e4SBarry Smith   PetscCall(SPARSEPACKgenqmd(&nrow, ia, ja, perm, iperm, deg, marker, rchset, nbrhd, qsize, qlink, &nofsub));
22*8be712e4SBarry Smith   PetscCall(MatRestoreRowIJ(mat, 1, PETSC_TRUE, PETSC_TRUE, NULL, &ia, &ja, &done));
23*8be712e4SBarry Smith 
24*8be712e4SBarry Smith   PetscCall(PetscFree2(qsize, qlink));
25*8be712e4SBarry Smith   PetscCall(PetscFree5(iperm, deg, marker, rchset, nbrhd));
26*8be712e4SBarry Smith   for (i = 0; i < nrow; i++) perm[i]--;
27*8be712e4SBarry Smith   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrow, perm, PETSC_COPY_VALUES, row));
28*8be712e4SBarry Smith   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrow, perm, PETSC_OWN_POINTER, col));
29*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
30*8be712e4SBarry Smith }
31