static char help[] = "Tests matrix ordering routines.\n\n"; #include extern PetscErrorCode MatGetOrdering_myordering(Mat,MatOrderingType,IS*,IS*); int main(int argc,char **args) { Mat C,Cperm; PetscInt i,j,m = 5,n = 5,Ii,J,ncols; PetscScalar v; PetscMPIInt size; IS rperm,cperm,icperm; const PetscInt *rperm_ptr,*cperm_ptr,*cols; const PetscScalar *vals; PetscBool TestMyorder=PETSC_FALSE; PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); /* create the matrix for the five point stencil, YET AGAIN */ PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C)); PetscCall(MatSetUp(C)); for (i=0; i0) {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} if (i0) {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} if (j /* This is modified from MatGetOrdering_Natural() */ PetscErrorCode MatGetOrdering_myordering(Mat mat,MatOrderingType type,IS *irow,IS *icol) { PetscInt n,i,*ii; PetscBool done; MPI_Comm comm; PetscFunctionBegin; PetscCall(PetscObjectGetComm((PetscObject)mat,&comm)); PetscCall(MatGetRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,&n,NULL,NULL,&done)); PetscCall(MatRestoreRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,NULL,NULL,NULL,&done)); if (done) { /* matrix may be "compressed" in symbolic factorization, due to i-nodes or block storage */ PetscCall(PetscMalloc1(n,&ii)); for (i=0; i