Actual source code: spqmd.c

  1: #define PETSCMAT_DLL

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

  7: /*
  8:     MatOrdering_QMD - Find the Quotient Minimum Degree ordering of a given matrix.
  9: */
 12: PetscErrorCode  MatOrdering_QMD(Mat mat,const MatOrderingType type,IS *row,IS *col)
 13: {
 14:   PetscInt       i,  *deg,*marker,*rchset,*nbrhd,*qsize,*qlink,nofsub,*iperm,nrow;
 16:   PetscInt       *ia,*ja,*perm;
 17:   PetscTruth      done;

 20:   MatGetRowIJ(mat,1,PETSC_TRUE,&nrow,&ia,&ja,&done);
 21:   if (!done) SETERRQ(PETSC_ERR_SUP,"Cannot get rows for matrix");

 23:   PetscMalloc(nrow * sizeof(PetscInt),&perm);
 24:   PetscMalloc(nrow * sizeof(PetscInt),&iperm);
 25:   PetscMalloc(nrow * sizeof(PetscInt),&deg);
 26:   PetscMalloc(nrow * sizeof(PetscInt),&marker);
 27:   PetscMalloc(nrow * sizeof(PetscInt),&rchset);
 28:   PetscMalloc(nrow * sizeof(PetscInt),&nbrhd);
 29:   PetscMalloc(nrow * sizeof(PetscInt),&qsize);
 30:   PetscMalloc(nrow * sizeof(PetscInt),&qlink);
 31:   /* WARNING - genqmd trashes ja */
 32:   SPARSEPACKgenqmd(&nrow,ia,ja,perm,iperm,deg,marker,rchset,nbrhd,qsize,qlink,&nofsub);
 33:   MatRestoreRowIJ(mat,1,PETSC_TRUE,&nrow,&ia,&ja,&done);

 35:   PetscFree(deg);
 36:   PetscFree(marker);
 37:   PetscFree(rchset);
 38:   PetscFree(nbrhd);
 39:   PetscFree(qsize);
 40:   PetscFree(qlink);
 41:   PetscFree(iperm);
 42:   for (i=0; i<nrow; i++) perm[i]--;
 43:   ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,row);
 44:   ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,col);
 45:   PetscFree(perm);
 46:   return(0);
 47: }