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