Actual source code: spqmd.c

  1: /*$Id: spqmd.c,v 1.40 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_QMD - Find the Quotient Minimum Degree ordering of a given matrix.
  9: */
 10: int MatOrdering_QMD(Mat mat,MatOrderingType type,IS *row,IS *col)
 11: {
 12:   int        i,  *deg,*marker,*rchset,*nbrhd,*qsize,*qlink,nofsub,*iperm,nrow;
 13:   int        ierr,*ia,*ja,*perm;
 14:   PetscTruth done;

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

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

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