Actual source code: spqmd.c
2: #include petscmat.h
3: #include src/mat/order/order.h
6: /*
7: MatOrdering_QMD - Find the Quotient Minimum Degree ordering of a given matrix.
8: */
11: PetscErrorCode MatOrdering_QMD(Mat mat,const MatOrderingType type,IS *row,IS *col)
12: {
13: PetscInt i, *deg,*marker,*rchset,*nbrhd,*qsize,*qlink,nofsub,*iperm,nrow;
15: PetscInt *ia,*ja,*perm;
16: PetscTruth done;
19: MatGetRowIJ(mat,1,PETSC_TRUE,&nrow,&ia,&ja,&done);
20: if (!done) SETERRQ(PETSC_ERR_SUP,"Cannot get rows for matrix");
22: PetscMalloc(nrow * sizeof(PetscInt),&perm);
23: PetscMalloc(nrow * sizeof(PetscInt),&iperm);
24: PetscMalloc(nrow * sizeof(PetscInt),°);
25: PetscMalloc(nrow * sizeof(PetscInt),&marker);
26: PetscMalloc(nrow * sizeof(PetscInt),&rchset);
27: PetscMalloc(nrow * sizeof(PetscInt),&nbrhd);
28: PetscMalloc(nrow * sizeof(PetscInt),&qsize);
29: PetscMalloc(nrow * sizeof(PetscInt),&qlink);
30: /* WARNING - genqmd trashes ja */
31: SPARSEPACKgenqmd(&nrow,ia,ja,perm,iperm,deg,marker,rchset,nbrhd,qsize,qlink,&nofsub);
32: MatRestoreRowIJ(mat,1,PETSC_TRUE,&nrow,&ia,&ja,&done);
34: PetscFree(deg);
35: PetscFree(marker);
36: PetscFree(rchset);
37: PetscFree(nbrhd);
38: PetscFree(qsize);
39: PetscFree(qlink);
40: PetscFree(iperm);
41: for (i=0; i<nrow; i++) perm[i]--;
42: ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,row);
43: ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,col);
44: PetscFree(perm);
45: return(0);
46: }