xref: /petsc/src/mat/tests/ex68.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatReorderForNonzeroDiagonal().\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **argv)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            mat,B,C;
9c4762a1bSJed Brown   PetscInt       i,j;
10c4762a1bSJed Brown   PetscMPIInt    size;
11c4762a1bSJed Brown   PetscScalar    v;
12c4762a1bSJed Brown   IS             isrow,iscol,identity;
13c4762a1bSJed Brown   PetscViewer    viewer;
14c4762a1bSJed Brown 
15*327415f7SBarry Smith   PetscFunctionBeginUser;
169566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   /* ------- Assemble matrix, --------- */
19c4762a1bSJed Brown 
209566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&mat));
219566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,4,4));
229566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(mat));
239566063dSJacob Faibussowitsch   PetscCall(MatSetUp(mat));
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   /* set anti-diagonal of matrix */
26c4762a1bSJed Brown   v    = 1.0;
27c4762a1bSJed Brown   i    = 0; j = 3;
289566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES));
29c4762a1bSJed Brown   v    = 2.0;
30c4762a1bSJed Brown   i    = 1; j = 2;
319566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES));
32c4762a1bSJed Brown   v    = 3.0;
33c4762a1bSJed Brown   i    = 2; j = 1;
349566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES));
35c4762a1bSJed Brown   v    = 4.0;
36c4762a1bSJed Brown   i    = 3; j = 0;
379566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES));
389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
40c4762a1bSJed Brown 
419566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
429566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_DENSE));
439566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix\n"));
449566063dSJacob Faibussowitsch   PetscCall(MatView(mat,viewer));
45c4762a1bSJed Brown 
469566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(mat,MATORDERINGNATURAL,&isrow,&iscol));
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(MatPermute(mat,isrow,iscol,&B));
499566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by identity\n"));
509566063dSJacob Faibussowitsch   PetscCall(MatView(B,viewer));
519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
52c4762a1bSJed Brown 
539566063dSJacob Faibussowitsch   PetscCall(MatReorderForNonzeroDiagonal(mat,1.e-8,isrow,iscol));
549566063dSJacob Faibussowitsch   PetscCall(MatPermute(mat,isrow,iscol,&B));
559566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by identity + NonzeroDiagonal()\n"));
569566063dSJacob Faibussowitsch   PetscCall(MatView(B,viewer));
579566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Row permutation\n"));
589566063dSJacob Faibussowitsch   PetscCall(ISView(isrow,viewer));
599566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Column permutation\n"));
609566063dSJacob Faibussowitsch   PetscCall(ISView(iscol,viewer));
619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
649566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(mat,MATORDERINGND,&isrow,&iscol));
679566063dSJacob Faibussowitsch   PetscCall(MatPermute(mat,isrow,iscol,&B));
689566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by ND\n"));
699566063dSJacob Faibussowitsch   PetscCall(MatView(B,viewer));
709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
719566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"ND row permutation\n"));
729566063dSJacob Faibussowitsch   PetscCall(ISView(isrow,viewer));
739566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"ND column permutation\n"));
749566063dSJacob Faibussowitsch   PetscCall(ISView(iscol,viewer));
75c4762a1bSJed Brown 
769566063dSJacob Faibussowitsch   PetscCall(MatReorderForNonzeroDiagonal(mat,1.e-8,isrow,iscol));
779566063dSJacob Faibussowitsch   PetscCall(MatPermute(mat,isrow,iscol,&B));
789566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by ND + NonzeroDiagonal()\n"));
799566063dSJacob Faibussowitsch   PetscCall(MatView(B,viewer));
809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
819566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"ND + NonzeroDiagonal() row permutation\n"));
829566063dSJacob Faibussowitsch   PetscCall(ISView(isrow,viewer));
839566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"ND + NonzeroDiagonal() column permutation\n"));
849566063dSJacob Faibussowitsch   PetscCall(ISView(iscol,viewer));
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
879566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
88c4762a1bSJed Brown 
899566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(mat,MATORDERINGRCM,&isrow,&iscol));
909566063dSJacob Faibussowitsch   PetscCall(MatPermute(mat,isrow,iscol,&B));
919566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by RCM\n"));
929566063dSJacob Faibussowitsch   PetscCall(MatView(B,viewer));
939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
949566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"RCM row permutation\n"));
959566063dSJacob Faibussowitsch   PetscCall(ISView(isrow,viewer));
969566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"RCM column permutation\n"));
979566063dSJacob Faibussowitsch   PetscCall(ISView(iscol,viewer));
98c4762a1bSJed Brown 
999566063dSJacob Faibussowitsch   PetscCall(MatReorderForNonzeroDiagonal(mat,1.e-8,isrow,iscol));
1009566063dSJacob Faibussowitsch   PetscCall(MatPermute(mat,isrow,iscol,&B));
1019566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by RCM + NonzeroDiagonal()\n"));
1029566063dSJacob Faibussowitsch   PetscCall(MatView(B,viewer));
1039566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"RCM + NonzeroDiagonal() row permutation\n"));
1049566063dSJacob Faibussowitsch   PetscCall(ISView(isrow,viewer));
1059566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"RCM + NonzeroDiagonal() column permutation\n"));
1069566063dSJacob Faibussowitsch   PetscCall(ISView(iscol,viewer));
1079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
108c4762a1bSJed Brown   if (size == 1) {
1099566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE));
1109566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,4,0,1,&identity));
1119566063dSJacob Faibussowitsch     PetscCall(MatPermute(B,identity,identity,&C));
1129566063dSJacob Faibussowitsch     PetscCall(MatConvert(C,MATSEQSBAIJ,MAT_INPLACE_MATRIX,&C));
1139566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
1149566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&identity));
115c4762a1bSJed Brown   }
1169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
117c4762a1bSJed Brown   /* Test MatLUFactor(); set diagonal as zeros as requested by PETSc matrix factorization */
118c4762a1bSJed Brown   for (i=0; i<4; i++) {
119c4762a1bSJed Brown     v = 0.0;
1209566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat,1,&i,1,&i,&v,INSERT_VALUES));
121c4762a1bSJed Brown   }
1229566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
1239566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
1249566063dSJacob Faibussowitsch   PetscCall(MatLUFactor(mat,isrow,iscol,NULL));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* Free data structures */
1279566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
1289566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
1299566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
130c4762a1bSJed Brown 
1319566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
132b122ec5aSJacob Faibussowitsch   return 0;
133c4762a1bSJed Brown }
134c4762a1bSJed Brown 
135c4762a1bSJed Brown /*TEST
136c4762a1bSJed Brown 
137c4762a1bSJed Brown    test:
138c4762a1bSJed Brown 
139c4762a1bSJed Brown TEST*/
140