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 PetscFunctionBeginUser; 19 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 20 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 21 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 22 23 /* create the matrix for the five point stencil, YET AGAIN */ 24 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * n, m * n, 5, NULL, &C)); 25 PetscCall(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 PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 33 } 34 if (i < m - 1) { 35 J = Ii + n; 36 PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 37 } 38 if (j > 0) { 39 J = Ii - 1; 40 PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 41 } 42 if (j < n - 1) { 43 J = Ii + 1; 44 PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 45 } 46 v = 4.0; 47 PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 48 } 49 } 50 PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 51 PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 52 53 PetscCall(MatGetOrdering(C, MATORDERINGND, &rperm, &cperm)); 54 PetscCall(ISView(rperm, PETSC_VIEWER_STDOUT_SELF)); 55 PetscCall(ISDestroy(&rperm)); 56 PetscCall(ISDestroy(&cperm)); 57 58 PetscCall(MatGetOrdering(C, MATORDERINGRCM, &rperm, &cperm)); 59 PetscCall(ISView(rperm, PETSC_VIEWER_STDOUT_SELF)); 60 PetscCall(ISDestroy(&rperm)); 61 PetscCall(ISDestroy(&cperm)); 62 63 PetscCall(MatGetOrdering(C, MATORDERINGQMD, &rperm, &cperm)); 64 PetscCall(ISView(rperm, PETSC_VIEWER_STDOUT_SELF)); 65 PetscCall(ISDestroy(&rperm)); 66 PetscCall(ISDestroy(&cperm)); 67 68 /* create Cperm = rperm*C*icperm */ 69 PetscCall(PetscOptionsGetBool(NULL, NULL, "-testmyordering", &TestMyorder, NULL)); 70 if (TestMyorder) { 71 PetscCall(MatGetOrdering_myordering(C, MATORDERINGQMD, &rperm, &cperm)); 72 printf("myordering's rperm:\n"); 73 PetscCall(ISView(rperm, PETSC_VIEWER_STDOUT_SELF)); 74 PetscCall(ISInvertPermutation(cperm, PETSC_DECIDE, &icperm)); 75 PetscCall(ISGetIndices(rperm, &rperm_ptr)); 76 PetscCall(ISGetIndices(icperm, &cperm_ptr)); 77 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * n, m * n, 5, NULL, &Cperm)); 78 for (i = 0; i < m * n; i++) { 79 PetscCall(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 PetscCall(MatSetValues(Cperm, 1, &i, 1, &cperm_ptr[cols[j]], &vals[j], INSERT_VALUES)); 83 } 84 } 85 PetscCall(MatAssemblyBegin(Cperm, MAT_FINAL_ASSEMBLY)); 86 PetscCall(MatAssemblyEnd(Cperm, MAT_FINAL_ASSEMBLY)); 87 PetscCall(ISRestoreIndices(rperm, &rperm_ptr)); 88 PetscCall(ISRestoreIndices(icperm, &cperm_ptr)); 89 90 PetscCall(ISDestroy(&rperm)); 91 PetscCall(ISDestroy(&cperm)); 92 PetscCall(ISDestroy(&icperm)); 93 PetscCall(MatDestroy(&Cperm)); 94 } 95 96 PetscCall(MatDestroy(&C)); 97 PetscCall(PetscFinalize()); 98 return 0; 99 } 100 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; 108 109 PetscFunctionBegin; 110 PetscCall(PetscObjectGetComm((PetscObject)mat, &comm)); 111 PetscCall(MatGetRowIJ(mat, 0, PETSC_FALSE, PETSC_TRUE, &n, NULL, NULL, &done)); 112 PetscCall(MatRestoreRowIJ(mat, 0, PETSC_FALSE, PETSC_TRUE, NULL, NULL, NULL, &done)); 113 if (done) { /* matrix may be "compressed" in symbolic factorization, due to i-nodes or block storage */ 114 PetscCall(PetscMalloc1(n, &ii)); 115 for (i = 0; i < n; i++) ii[i] = n - i - 1; /* replace your index here */ 116 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, ii, PETSC_COPY_VALUES, irow)); 117 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, ii, PETSC_OWN_POINTER, icol)); 118 } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "MatRestoreRowIJ fails!"); 119 PetscCall(ISSetPermutation(*irow)); 120 PetscCall(ISSetPermutation(*icol)); 121 PetscFunctionReturn(PETSC_SUCCESS); 122 } 123 124 /*TEST 125 126 test: 127 128 TEST*/ 129