xref: /petsc/src/mat/tests/ex166.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static char help[] = "Tests MatPermute() for a square matrix in parallel.\n\n";
2c4762a1bSJed Brown #include <petscmat.h>
3c4762a1bSJed Brown /* Results:
4c4762a1bSJed Brown    Sequential:
5c4762a1bSJed Brown    - seqaij:   correct permutation
6c4762a1bSJed Brown    - seqbaij:  permutation not supported for this MATTYPE
7c4762a1bSJed Brown    - seqsbaij: permutation not supported for this MATTYPE
8c4762a1bSJed Brown    Parallel:
9c4762a1bSJed Brown    - mpiaij:   incorrect permutation (disable this op for now)
10c4762a1bSJed Brown    - seqbaij:  correct permutation, but all manner of memory leaks
11c4762a1bSJed Brown                (disable this op for now, because it appears broken for nonsquare matrices; see ex151)
12c4762a1bSJed Brown    - seqsbaij: permutation not supported for this MATTYPE
13c4762a1bSJed Brown 
14c4762a1bSJed Brown  */
159371c9d4SSatish Balay int main(int argc, char **argv) {
169371c9d4SSatish Balay   const struct {
179371c9d4SSatish Balay     PetscInt    i, j;
189371c9d4SSatish Balay     PetscScalar v;
199371c9d4SSatish Balay   } entries[] = {
209371c9d4SSatish Balay     {0, 3, 1.},
219371c9d4SSatish Balay     {1, 2, 2.},
229371c9d4SSatish Balay     {2, 1, 3.},
239371c9d4SSatish Balay     {2, 4, 4.},
249371c9d4SSatish Balay     {3, 0, 5.},
259371c9d4SSatish Balay     {3, 3, 6.},
269371c9d4SSatish Balay     {4, 1, 7.},
279371c9d4SSatish Balay     {4, 4, 8.}
289371c9d4SSatish Balay   };
29c4762a1bSJed Brown   const PetscInt ixrow[5] = {4, 2, 1, 3, 0}, ixcol[5] = {3, 2, 1, 4, 0};
30c4762a1bSJed Brown   Mat            A, B;
31c4762a1bSJed Brown   PetscInt       i, rstart, rend, cstart, cend;
32c4762a1bSJed Brown   IS             isrow, iscol;
33c4762a1bSJed Brown   PetscViewer    viewer, sviewer;
34c4762a1bSJed Brown   PetscBool      view_sparse;
35c4762a1bSJed Brown 
36327415f7SBarry Smith   PetscFunctionBeginUser;
379566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
38c4762a1bSJed Brown   /* ------- Assemble matrix, --------- */
399566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
409566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 5, 5));
419566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
429566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
439566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
449566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(A, &cstart, &cend));
45c4762a1bSJed Brown 
46*48a46eb9SPierre Jolivet   for (i = 0; i < (PetscInt)PETSC_STATIC_ARRAY_LENGTH(entries); i++) PetscCall(MatSetValue(A, entries[i].i, entries[i].j, entries[i].v, INSERT_VALUES));
479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* ------ Prepare index sets ------ */
519566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, rend - rstart, ixrow + rstart, PETSC_USE_POINTER, &isrow));
529566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 5, ixcol, PETSC_USE_POINTER, &iscol));
539566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(isrow));
549566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(iscol));
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD, &viewer));
57c4762a1bSJed Brown   view_sparse = PETSC_FALSE;
589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_sparse", &view_sparse, NULL));
59*48a46eb9SPierre Jolivet   if (!view_sparse) PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_DENSE));
609566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Original matrix\n"));
619566063dSJacob Faibussowitsch   PetscCall(MatView(A, viewer));
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(MatPermute(A, isrow, iscol, &B));
649566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Permuted matrix\n"));
659566063dSJacob Faibussowitsch   PetscCall(MatView(B, viewer));
66*48a46eb9SPierre Jolivet   if (!view_sparse) PetscCall(PetscViewerPopFormat(viewer));
67c4762a1bSJed Brown 
689566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Row permutation\n"));
699566063dSJacob Faibussowitsch   PetscCall(ISView(isrow, viewer));
709566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Column permutation\n"));
719566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
729566063dSJacob Faibussowitsch   PetscCall(ISView(iscol, sviewer));
739566063dSJacob Faibussowitsch   PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   /* Free data structures */
769566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
779566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
80c4762a1bSJed Brown 
819566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
82b122ec5aSJacob Faibussowitsch   return 0;
83c4762a1bSJed Brown }
84