1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests matrix ordering routines.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown extern PetscErrorCode MatGetOrdering_myordering(Mat, MatOrderingType, IS *, IS *); 6c4762a1bSJed Brown 7d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 8d71ae5a4SJacob Faibussowitsch { 9c4762a1bSJed Brown Mat C, Cperm; 10c4762a1bSJed Brown PetscInt i, j, m = 5, n = 5, Ii, J, ncols; 11c4762a1bSJed Brown PetscScalar v; 12c4762a1bSJed Brown PetscMPIInt size; 13c4762a1bSJed Brown IS rperm, cperm, icperm; 14c4762a1bSJed Brown const PetscInt *rperm_ptr, *cperm_ptr, *cols; 15c4762a1bSJed Brown const PetscScalar *vals; 16c4762a1bSJed Brown PetscBool TestMyorder = PETSC_FALSE; 17c4762a1bSJed Brown 18327415f7SBarry Smith PetscFunctionBeginUser; 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 21be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* create the matrix for the five point stencil, YET AGAIN */ 249566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * n, m * n, 5, NULL, &C)); 259566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 26c4762a1bSJed Brown for (i = 0; i < m; i++) { 27c4762a1bSJed Brown for (j = 0; j < n; j++) { 289371c9d4SSatish Balay v = -1.0; 299371c9d4SSatish Balay Ii = j + n * i; 309371c9d4SSatish Balay if (i > 0) { 319371c9d4SSatish Balay J = Ii - n; 329371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 339371c9d4SSatish Balay } 349371c9d4SSatish Balay if (i < m - 1) { 359371c9d4SSatish Balay J = Ii + n; 369371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 379371c9d4SSatish Balay } 389371c9d4SSatish Balay if (j > 0) { 399371c9d4SSatish Balay J = Ii - 1; 409371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 419371c9d4SSatish Balay } 429371c9d4SSatish Balay if (j < n - 1) { 439371c9d4SSatish Balay J = Ii + 1; 449371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 459371c9d4SSatish Balay } 469371c9d4SSatish Balay v = 4.0; 479371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 48c4762a1bSJed Brown } 49c4762a1bSJed Brown } 509566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 519566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 52c4762a1bSJed Brown 539566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C, MATORDERINGND, &rperm, &cperm)); 549566063dSJacob Faibussowitsch PetscCall(ISView(rperm, PETSC_VIEWER_STDOUT_SELF)); 559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&rperm)); 569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cperm)); 57c4762a1bSJed Brown 589566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C, MATORDERINGRCM, &rperm, &cperm)); 599566063dSJacob Faibussowitsch PetscCall(ISView(rperm, PETSC_VIEWER_STDOUT_SELF)); 609566063dSJacob Faibussowitsch PetscCall(ISDestroy(&rperm)); 619566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cperm)); 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C, MATORDERINGQMD, &rperm, &cperm)); 649566063dSJacob Faibussowitsch PetscCall(ISView(rperm, PETSC_VIEWER_STDOUT_SELF)); 659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&rperm)); 669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cperm)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* create Cperm = rperm*C*icperm */ 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-testmyordering", &TestMyorder, NULL)); 70c4762a1bSJed Brown if (TestMyorder) { 719566063dSJacob Faibussowitsch PetscCall(MatGetOrdering_myordering(C, MATORDERINGQMD, &rperm, &cperm)); 72c4762a1bSJed Brown printf("myordering's rperm:\n"); 739566063dSJacob Faibussowitsch PetscCall(ISView(rperm, PETSC_VIEWER_STDOUT_SELF)); 749566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(cperm, PETSC_DECIDE, &icperm)); 759566063dSJacob Faibussowitsch PetscCall(ISGetIndices(rperm, &rperm_ptr)); 769566063dSJacob Faibussowitsch PetscCall(ISGetIndices(icperm, &cperm_ptr)); 779566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * n, m * n, 5, NULL, &Cperm)); 78c4762a1bSJed Brown for (i = 0; i < m * n; i++) { 799566063dSJacob Faibussowitsch PetscCall(MatGetRow(C, rperm_ptr[i], &ncols, &cols, &vals)); 80c4762a1bSJed Brown for (j = 0; j < ncols; j++) { 81c4762a1bSJed Brown /* printf(" (%d %d %g)\n",i,cperm_ptr[cols[j]],vals[j]); */ 829566063dSJacob Faibussowitsch PetscCall(MatSetValues(Cperm, 1, &i, 1, &cperm_ptr[cols[j]], &vals[j], INSERT_VALUES)); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown } 859566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Cperm, MAT_FINAL_ASSEMBLY)); 869566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Cperm, MAT_FINAL_ASSEMBLY)); 879566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(rperm, &rperm_ptr)); 889566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(icperm, &cperm_ptr)); 89c4762a1bSJed Brown 909566063dSJacob Faibussowitsch PetscCall(ISDestroy(&rperm)); 919566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cperm)); 929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icperm)); 939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cperm)); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 979566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 98b122ec5aSJacob Faibussowitsch return 0; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown #include <petsc/private/matimpl.h> 102c4762a1bSJed Brown /* This is modified from MatGetOrdering_Natural() */ 103d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetOrdering_myordering(Mat mat, MatOrderingType type, IS *irow, IS *icol) 104d71ae5a4SJacob Faibussowitsch { 105c4762a1bSJed Brown PetscInt n, i, *ii; 106c4762a1bSJed Brown PetscBool done; 107c4762a1bSJed Brown MPI_Comm comm; 108c4762a1bSJed Brown 109c4762a1bSJed Brown PetscFunctionBegin; 1109566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)mat, &comm)); 1119566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(mat, 0, PETSC_FALSE, PETSC_TRUE, &n, NULL, NULL, &done)); 1129566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(mat, 0, PETSC_FALSE, PETSC_TRUE, NULL, NULL, NULL, &done)); 113c4762a1bSJed Brown if (done) { /* matrix may be "compressed" in symbolic factorization, due to i-nodes or block storage */ 1149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &ii)); 115c4762a1bSJed Brown for (i = 0; i < n; i++) ii[i] = n - i - 1; /* replace your index here */ 1169566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, ii, PETSC_COPY_VALUES, irow)); 1179566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, n, ii, PETSC_OWN_POINTER, icol)); 118c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "MatRestoreRowIJ fails!"); 1199566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*irow)); 1209566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(*icol)); 121*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown 124c4762a1bSJed Brown /*TEST 125c4762a1bSJed Brown 126c4762a1bSJed Brown test: 127c4762a1bSJed Brown 128c4762a1bSJed Brown TEST*/ 129