xref: /petsc/src/mat/tests/ex68.c (revision 9566063d113dddea24716c546802770db7481bc0)
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*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   /* ------- Assemble matrix, --------- */
18c4762a1bSJed Brown 
19*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&mat));
20*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,4,4));
21*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(mat));
22*9566063dSJacob Faibussowitsch   PetscCall(MatSetUp(mat));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /* set anti-diagonal of matrix */
25c4762a1bSJed Brown   v    = 1.0;
26c4762a1bSJed Brown   i    = 0; j = 3;
27*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES));
28c4762a1bSJed Brown   v    = 2.0;
29c4762a1bSJed Brown   i    = 1; j = 2;
30*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES));
31c4762a1bSJed Brown   v    = 3.0;
32c4762a1bSJed Brown   i    = 2; j = 1;
33*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES));
34c4762a1bSJed Brown   v    = 4.0;
35c4762a1bSJed Brown   i    = 3; j = 0;
36*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES));
37*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
38*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
39c4762a1bSJed Brown 
40*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
41*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_DENSE));
42*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix\n"));
43*9566063dSJacob Faibussowitsch   PetscCall(MatView(mat,viewer));
44c4762a1bSJed Brown 
45*9566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(mat,MATORDERINGNATURAL,&isrow,&iscol));
46c4762a1bSJed Brown 
47*9566063dSJacob Faibussowitsch   PetscCall(MatPermute(mat,isrow,iscol,&B));
48*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by identity\n"));
49*9566063dSJacob Faibussowitsch   PetscCall(MatView(B,viewer));
50*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
51c4762a1bSJed Brown 
52*9566063dSJacob Faibussowitsch   PetscCall(MatReorderForNonzeroDiagonal(mat,1.e-8,isrow,iscol));
53*9566063dSJacob Faibussowitsch   PetscCall(MatPermute(mat,isrow,iscol,&B));
54*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by identity + NonzeroDiagonal()\n"));
55*9566063dSJacob Faibussowitsch   PetscCall(MatView(B,viewer));
56*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Row permutation\n"));
57*9566063dSJacob Faibussowitsch   PetscCall(ISView(isrow,viewer));
58*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Column permutation\n"));
59*9566063dSJacob Faibussowitsch   PetscCall(ISView(iscol,viewer));
60*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
61c4762a1bSJed Brown 
62*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
63*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
64c4762a1bSJed Brown 
65*9566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(mat,MATORDERINGND,&isrow,&iscol));
66*9566063dSJacob Faibussowitsch   PetscCall(MatPermute(mat,isrow,iscol,&B));
67*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by ND\n"));
68*9566063dSJacob Faibussowitsch   PetscCall(MatView(B,viewer));
69*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
70*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"ND row permutation\n"));
71*9566063dSJacob Faibussowitsch   PetscCall(ISView(isrow,viewer));
72*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"ND column permutation\n"));
73*9566063dSJacob Faibussowitsch   PetscCall(ISView(iscol,viewer));
74c4762a1bSJed Brown 
75*9566063dSJacob Faibussowitsch   PetscCall(MatReorderForNonzeroDiagonal(mat,1.e-8,isrow,iscol));
76*9566063dSJacob Faibussowitsch   PetscCall(MatPermute(mat,isrow,iscol,&B));
77*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by ND + NonzeroDiagonal()\n"));
78*9566063dSJacob Faibussowitsch   PetscCall(MatView(B,viewer));
79*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
80*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"ND + NonzeroDiagonal() row permutation\n"));
81*9566063dSJacob Faibussowitsch   PetscCall(ISView(isrow,viewer));
82*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"ND + NonzeroDiagonal() column permutation\n"));
83*9566063dSJacob Faibussowitsch   PetscCall(ISView(iscol,viewer));
84c4762a1bSJed Brown 
85*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
86*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
87c4762a1bSJed Brown 
88*9566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(mat,MATORDERINGRCM,&isrow,&iscol));
89*9566063dSJacob Faibussowitsch   PetscCall(MatPermute(mat,isrow,iscol,&B));
90*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by RCM\n"));
91*9566063dSJacob Faibussowitsch   PetscCall(MatView(B,viewer));
92*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
93*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"RCM row permutation\n"));
94*9566063dSJacob Faibussowitsch   PetscCall(ISView(isrow,viewer));
95*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"RCM column permutation\n"));
96*9566063dSJacob Faibussowitsch   PetscCall(ISView(iscol,viewer));
97c4762a1bSJed Brown 
98*9566063dSJacob Faibussowitsch   PetscCall(MatReorderForNonzeroDiagonal(mat,1.e-8,isrow,iscol));
99*9566063dSJacob Faibussowitsch   PetscCall(MatPermute(mat,isrow,iscol,&B));
100*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix permuted by RCM + NonzeroDiagonal()\n"));
101*9566063dSJacob Faibussowitsch   PetscCall(MatView(B,viewer));
102*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"RCM + NonzeroDiagonal() row permutation\n"));
103*9566063dSJacob Faibussowitsch   PetscCall(ISView(isrow,viewer));
104*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"RCM + NonzeroDiagonal() column permutation\n"));
105*9566063dSJacob Faibussowitsch   PetscCall(ISView(iscol,viewer));
106*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
107c4762a1bSJed Brown   if (size == 1) {
108*9566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE));
109*9566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,4,0,1,&identity));
110*9566063dSJacob Faibussowitsch     PetscCall(MatPermute(B,identity,identity,&C));
111*9566063dSJacob Faibussowitsch     PetscCall(MatConvert(C,MATSEQSBAIJ,MAT_INPLACE_MATRIX,&C));
112*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&C));
113*9566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&identity));
114c4762a1bSJed Brown   }
115*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
116c4762a1bSJed Brown   /* Test MatLUFactor(); set diagonal as zeros as requested by PETSc matrix factorization */
117c4762a1bSJed Brown   for (i=0; i<4; i++) {
118c4762a1bSJed Brown     v = 0.0;
119*9566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat,1,&i,1,&i,&v,INSERT_VALUES));
120c4762a1bSJed Brown   }
121*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
122*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
123*9566063dSJacob Faibussowitsch   PetscCall(MatLUFactor(mat,isrow,iscol,NULL));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* Free data structures */
126*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
127*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
128*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
129c4762a1bSJed Brown 
130*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
131b122ec5aSJacob Faibussowitsch   return 0;
132c4762a1bSJed Brown }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown /*TEST
135c4762a1bSJed Brown 
136c4762a1bSJed Brown    test:
137c4762a1bSJed Brown 
138c4762a1bSJed Brown TEST*/
139