xref: /petsc/src/mat/tests/ex166.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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  */
15c4762a1bSJed Brown int main(int argc,char **argv)
16c4762a1bSJed Brown {
17c4762a1bSJed Brown   const struct {PetscInt i,j; PetscScalar v;} entries[] = {{0,3,1.},{1,2,2.},{2,1,3.},{2,4,4.},{3,0,5.},{3,3,6.},{4,1,7.},{4,4,8.}};
18c4762a1bSJed Brown   const PetscInt ixrow[5]                               = {4,2,1,3,0},ixcol[5] = {3,2,1,4,0};
19c4762a1bSJed Brown   Mat            A,B;
20c4762a1bSJed Brown   PetscErrorCode ierr;
21c4762a1bSJed Brown   PetscInt       i,rstart,rend,cstart,cend;
22c4762a1bSJed Brown   IS             isrow,iscol;
23c4762a1bSJed Brown   PetscViewer    viewer,sviewer;
24c4762a1bSJed Brown   PetscBool      view_sparse;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
27c4762a1bSJed Brown   /* ------- Assemble matrix, --------- */
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,5,5));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRangeColumn(A,&cstart,&cend));
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   for (i=0; i<(PetscInt)(sizeof(entries)/sizeof(entries[0])); i++) {
36*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValue(A,entries[i].i,entries[i].j,entries[i].v,INSERT_VALUES));
37c4762a1bSJed Brown   }
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* ------ Prepare index sets ------ */
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,rend-rstart,ixrow+rstart,PETSC_USE_POINTER,&isrow));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,5,ixcol,PETSC_USE_POINTER,&iscol));
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISSetPermutation(isrow));
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISSetPermutation(iscol));
46c4762a1bSJed Brown 
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
48c4762a1bSJed Brown   view_sparse = PETSC_FALSE;
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL, "-view_sparse", &view_sparse, NULL));
50c4762a1bSJed Brown   if (!view_sparse) {
51*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_DENSE));
52c4762a1bSJed Brown   }
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Original matrix\n"));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,viewer));
55c4762a1bSJed Brown 
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPermute(A,isrow,iscol,&B));
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Permuted matrix\n"));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,viewer));
59c4762a1bSJed Brown   if (!view_sparse) {
60*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPopFormat(viewer));
61c4762a1bSJed Brown   }
62c4762a1bSJed Brown 
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Row permutation\n"));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(isrow,viewer));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Column permutation\n"));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(iscol,sviewer));
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   /* Free data structures */
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrow));
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iscol));
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   ierr = PetscFinalize();
77c4762a1bSJed Brown   return ierr;
78c4762a1bSJed Brown }
79