Actual source code: spqmd.c

  1: #include <petscmat.h>
  2: #include <petsc/private/matorderimpl.h>

  4: /*
  5:     MatGetOrdering_QMD - Find the Quotient Minimum Degree ordering of a given matrix.
  6: */
  7: PETSC_INTERN PetscErrorCode MatGetOrdering_QMD(Mat mat, MatOrderingType type, IS *row, IS *col)
  8: {
  9:   PetscInt        i, *deg, *marker, *rchset, *nbrhd, *qsize, *qlink, nofsub, *iperm, nrow, *perm;
 10:   const PetscInt *ia, *ja;
 11:   PetscBool       done;

 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");

 17:   PetscCall(PetscMalloc1(nrow, &perm));
 18:   PetscCall(PetscMalloc5(nrow, &iperm, nrow, &deg, nrow, &marker, nrow, &rchset, nrow, &nbrhd));
 19:   PetscCall(PetscMalloc2(nrow, &qsize, nrow, &qlink));
 20:   /* WARNING - genqmd trashes ja */
 21:   PetscCall(SPARSEPACKgenqmd(&nrow, ia, ja, perm, iperm, deg, marker, rchset, nbrhd, qsize, qlink, &nofsub));
 22:   PetscCall(MatRestoreRowIJ(mat, 1, PETSC_TRUE, PETSC_TRUE, NULL, &ia, &ja, &done));

 24:   PetscCall(PetscFree2(qsize, qlink));
 25:   PetscCall(PetscFree5(iperm, deg, marker, rchset, nbrhd));
 26:   for (i = 0; i < nrow; i++) perm[i]--;
 27:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrow, perm, PETSC_COPY_VALUES, row));
 28:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrow, perm, PETSC_OWN_POINTER, col));
 29:   PetscFunctionReturn(PETSC_SUCCESS);
 30: }