Actual source code: spnd.c

  1: /*$Id: spnd.c,v 1.41 2001/03/23 23:22:51 balay Exp $*/

 3:  #include petscmat.h
 4:  #include src/mat/order/order.h

  6: EXTERN_C_BEGIN
  7: /*
  8:     MatOrdering_ND - Find the nested dissection ordering of a given matrix.
  9: */
 10: int MatOrdering_ND(Mat mat,MatOrderingType type,IS *row,IS *col)
 11: {
 12:   int        ierr,i, *mask,*xls,*ls,nrow,*ia,*ja,*perm;
 13:   PetscTruth done;

 16:   MatGetRowIJ(mat,1,PETSC_TRUE,&nrow,&ia,&ja,&done);
 17:   if (!done) SETERRQ1(PETSC_ERR_SUP,"Cannot get rows for matrix type %s",((PetscObject)mat)->type_name);

 19:   PetscMalloc((4*nrow +1) * sizeof(int),&mask);
 20:   perm = mask + nrow;
 21:   xls  = perm + nrow;
 22:   ls   = xls  + nrow + 1;

 24:   SPARSEPACKgennd(&nrow,ia,ja,mask,perm,xls,ls);
 25:   MatRestoreRowIJ(mat,1,PETSC_TRUE,&nrow,&ia,&ja,&done);

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

 30:   ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,row);
 31:   ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,col);
 32:   PetscFree(mask);

 34:   return(0);
 35: }
 36: EXTERN_C_END