1 2 static char help[] = "Tests matrix ordering routines.\n\n"; 3 4 #include <petscmat.h> 5 extern PetscErrorCode MatGetOrdering_myordering(Mat,MatOrderingType,IS*,IS*); 6 7 int main(int argc,char **args) 8 { 9 Mat C,Cperm; 10 PetscInt i,j,m = 5,n = 5,Ii,J,ncols; 11 PetscScalar v; 12 PetscMPIInt size; 13 IS rperm,cperm,icperm; 14 const PetscInt *rperm_ptr,*cperm_ptr,*cols; 15 const PetscScalar *vals; 16 PetscBool TestMyorder=PETSC_FALSE; 17 18 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 19 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 20 PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 21 22 /* create the matrix for the five point stencil, YET AGAIN */ 23 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C)); 24 PetscCall(MatSetUp(C)); 25 for (i=0; i<m; i++) { 26 for (j=0; j<n; j++) { 27 v = -1.0; Ii = j + n*i; 28 if (i>0) {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 29 if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 30 if (j>0) {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 31 if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 32 v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 33 } 34 } 35 PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 36 PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 37 38 PetscCall(MatGetOrdering(C,MATORDERINGND,&rperm,&cperm)); 39 PetscCall(ISView(rperm,PETSC_VIEWER_STDOUT_SELF)); 40 PetscCall(ISDestroy(&rperm)); 41 PetscCall(ISDestroy(&cperm)); 42 43 PetscCall(MatGetOrdering(C,MATORDERINGRCM,&rperm,&cperm)); 44 PetscCall(ISView(rperm,PETSC_VIEWER_STDOUT_SELF)); 45 PetscCall(ISDestroy(&rperm)); 46 PetscCall(ISDestroy(&cperm)); 47 48 PetscCall(MatGetOrdering(C,MATORDERINGQMD,&rperm,&cperm)); 49 PetscCall(ISView(rperm,PETSC_VIEWER_STDOUT_SELF)); 50 PetscCall(ISDestroy(&rperm)); 51 PetscCall(ISDestroy(&cperm)); 52 53 /* create Cperm = rperm*C*icperm */ 54 PetscCall(PetscOptionsGetBool(NULL,NULL,"-testmyordering",&TestMyorder,NULL)); 55 if (TestMyorder) { 56 PetscCall(MatGetOrdering_myordering(C,MATORDERINGQMD,&rperm,&cperm)); 57 printf("myordering's rperm:\n"); 58 PetscCall(ISView(rperm,PETSC_VIEWER_STDOUT_SELF)); 59 PetscCall(ISInvertPermutation(cperm,PETSC_DECIDE,&icperm)); 60 PetscCall(ISGetIndices(rperm,&rperm_ptr)); 61 PetscCall(ISGetIndices(icperm,&cperm_ptr)); 62 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&Cperm)); 63 for (i=0; i<m*n; i++) { 64 PetscCall(MatGetRow(C,rperm_ptr[i],&ncols,&cols,&vals)); 65 for (j=0; j<ncols; j++) { 66 /* printf(" (%d %d %g)\n",i,cperm_ptr[cols[j]],vals[j]); */ 67 PetscCall(MatSetValues(Cperm,1,&i,1,&cperm_ptr[cols[j]],&vals[j],INSERT_VALUES)); 68 } 69 } 70 PetscCall(MatAssemblyBegin(Cperm,MAT_FINAL_ASSEMBLY)); 71 PetscCall(MatAssemblyEnd(Cperm,MAT_FINAL_ASSEMBLY)); 72 PetscCall(ISRestoreIndices(rperm,&rperm_ptr)); 73 PetscCall(ISRestoreIndices(icperm,&cperm_ptr)); 74 75 PetscCall(ISDestroy(&rperm)); 76 PetscCall(ISDestroy(&cperm)); 77 PetscCall(ISDestroy(&icperm)); 78 PetscCall(MatDestroy(&Cperm)); 79 } 80 81 PetscCall(MatDestroy(&C)); 82 PetscCall(PetscFinalize()); 83 return 0; 84 } 85 86 #include <petsc/private/matimpl.h> 87 /* This is modified from MatGetOrdering_Natural() */ 88 PetscErrorCode MatGetOrdering_myordering(Mat mat,MatOrderingType type,IS *irow,IS *icol) 89 { 90 PetscInt n,i,*ii; 91 PetscBool done; 92 MPI_Comm comm; 93 94 PetscFunctionBegin; 95 PetscCall(PetscObjectGetComm((PetscObject)mat,&comm)); 96 PetscCall(MatGetRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,&n,NULL,NULL,&done)); 97 PetscCall(MatRestoreRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,NULL,NULL,NULL,&done)); 98 if (done) { /* matrix may be "compressed" in symbolic factorization, due to i-nodes or block storage */ 99 PetscCall(PetscMalloc1(n,&ii)); 100 for (i=0; i<n; i++) ii[i] = n-i-1; /* replace your index here */ 101 PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_COPY_VALUES,irow)); 102 PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_OWN_POINTER,icol)); 103 } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"MatRestoreRowIJ fails!"); 104 PetscCall(ISSetPermutation(*irow)); 105 PetscCall(ISSetPermutation(*icol)); 106 PetscFunctionReturn(0); 107 } 108 109 /*TEST 110 111 test: 112 113 TEST*/ 114