1*8be712e4SBarry Smith #include <petscmat.h> 2*8be712e4SBarry Smith #include <petsc/private/matorderimpl.h> 3*8be712e4SBarry Smith 4*8be712e4SBarry Smith /* 5*8be712e4SBarry Smith MatGetOrdering_ND - Find the nested dissection ordering of a given matrix. 6*8be712e4SBarry Smith */ 7*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatGetOrdering_ND(Mat mat, MatOrderingType type, IS *row, IS *col) 8*8be712e4SBarry Smith { 9*8be712e4SBarry Smith PetscInt i, *mask, *xls, *ls, nrow, *perm; 10*8be712e4SBarry Smith const PetscInt *ia, *ja; 11*8be712e4SBarry Smith PetscBool done; 12*8be712e4SBarry Smith Mat B = NULL; 13*8be712e4SBarry Smith 14*8be712e4SBarry Smith PetscFunctionBegin; 15*8be712e4SBarry Smith PetscCall(MatGetRowIJ(mat, 1, PETSC_TRUE, PETSC_TRUE, &nrow, &ia, &ja, &done)); 16*8be712e4SBarry Smith if (!done) { 17*8be712e4SBarry Smith PetscCall(MatConvert(mat, MATAIJ, MAT_INITIAL_MATRIX, &B)); 18*8be712e4SBarry Smith PetscCall(MatGetRowIJ(B, 1, PETSC_TRUE, PETSC_TRUE, &nrow, &ia, &ja, &done)); 19*8be712e4SBarry Smith } 20*8be712e4SBarry Smith 21*8be712e4SBarry Smith PetscCall(PetscMalloc4(nrow, &mask, nrow, &perm, nrow + 1, &xls, nrow, &ls)); 22*8be712e4SBarry Smith PetscCall(SPARSEPACKgennd(&nrow, ia, ja, mask, perm, xls, ls)); 23*8be712e4SBarry Smith if (B) { 24*8be712e4SBarry Smith PetscCall(MatRestoreRowIJ(B, 1, PETSC_TRUE, PETSC_TRUE, NULL, &ia, &ja, &done)); 25*8be712e4SBarry Smith PetscCall(MatDestroy(&B)); 26*8be712e4SBarry Smith } else { 27*8be712e4SBarry Smith PetscCall(MatRestoreRowIJ(mat, 1, PETSC_TRUE, PETSC_TRUE, NULL, &ia, &ja, &done)); 28*8be712e4SBarry Smith } 29*8be712e4SBarry Smith 30*8be712e4SBarry Smith /* shift because Sparsepack indices start at one */ 31*8be712e4SBarry Smith for (i = 0; i < nrow; i++) perm[i]--; 32*8be712e4SBarry Smith 33*8be712e4SBarry Smith PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrow, perm, PETSC_COPY_VALUES, row)); 34*8be712e4SBarry Smith PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrow, perm, PETSC_COPY_VALUES, col)); 35*8be712e4SBarry Smith PetscCall(PetscFree4(mask, perm, xls, ls)); 36*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 37*8be712e4SBarry Smith } 38