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