Actual source code: spnd.c


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

  5: /*
  6:     MatGetOrdering_ND - Find the nested dissection ordering of a given matrix.
  7: */
  8: PETSC_INTERN PetscErrorCode MatGetOrdering_ND(Mat mat, MatOrderingType type, IS *row, IS *col)
  9: {
 10:   PetscInt        i, *mask, *xls, *ls, nrow, *perm;
 11:   const PetscInt *ia, *ja;
 12:   PetscBool       done;
 13:   Mat             B = NULL;

 15:   PetscFunctionBegin;
 16:   PetscCall(MatGetRowIJ(mat, 1, PETSC_TRUE, PETSC_TRUE, &nrow, &ia, &ja, &done));
 17:   if (!done) {
 18:     PetscCall(MatConvert(mat, MATAIJ, MAT_INITIAL_MATRIX, &B));
 19:     PetscCall(MatGetRowIJ(B, 1, PETSC_TRUE, PETSC_TRUE, &nrow, &ia, &ja, &done));
 20:   }

 22:   PetscCall(PetscMalloc4(nrow, &mask, nrow, &perm, nrow + 1, &xls, nrow, &ls));
 23:   PetscCall(SPARSEPACKgennd(&nrow, ia, ja, mask, perm, xls, ls));
 24:   if (B) {
 25:     PetscCall(MatRestoreRowIJ(B, 1, PETSC_TRUE, PETSC_TRUE, NULL, &ia, &ja, &done));
 26:     PetscCall(MatDestroy(&B));
 27:   } else {
 28:     PetscCall(MatRestoreRowIJ(mat, 1, PETSC_TRUE, PETSC_TRUE, NULL, &ia, &ja, &done));
 29:   }

 31:   /* shift because Sparsepack indices start at one */
 32:   for (i = 0; i < nrow; i++) perm[i]--;

 34:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrow, perm, PETSC_COPY_VALUES, row));
 35:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrow, perm, PETSC_COPY_VALUES, col));
 36:   PetscCall(PetscFree4(mask, perm, xls, ls));
 37:   PetscFunctionReturn(PETSC_SUCCESS);
 38: }