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: {
 11:   PetscInt       i, *mask,*xls,*ls,nrow,*perm;
 12:   const PetscInt *ia,*ja;
 13:   PetscBool      done;
 14:   Mat            B = NULL;

 17:   MatGetRowIJ(mat,1,PETSC_TRUE,PETSC_TRUE,&nrow,&ia,&ja,&done);
 18:   if (!done) {
 19:     MatConvert(mat,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
 20:     MatGetRowIJ(B,1,PETSC_TRUE,PETSC_TRUE,&nrow,&ia,&ja,&done);
 21:   }

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

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

 35:   ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,PETSC_COPY_VALUES,row);
 36:   ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,PETSC_COPY_VALUES,col);
 37:   PetscFree4(mask,perm,xls,ls);
 38:   return(0);
 39: }