Actual source code: spqmd.c
petsc-3.7.3 2016-08-01
2: #include <petscmat.h>
3: #include <petsc/private/matorderimpl.h>
5: /*
6: MatGetOrdering_QMD - Find the Quotient Minimum Degree ordering of a given matrix.
7: */
10: PETSC_INTERN PetscErrorCode MatGetOrdering_QMD(Mat mat,MatOrderingType type,IS *row,IS *col)
11: {
12: PetscInt i, *deg,*marker,*rchset,*nbrhd,*qsize,*qlink,nofsub,*iperm,nrow,*perm;
14: const PetscInt *ia,*ja;
15: PetscBool done;
18: MatGetRowIJ(mat,1,PETSC_TRUE,PETSC_TRUE,&nrow,&ia,&ja,&done);
19: if (!done) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get rows for matrix");
21: PetscMalloc1(nrow,&perm);
22: PetscMalloc5(nrow,&iperm,nrow,°,nrow,&marker,nrow,&rchset,nrow,&nbrhd);
23: PetscMalloc2(nrow,&qsize,nrow,&qlink);
24: /* WARNING - genqmd trashes ja */
25: SPARSEPACKgenqmd(&nrow,ia,ja,perm,iperm,deg,marker,rchset,nbrhd,qsize,qlink,&nofsub);
26: MatRestoreRowIJ(mat,1,PETSC_TRUE,PETSC_TRUE,NULL,&ia,&ja,&done);
28: PetscFree2(qsize,qlink);
29: PetscFree5(iperm,deg,marker,rchset,nbrhd);
30: for (i=0; i<nrow; i++) perm[i]--;
31: ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,PETSC_COPY_VALUES,row);
32: ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,PETSC_OWN_POINTER,col);
33: return(0);
34: }