xref: /petsc/src/mat/graphops/order/sprcm.c (revision 98d129c30f3ee9fdddc40fdbc5a989b7be64f888)
1 #include <petscmat.h>
2 #include <petsc/private/matorderimpl.h>
3 
4 /*
5     MatGetOrdering_RCM - Find the Reverse Cuthill-McKee ordering of a given matrix.
6 */
7 PETSC_INTERN PetscErrorCode MatGetOrdering_RCM(Mat mat, MatOrderingType type, IS *row, IS *col)
8 {
9   PetscInt        i, *mask, *xls, nrow, *perm;
10   const PetscInt *ia, *ja;
11   PetscBool       done;
12 
13   PetscFunctionBegin;
14   PetscCall(MatGetRowIJ(mat, 1, PETSC_TRUE, PETSC_TRUE, &nrow, &ia, &ja, &done));
15   PetscCheck(done, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot get rows for matrix");
16 
17   PetscCall(PetscMalloc3(nrow, &mask, nrow, &perm, 2 * nrow, &xls));
18   PetscCall(SPARSEPACKgenrcm(&nrow, ia, ja, perm, mask, xls));
19   PetscCall(MatRestoreRowIJ(mat, 1, PETSC_TRUE, PETSC_TRUE, NULL, &ia, &ja, &done));
20 
21   /* shift because Sparsepack indices start at one */
22   for (i = 0; i < nrow; i++) perm[i]--;
23   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrow, perm, PETSC_COPY_VALUES, row));
24   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrow, perm, PETSC_COPY_VALUES, col));
25   PetscCall(PetscFree3(mask, perm, xls));
26   PetscFunctionReturn(PETSC_SUCCESS);
27 }
28