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: MatGetRowIJ(mat,1,PETSC_TRUE,PETSC_TRUE,&nrow,&ia,&ja,&done);
16: if (!done) {
17: MatConvert(mat,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
18: MatGetRowIJ(B,1,PETSC_TRUE,PETSC_TRUE,&nrow,&ia,&ja,&done);
19: }
21: PetscMalloc4(nrow,&mask,nrow,&perm,nrow+1,&xls,nrow,&ls);
22: SPARSEPACKgennd(&nrow,ia,ja,mask,perm,xls,ls);
23: if (B) {
24: MatRestoreRowIJ(B,1,PETSC_TRUE,PETSC_TRUE,NULL,&ia,&ja,&done);
25: MatDestroy(&B);
26: } else {
27: MatRestoreRowIJ(mat,1,PETSC_TRUE,PETSC_TRUE,NULL,&ia,&ja,&done);
28: }
30: /* shift because Sparsepack indices start at one */
31: for (i=0; i<nrow; i++) perm[i]--;
33: ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,PETSC_COPY_VALUES,row);
34: ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,PETSC_COPY_VALUES,col);
35: PetscFree4(mask,perm,xls,ls);
36: return 0;
37: }