Actual source code: ex22.c


  2: static char help[] = "Tests matrix ordering routines.\n\n";

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

  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;

 19:   PetscInitialize(&argc, &args, (char *)0, help);
 20:   MPI_Comm_size(PETSC_COMM_WORLD, &size);

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

 53:   MatGetOrdering(C, MATORDERINGND, &rperm, &cperm);
 54:   ISView(rperm, PETSC_VIEWER_STDOUT_SELF);
 55:   ISDestroy(&rperm);
 56:   ISDestroy(&cperm);

 58:   MatGetOrdering(C, MATORDERINGRCM, &rperm, &cperm);
 59:   ISView(rperm, PETSC_VIEWER_STDOUT_SELF);
 60:   ISDestroy(&rperm);
 61:   ISDestroy(&cperm);

 63:   MatGetOrdering(C, MATORDERINGQMD, &rperm, &cperm);
 64:   ISView(rperm, PETSC_VIEWER_STDOUT_SELF);
 65:   ISDestroy(&rperm);
 66:   ISDestroy(&cperm);

 68:   /* create Cperm = rperm*C*icperm */
 69:   PetscOptionsGetBool(NULL, NULL, "-testmyordering", &TestMyorder, NULL);
 70:   if (TestMyorder) {
 71:     MatGetOrdering_myordering(C, MATORDERINGQMD, &rperm, &cperm);
 72:     printf("myordering's rperm:\n");
 73:     ISView(rperm, PETSC_VIEWER_STDOUT_SELF);
 74:     ISInvertPermutation(cperm, PETSC_DECIDE, &icperm);
 75:     ISGetIndices(rperm, &rperm_ptr);
 76:     ISGetIndices(icperm, &cperm_ptr);
 77:     MatCreateSeqAIJ(PETSC_COMM_SELF, m * n, m * n, 5, NULL, &Cperm);
 78:     for (i = 0; i < m * n; i++) {
 79:       MatGetRow(C, rperm_ptr[i], &ncols, &cols, &vals);
 80:       for (j = 0; j < ncols; j++) {
 81:         /* printf(" (%d %d %g)\n",i,cperm_ptr[cols[j]],vals[j]); */
 82:         MatSetValues(Cperm, 1, &i, 1, &cperm_ptr[cols[j]], &vals[j], INSERT_VALUES);
 83:       }
 84:     }
 85:     MatAssemblyBegin(Cperm, MAT_FINAL_ASSEMBLY);
 86:     MatAssemblyEnd(Cperm, MAT_FINAL_ASSEMBLY);
 87:     ISRestoreIndices(rperm, &rperm_ptr);
 88:     ISRestoreIndices(icperm, &cperm_ptr);

 90:     ISDestroy(&rperm);
 91:     ISDestroy(&cperm);
 92:     ISDestroy(&icperm);
 93:     MatDestroy(&Cperm);
 94:   }

 96:   MatDestroy(&C);
 97:   PetscFinalize();
 98:   return 0;
 99: }

101: #include <petsc/private/matimpl.h>
102: /* This is modified from MatGetOrdering_Natural() */
103: PetscErrorCode MatGetOrdering_myordering(Mat mat, MatOrderingType type, IS *irow, IS *icol)
104: {
105:   PetscInt  n, i, *ii;
106:   PetscBool done;
107:   MPI_Comm  comm;

109:   PetscObjectGetComm((PetscObject)mat, &comm);
110:   MatGetRowIJ(mat, 0, PETSC_FALSE, PETSC_TRUE, &n, NULL, NULL, &done);
111:   MatRestoreRowIJ(mat, 0, PETSC_FALSE, PETSC_TRUE, NULL, NULL, NULL, &done);
112:   if (done) { /* matrix may be "compressed" in symbolic factorization, due to i-nodes or block storage */
113:     PetscMalloc1(n, &ii);
114:     for (i = 0; i < n; i++) ii[i] = n - i - 1; /* replace your index here */
115:     ISCreateGeneral(PETSC_COMM_SELF, n, ii, PETSC_COPY_VALUES, irow);
116:     ISCreateGeneral(PETSC_COMM_SELF, n, ii, PETSC_OWN_POINTER, icol);
117:   } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "MatRestoreRowIJ fails!");
118:   ISSetPermutation(*irow);
119:   ISSetPermutation(*icol);
120:   return 0;
121: }

123: /*TEST

125:    test:

127: TEST*/