xref: /petsc/src/mat/tests/ex22.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
1c4762a1bSJed Brown static char help[] = "Tests matrix ordering routines.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown extern PetscErrorCode MatGetOrdering_myordering(Mat, MatOrderingType, IS *, IS *);
5c4762a1bSJed Brown 
main(int argc,char ** args)6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat                C, Cperm;
9c4762a1bSJed Brown   PetscInt           i, j, m = 5, n = 5, Ii, J, ncols;
10c4762a1bSJed Brown   PetscScalar        v;
11c4762a1bSJed Brown   PetscMPIInt        size;
12c4762a1bSJed Brown   IS                 rperm, cperm, icperm;
13c4762a1bSJed Brown   const PetscInt    *rperm_ptr, *cperm_ptr, *cols;
14c4762a1bSJed Brown   const PetscScalar *vals;
15c4762a1bSJed Brown   PetscBool          TestMyorder = PETSC_FALSE;
16c4762a1bSJed Brown 
17327415f7SBarry Smith   PetscFunctionBeginUser;
18c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
20be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN */
239566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * n, m * n, 5, NULL, &C));
249566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
25c4762a1bSJed Brown   for (i = 0; i < m; i++) {
26c4762a1bSJed Brown     for (j = 0; j < n; j++) {
279371c9d4SSatish Balay       v  = -1.0;
289371c9d4SSatish Balay       Ii = j + n * i;
299371c9d4SSatish Balay       if (i > 0) {
309371c9d4SSatish Balay         J = Ii - n;
319371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
329371c9d4SSatish Balay       }
339371c9d4SSatish Balay       if (i < m - 1) {
349371c9d4SSatish Balay         J = Ii + n;
359371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
369371c9d4SSatish Balay       }
379371c9d4SSatish Balay       if (j > 0) {
389371c9d4SSatish Balay         J = Ii - 1;
399371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
409371c9d4SSatish Balay       }
419371c9d4SSatish Balay       if (j < n - 1) {
429371c9d4SSatish Balay         J = Ii + 1;
439371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
449371c9d4SSatish Balay       }
459371c9d4SSatish Balay       v = 4.0;
469371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
47c4762a1bSJed Brown     }
48c4762a1bSJed Brown   }
499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(C, MATORDERINGND, &rperm, &cperm));
539566063dSJacob Faibussowitsch   PetscCall(ISView(rperm, PETSC_VIEWER_STDOUT_SELF));
549566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&rperm));
559566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cperm));
56c4762a1bSJed Brown 
579566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(C, MATORDERINGRCM, &rperm, &cperm));
589566063dSJacob Faibussowitsch   PetscCall(ISView(rperm, PETSC_VIEWER_STDOUT_SELF));
599566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&rperm));
609566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cperm));
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(C, MATORDERINGQMD, &rperm, &cperm));
639566063dSJacob Faibussowitsch   PetscCall(ISView(rperm, PETSC_VIEWER_STDOUT_SELF));
649566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&rperm));
659566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cperm));
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   /* create Cperm = rperm*C*icperm */
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-testmyordering", &TestMyorder, NULL));
69c4762a1bSJed Brown   if (TestMyorder) {
709566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering_myordering(C, MATORDERINGQMD, &rperm, &cperm));
71c4762a1bSJed Brown     printf("myordering's rperm:\n");
729566063dSJacob Faibussowitsch     PetscCall(ISView(rperm, PETSC_VIEWER_STDOUT_SELF));
739566063dSJacob Faibussowitsch     PetscCall(ISInvertPermutation(cperm, PETSC_DECIDE, &icperm));
749566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(rperm, &rperm_ptr));
759566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(icperm, &cperm_ptr));
769566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * n, m * n, 5, NULL, &Cperm));
77c4762a1bSJed Brown     for (i = 0; i < m * n; i++) {
789566063dSJacob Faibussowitsch       PetscCall(MatGetRow(C, rperm_ptr[i], &ncols, &cols, &vals));
79c4762a1bSJed Brown       for (j = 0; j < ncols; j++) {
80c4762a1bSJed Brown         /* printf(" (%d %d %g)\n",i,cperm_ptr[cols[j]],vals[j]); */
819566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Cperm, 1, &i, 1, &cperm_ptr[cols[j]], &vals[j], INSERT_VALUES));
82c4762a1bSJed Brown       }
83c4762a1bSJed Brown     }
849566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Cperm, MAT_FINAL_ASSEMBLY));
859566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Cperm, MAT_FINAL_ASSEMBLY));
869566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(rperm, &rperm_ptr));
879566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(icperm, &cperm_ptr));
88c4762a1bSJed Brown 
899566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&rperm));
909566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cperm));
919566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&icperm));
929566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Cperm));
93c4762a1bSJed Brown   }
94c4762a1bSJed Brown 
959566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
969566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
97b122ec5aSJacob Faibussowitsch   return 0;
98c4762a1bSJed Brown }
99c4762a1bSJed Brown 
100c4762a1bSJed Brown #include <petsc/private/matimpl.h>
101c4762a1bSJed Brown /* This is modified from MatGetOrdering_Natural() */
MatGetOrdering_myordering(Mat mat,MatOrderingType type,IS * irow,IS * icol)102d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetOrdering_myordering(Mat mat, MatOrderingType type, IS *irow, IS *icol)
103d71ae5a4SJacob Faibussowitsch {
104c4762a1bSJed Brown   PetscInt  n, i, *ii;
105c4762a1bSJed Brown   PetscBool done;
106c4762a1bSJed Brown   MPI_Comm  comm;
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   PetscFunctionBegin;
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
1109566063dSJacob Faibussowitsch   PetscCall(MatGetRowIJ(mat, 0, PETSC_FALSE, PETSC_TRUE, &n, NULL, NULL, &done));
1119566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowIJ(mat, 0, PETSC_FALSE, PETSC_TRUE, NULL, NULL, NULL, &done));
112*966bd95aSPierre Jolivet   PetscCheck(done, PETSC_COMM_WORLD, PETSC_ERR_SUP, "MatRestoreRowIJ fails!");
113*966bd95aSPierre Jolivet   /* 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));
1189566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*irow));
1199566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*icol));
1203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
121c4762a1bSJed Brown }
122c4762a1bSJed Brown 
123c4762a1bSJed Brown /*TEST
124c4762a1bSJed Brown 
125c4762a1bSJed Brown    test:
126c4762a1bSJed Brown 
127c4762a1bSJed Brown TEST*/
128