Actual source code: ex22.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: static char help[] = "Tests matrix ordering routines.\n\n";

  4: #include <petscmat.h>
  5: extern PetscErrorCode MatGetOrdering_myordering(Mat,MatOrderingType,IS*,IS*);

  9: int main(int argc,char **args)
 10: {
 11:   Mat               C,Cperm;
 12:   PetscInt          i,j,m = 5,n = 5,Ii,J,ncols;
 13:   PetscErrorCode    ierr;
 14:   PetscScalar       v;
 15:   PetscMPIInt       size;
 16:   IS                rperm,cperm,icperm;
 17:   const PetscInt    *rperm_ptr,*cperm_ptr,*cols;
 18:   const PetscScalar *vals;
 19:   PetscBool         TestMyorder=PETSC_FALSE;

 21:   PetscInitialize(&argc,&args,(char*)0,help);
 22:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 23:   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");

 25:   /* create the matrix for the five point stencil, YET AGAIN */
 26:   MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C);
 27:   MatSetUp(C);
 28:   for (i=0; i<m; i++) {
 29:     for (j=0; j<n; j++) {
 30:       v = -1.0;  Ii = j + n*i;
 31:       if (i>0)   {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 32:       if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 33:       if (j>0)   {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 34:       if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 35:       v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 36:     }
 37:   }
 38:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 39:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 41:   MatGetOrdering(C,MATORDERINGND,&rperm,&cperm);
 42:   ISView(rperm,PETSC_VIEWER_STDOUT_SELF);
 43:   ISDestroy(&rperm);
 44:   ISDestroy(&cperm);

 46:   MatGetOrdering(C,MATORDERINGRCM,&rperm,&cperm);
 47:   ISView(rperm,PETSC_VIEWER_STDOUT_SELF);
 48:   ISDestroy(&rperm);
 49:   ISDestroy(&cperm);

 51:   MatGetOrdering(C,MATORDERINGQMD,&rperm,&cperm);
 52:   ISView(rperm,PETSC_VIEWER_STDOUT_SELF);
 53:   ISDestroy(&rperm);
 54:   ISDestroy(&cperm);

 56:   /* create Cperm = rperm*C*icperm */
 57:   PetscOptionsGetBool(NULL,"-testmyordering",&TestMyorder,NULL);
 58:   if (TestMyorder) {
 59:     MatGetOrdering_myordering(C,MATORDERINGQMD,&rperm,&cperm);
 60:     printf("myordering's rperm:\n");
 61:     ISView(rperm,PETSC_VIEWER_STDOUT_SELF);
 62:     ISInvertPermutation(cperm,PETSC_DECIDE,&icperm);
 63:     ISGetIndices(rperm,&rperm_ptr);
 64:     ISGetIndices(icperm,&cperm_ptr);
 65:     MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&Cperm);
 66:     for (i=0; i<m*n; i++) {
 67:       MatGetRow(C,rperm_ptr[i],&ncols,&cols,&vals);
 68:       for (j=0; j<ncols; j++) {
 69:         /* printf(" (%d %d %g)\n",i,cperm_ptr[cols[j]],vals[j]); */
 70:         MatSetValues(Cperm,1,&i,1,&cperm_ptr[cols[j]],&vals[j],INSERT_VALUES);
 71:       }
 72:     }
 73:     MatAssemblyBegin(Cperm,MAT_FINAL_ASSEMBLY);
 74:     MatAssemblyEnd(Cperm,MAT_FINAL_ASSEMBLY);
 75:     ISRestoreIndices(rperm,&rperm_ptr);
 76:     ISRestoreIndices(icperm,&cperm_ptr);

 78:     ISDestroy(&rperm);
 79:     ISDestroy(&cperm);
 80:     ISDestroy(&icperm);
 81:     MatDestroy(&Cperm);
 82:   }

 84:   MatDestroy(&C);
 85:   PetscFinalize();
 86:   return 0;
 87: }

 89: #include <petsc-private/matimpl.h>
 90: /* This is modified from MatGetOrdering_Natural() */
 93: PetscErrorCode MatGetOrdering_myordering(Mat mat,MatOrderingType type,IS *irow,IS *icol)
 94: {
 96:   PetscInt       n,i,*ii;
 97:   PetscBool      done;
 98:   MPI_Comm       comm;

101:   PetscObjectGetComm((PetscObject)mat,&comm);
102:   MatGetRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,&n,NULL,NULL,&done);
103:   MatRestoreRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,NULL,NULL,NULL,&done);
104:   if (done) { /* matrix may be "compressed" in symbolic factorization, due to i-nodes or block storage */
105:     PetscMalloc1(n,&ii);
106:     for (i=0; i<n; i++) ii[i] = n-i-1; /* replace your index here */
107:     ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_COPY_VALUES,irow);
108:     ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_OWN_POINTER,icol);
109:   } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"MatRestoreRowIJ fails!");
110:   ISSetIdentity(*irow);
111:   ISSetIdentity(*icol);

113:   ISSetPermutation(*irow);
114:   ISSetPermutation(*icol);
115:   return(0);
116: }