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 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 21 PetscCheckFalse(size != 1,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 CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C)); 25 CHKERRQ(MatSetUp(C)); 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; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 30 if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 31 if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 32 if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 33 v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 34 } 35 } 36 CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 37 CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 38 39 CHKERRQ(MatGetOrdering(C,MATORDERINGND,&rperm,&cperm)); 40 CHKERRQ(ISView(rperm,PETSC_VIEWER_STDOUT_SELF)); 41 CHKERRQ(ISDestroy(&rperm)); 42 CHKERRQ(ISDestroy(&cperm)); 43 44 CHKERRQ(MatGetOrdering(C,MATORDERINGRCM,&rperm,&cperm)); 45 CHKERRQ(ISView(rperm,PETSC_VIEWER_STDOUT_SELF)); 46 CHKERRQ(ISDestroy(&rperm)); 47 CHKERRQ(ISDestroy(&cperm)); 48 49 CHKERRQ(MatGetOrdering(C,MATORDERINGQMD,&rperm,&cperm)); 50 CHKERRQ(ISView(rperm,PETSC_VIEWER_STDOUT_SELF)); 51 CHKERRQ(ISDestroy(&rperm)); 52 CHKERRQ(ISDestroy(&cperm)); 53 54 /* create Cperm = rperm*C*icperm */ 55 CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-testmyordering",&TestMyorder,NULL)); 56 if (TestMyorder) { 57 CHKERRQ(MatGetOrdering_myordering(C,MATORDERINGQMD,&rperm,&cperm)); 58 printf("myordering's rperm:\n"); 59 CHKERRQ(ISView(rperm,PETSC_VIEWER_STDOUT_SELF)); 60 CHKERRQ(ISInvertPermutation(cperm,PETSC_DECIDE,&icperm)); 61 CHKERRQ(ISGetIndices(rperm,&rperm_ptr)); 62 CHKERRQ(ISGetIndices(icperm,&cperm_ptr)); 63 CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&Cperm)); 64 for (i=0; i<m*n; i++) { 65 CHKERRQ(MatGetRow(C,rperm_ptr[i],&ncols,&cols,&vals)); 66 for (j=0; j<ncols; j++) { 67 /* printf(" (%d %d %g)\n",i,cperm_ptr[cols[j]],vals[j]); */ 68 CHKERRQ(MatSetValues(Cperm,1,&i,1,&cperm_ptr[cols[j]],&vals[j],INSERT_VALUES)); 69 } 70 } 71 CHKERRQ(MatAssemblyBegin(Cperm,MAT_FINAL_ASSEMBLY)); 72 CHKERRQ(MatAssemblyEnd(Cperm,MAT_FINAL_ASSEMBLY)); 73 CHKERRQ(ISRestoreIndices(rperm,&rperm_ptr)); 74 CHKERRQ(ISRestoreIndices(icperm,&cperm_ptr)); 75 76 CHKERRQ(ISDestroy(&rperm)); 77 CHKERRQ(ISDestroy(&cperm)); 78 CHKERRQ(ISDestroy(&icperm)); 79 CHKERRQ(MatDestroy(&Cperm)); 80 } 81 82 CHKERRQ(MatDestroy(&C)); 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 PetscInt n,i,*ii; 92 PetscBool done; 93 MPI_Comm comm; 94 95 PetscFunctionBegin; 96 CHKERRQ(PetscObjectGetComm((PetscObject)mat,&comm)); 97 CHKERRQ(MatGetRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,&n,NULL,NULL,&done)); 98 CHKERRQ(MatRestoreRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,NULL,NULL,NULL,&done)); 99 if (done) { /* matrix may be "compressed" in symbolic factorization, due to i-nodes or block storage */ 100 CHKERRQ(PetscMalloc1(n,&ii)); 101 for (i=0; i<n; i++) ii[i] = n-i-1; /* replace your index here */ 102 CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_COPY_VALUES,irow)); 103 CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_OWN_POINTER,icol)); 104 } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"MatRestoreRowIJ fails!"); 105 CHKERRQ(ISSetPermutation(*irow)); 106 CHKERRQ(ISSetPermutation(*icol)); 107 PetscFunctionReturn(0); 108 } 109 110 /*TEST 111 112 test: 113 114 TEST*/ 115