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