xref: /petsc/src/mat/tests/ex30.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1 
2 static char help[] = "Tests ILU and ICC factorization with and without matrix ordering on seqaij format, and illustrates drawing of matrix sparsity structure with MatView().\n\
3   Input parameters are:\n\
4   -lf <level> : level of fill for ILU (default is 0)\n\
5   -lu : use full LU or Cholesky factorization\n\
6   -m <value>,-n <value> : grid dimensions\n\
7 Note that most users should employ the KSP interface to the\n\
8 linear solvers instead of using the factorization routines\n\
9 directly.\n\n";
10 
11 #include <petscmat.h>
12 
13 int main(int argc,char **args)
14 {
15   Mat            C,A;
16   PetscInt       i,j,m = 5,n = 5,Ii,J,lf = 0;
17   PetscErrorCode ierr;
18   PetscBool      LU=PETSC_FALSE,CHOLESKY,TRIANGULAR=PETSC_FALSE,MATDSPL=PETSC_FALSE,flg,matordering;
19   PetscScalar    v;
20   IS             row,col;
21   PetscViewer    viewer1,viewer2;
22   MatFactorInfo  info;
23   Vec            x,y,b,ytmp;
24   PetscReal      norm2,norm2_inplace, tol = 100.*PETSC_MACHINE_EPSILON;
25   PetscRandom    rdm;
26   PetscMPIInt    size;
27 
28   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
29   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
30   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
31   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
32   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
33   ierr = PetscOptionsGetInt(NULL,NULL,"-lf",&lf,NULL);CHKERRQ(ierr);
34 
35   ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,400,400,&viewer1);CHKERRQ(ierr);
36   ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,400,0,400,400,&viewer2);CHKERRQ(ierr);
37 
38   ierr = MatCreate(PETSC_COMM_SELF,&C);CHKERRQ(ierr);
39   ierr = MatSetSizes(C,m*n,m*n,m*n,m*n);CHKERRQ(ierr);
40   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
41   ierr = MatSetUp(C);CHKERRQ(ierr);
42 
43   /* Create matrix C in seqaij format and sC in seqsbaij. (This is five-point stencil with some extra elements) */
44   for (i=0; i<m; i++) {
45     for (j=0; j<n; j++) {
46       v = -1.0;  Ii = j + n*i;
47       J = Ii - n; if (J>=0)  {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
48       J = Ii + n; if (J<m*n) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
49       J = Ii - 1; if (J>=0)  {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
50       J = Ii + 1; if (J<m*n) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
51       v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
52     }
53   }
54   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
55   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
56 
57   ierr = MatIsSymmetric(C,0.0,&flg);CHKERRQ(ierr);
58   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"C is non-symmetric");
59 
60   /* Create vectors for error checking */
61   ierr = MatCreateVecs(C,&x,&b);CHKERRQ(ierr);
62   ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
63   ierr = VecDuplicate(x,&ytmp);CHKERRQ(ierr);
64   ierr = PetscRandomCreate(PETSC_COMM_SELF,&rdm);CHKERRQ(ierr);
65   ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr);
66   ierr = VecSetRandom(x,rdm);CHKERRQ(ierr);
67   ierr = MatMult(C,x,b);CHKERRQ(ierr);
68 
69   ierr = PetscOptionsHasName(NULL,NULL,"-mat_ordering",&matordering);CHKERRQ(ierr);
70   if (matordering) {
71     ierr = MatGetOrdering(C,MATORDERINGRCM,&row,&col);CHKERRQ(ierr);
72   } else {
73     ierr = MatGetOrdering(C,MATORDERINGNATURAL,&row,&col);CHKERRQ(ierr);
74   }
75 
76   ierr = PetscOptionsHasName(NULL,NULL,"-display_matrices",&MATDSPL);CHKERRQ(ierr);
77   if (MATDSPL) {
78     printf("original matrix:\n");
79     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
80     ierr = MatView(C,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
81     ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
82     ierr = MatView(C,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
83     ierr = MatView(C,viewer1);CHKERRQ(ierr);
84   }
85 
86   /* Compute LU or ILU factor A */
87   ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr);
88 
89   info.fill          = 1.0;
90   info.diagonal_fill = 0;
91   info.zeropivot     = 0.0;
92 
93   ierr = PetscOptionsHasName(NULL,NULL,"-lu",&LU);CHKERRQ(ierr);
94   if (LU) {
95     printf("Test LU...\n");
96     ierr = MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&A);CHKERRQ(ierr);
97     ierr = MatLUFactorSymbolic(A,C,row,col,&info);CHKERRQ(ierr);
98   } else {
99     printf("Test ILU...\n");
100     info.levels = lf;
101 
102     ierr = MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_ILU,&A);CHKERRQ(ierr);
103     ierr = MatILUFactorSymbolic(A,C,row,col,&info);CHKERRQ(ierr);
104   }
105   ierr = MatLUFactorNumeric(A,C,&info);CHKERRQ(ierr);
106 
107   /* Solve A*y = b, then check the error */
108   ierr = MatSolve(A,b,y);CHKERRQ(ierr);
109   ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr);
110   ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr);
111   ierr = MatDestroy(&A);CHKERRQ(ierr);
112 
113   /* Test in-place ILU(0) and compare it with the out-place ILU(0) */
114   if (!LU && lf==0) {
115     ierr = MatDuplicate(C,MAT_COPY_VALUES,&A);CHKERRQ(ierr);
116     ierr = MatILUFactor(A,row,col,&info);CHKERRQ(ierr);
117     /*
118     printf("In-place factored matrix:\n");
119     ierr = MatView(C,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
120     */
121     ierr = MatSolve(A,b,y);CHKERRQ(ierr);
122     ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr);
123     ierr = VecNorm(y,NORM_2,&norm2_inplace);CHKERRQ(ierr);
124     if (PetscAbs(norm2 - norm2_inplace) > tol) SETERRQ2(PETSC_COMM_SELF,1,"ILU(0) %g and in-place ILU(0) %g give different residuals",(double)norm2,(double)norm2_inplace);
125     ierr = MatDestroy(&A);CHKERRQ(ierr);
126   }
127 
128   /* Test Cholesky and ICC on seqaij matrix with matrix reordering on aij matrix C */
129   CHOLESKY = LU;
130   if (CHOLESKY) {
131     printf("Test Cholesky...\n");
132     lf   = -1;
133     ierr = MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&A);CHKERRQ(ierr);
134     ierr = MatCholeskyFactorSymbolic(A,C,row,&info);CHKERRQ(ierr);
135   } else {
136     printf("Test ICC...\n");
137     info.levels        = lf;
138     info.fill          = 1.0;
139     info.diagonal_fill = 0;
140     info.zeropivot     = 0.0;
141 
142     ierr = MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_ICC,&A);CHKERRQ(ierr);
143     ierr = MatICCFactorSymbolic(A,C,row,&info);CHKERRQ(ierr);
144   }
145   ierr = MatCholeskyFactorNumeric(A,C,&info);CHKERRQ(ierr);
146 
147   /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering on aij matrix C */
148   if (lf == -1) {
149     ierr = PetscOptionsHasName(NULL,NULL,"-triangular_solve",&TRIANGULAR);CHKERRQ(ierr);
150     if (TRIANGULAR) {
151       printf("Test MatForwardSolve...\n");
152       ierr = MatForwardSolve(A,b,ytmp);CHKERRQ(ierr);
153       printf("Test MatBackwardSolve...\n");
154       ierr = MatBackwardSolve(A,ytmp,y);CHKERRQ(ierr);
155       ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr);
156       ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr);
157       if (norm2 > tol) {
158         ierr = PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g\n",(double)norm2);CHKERRQ(ierr);
159       }
160     }
161   }
162 
163   ierr = MatSolve(A,b,y);CHKERRQ(ierr);
164   ierr = MatDestroy(&A);CHKERRQ(ierr);
165   ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr);
166   ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr);
167   if (lf == -1 && norm2 > tol) {
168     PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ:   Cholesky/ICC levels %d, residual %g\n",lf,norm2);CHKERRQ(ierr);
169   }
170 
171   /* Test in-place ICC(0) and compare it with the out-place ICC(0) */
172   if (!CHOLESKY && lf==0 && !matordering) {
173     ierr = MatConvert(C,MATSBAIJ,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
174     ierr = MatICCFactor(A,row,&info);CHKERRQ(ierr);
175     /*
176     printf("In-place factored matrix:\n");
177     ierr = MatView(A,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
178     */
179     ierr = MatSolve(A,b,y);CHKERRQ(ierr);
180     ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr);
181     ierr = VecNorm(y,NORM_2,&norm2_inplace);CHKERRQ(ierr);
182     if (PetscAbs(norm2 - norm2_inplace) > tol) SETERRQ2(PETSC_COMM_SELF,1,"ICC(0) %g and in-place ICC(0) %g give different residuals",(double)norm2,(double)norm2_inplace);
183     ierr = MatDestroy(&A);CHKERRQ(ierr);
184   }
185 
186   /* Free data structures */
187   ierr = ISDestroy(&row);CHKERRQ(ierr);
188   ierr = ISDestroy(&col);CHKERRQ(ierr);
189   ierr = MatDestroy(&C);CHKERRQ(ierr);
190   ierr = PetscViewerDestroy(&viewer1);CHKERRQ(ierr);
191   ierr = PetscViewerDestroy(&viewer2);CHKERRQ(ierr);
192   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
193   ierr = VecDestroy(&x);CHKERRQ(ierr);
194   ierr = VecDestroy(&y);CHKERRQ(ierr);
195   ierr = VecDestroy(&ytmp);CHKERRQ(ierr);
196   ierr = VecDestroy(&b);CHKERRQ(ierr);
197   ierr = PetscFinalize();
198   return ierr;
199 }
200 
201 
202 /*TEST
203 
204    test:
205       args: -mat_ordering -display_matrices -nox
206       filter: grep -v "MPI processes"
207 
208    test:
209       suffix: 2
210       args: -mat_ordering -display_matrices -nox -lu
211 
212    test:
213       suffix: 3
214       args: -mat_ordering -lu -triangular_solve
215 
216    test:
217       suffix: 4
218 
219    test:
220       suffix: 5
221       args: -lu
222 
223    test:
224       suffix: 6
225       args: -lu -triangular_solve
226       output_file: output/ex30_3.out
227 
228 TEST*/
229