xref: /petsc/src/mat/tests/ex30.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
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   PetscBool      LU=PETSC_FALSE,CHOLESKY,TRIANGULAR=PETSC_FALSE,MATDSPL=PETSC_FALSE,flg,matordering;
18   PetscScalar    v;
19   IS             row,col;
20   PetscViewer    viewer1,viewer2;
21   MatFactorInfo  info;
22   Vec            x,y,b,ytmp;
23   PetscReal      norm2,norm2_inplace, tol = 100.*PETSC_MACHINE_EPSILON;
24   PetscRandom    rdm;
25   PetscMPIInt    size;
26 
27   PetscFunctionBeginUser;
28   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
29   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
30   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
31   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
32   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
33   PetscCall(PetscOptionsGetInt(NULL,NULL,"-lf",&lf,NULL));
34 
35   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,400,400,&viewer1));
36   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,400,0,400,400,&viewer2));
37 
38   PetscCall(MatCreate(PETSC_COMM_SELF,&C));
39   PetscCall(MatSetSizes(C,m*n,m*n,m*n,m*n));
40   PetscCall(MatSetFromOptions(C));
41   PetscCall(MatSetUp(C));
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)  PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
48       J = Ii + n; if (J<m*n) PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
49       J = Ii - 1; if (J>=0)  PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
50       J = Ii + 1; if (J<m*n) PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
51       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
52     }
53   }
54   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
55   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
56 
57   PetscCall(MatIsSymmetric(C,0.0,&flg));
58   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"C is non-symmetric");
59 
60   /* Create vectors for error checking */
61   PetscCall(MatCreateVecs(C,&x,&b));
62   PetscCall(VecDuplicate(x,&y));
63   PetscCall(VecDuplicate(x,&ytmp));
64   PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm));
65   PetscCall(PetscRandomSetFromOptions(rdm));
66   PetscCall(VecSetRandom(x,rdm));
67   PetscCall(MatMult(C,x,b));
68 
69   PetscCall(PetscOptionsHasName(NULL,NULL,"-mat_ordering",&matordering));
70   if (matordering) {
71     PetscCall(MatGetOrdering(C,MATORDERINGRCM,&row,&col));
72   } else {
73     PetscCall(MatGetOrdering(C,MATORDERINGNATURAL,&row,&col));
74   }
75 
76   PetscCall(PetscOptionsHasName(NULL,NULL,"-display_matrices",&MATDSPL));
77   if (MATDSPL) {
78     printf("original matrix:\n");
79     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_INFO));
80     PetscCall(MatView(C,PETSC_VIEWER_STDOUT_SELF));
81     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF));
82     PetscCall(MatView(C,PETSC_VIEWER_STDOUT_SELF));
83     PetscCall(MatView(C,viewer1));
84   }
85 
86   /* Compute LU or ILU factor A */
87   PetscCall(MatFactorInfoInitialize(&info));
88 
89   info.fill          = 1.0;
90   info.diagonal_fill = 0;
91   info.zeropivot     = 0.0;
92 
93   PetscCall(PetscOptionsHasName(NULL,NULL,"-lu",&LU));
94   if (LU) {
95     printf("Test LU...\n");
96     PetscCall(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&A));
97     PetscCall(MatLUFactorSymbolic(A,C,row,col,&info));
98   } else {
99     printf("Test ILU...\n");
100     info.levels = lf;
101 
102     PetscCall(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_ILU,&A));
103     PetscCall(MatILUFactorSymbolic(A,C,row,col,&info));
104   }
105   PetscCall(MatLUFactorNumeric(A,C,&info));
106 
107   /* Solve A*y = b, then check the error */
108   PetscCall(MatSolve(A,b,y));
109   PetscCall(VecAXPY(y,-1.0,x));
110   PetscCall(VecNorm(y,NORM_2,&norm2));
111   PetscCall(MatDestroy(&A));
112 
113   /* Test in-place ILU(0) and compare it with the out-place ILU(0) */
114   if (!LU && lf==0) {
115     PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&A));
116     PetscCall(MatILUFactor(A,row,col,&info));
117     /*
118     printf("In-place factored matrix:\n");
119     PetscCall(MatView(C,PETSC_VIEWER_STDOUT_SELF));
120     */
121     PetscCall(MatSolve(A,b,y));
122     PetscCall(VecAXPY(y,-1.0,x));
123     PetscCall(VecNorm(y,NORM_2,&norm2_inplace));
124     PetscCheck(PetscAbs(norm2 - norm2_inplace) <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ILU(0) %g and in-place ILU(0) %g give different residuals",(double)norm2,(double)norm2_inplace);
125     PetscCall(MatDestroy(&A));
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     PetscCall(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&A));
134     PetscCall(MatCholeskyFactorSymbolic(A,C,row,&info));
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     PetscCall(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_ICC,&A));
143     PetscCall(MatICCFactorSymbolic(A,C,row,&info));
144   }
145   PetscCall(MatCholeskyFactorNumeric(A,C,&info));
146 
147   /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering on aij matrix C */
148   if (lf == -1) {
149     PetscCall(PetscOptionsHasName(NULL,NULL,"-triangular_solve",&TRIANGULAR));
150     if (TRIANGULAR) {
151       printf("Test MatForwardSolve...\n");
152       PetscCall(MatForwardSolve(A,b,ytmp));
153       printf("Test MatBackwardSolve...\n");
154       PetscCall(MatBackwardSolve(A,ytmp,y));
155       PetscCall(VecAXPY(y,-1.0,x));
156       PetscCall(VecNorm(y,NORM_2,&norm2));
157       if (norm2 > tol) {
158         PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g\n",(double)norm2));
159       }
160     }
161   }
162 
163   PetscCall(MatSolve(A,b,y));
164   PetscCall(MatDestroy(&A));
165   PetscCall(VecAXPY(y,-1.0,x));
166   PetscCall(VecNorm(y,NORM_2,&norm2));
167   if (lf == -1 && norm2 > tol) {
168     PetscCall(PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ:   Cholesky/ICC levels %" PetscInt_FMT ", residual %g\n",lf,(double)norm2));
169   }
170 
171   /* Test in-place ICC(0) and compare it with the out-place ICC(0) */
172   if (!CHOLESKY && lf==0 && !matordering) {
173     PetscCall(MatConvert(C,MATSBAIJ,MAT_INITIAL_MATRIX,&A));
174     PetscCall(MatICCFactor(A,row,&info));
175     /*
176     printf("In-place factored matrix:\n");
177     PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF));
178     */
179     PetscCall(MatSolve(A,b,y));
180     PetscCall(VecAXPY(y,-1.0,x));
181     PetscCall(VecNorm(y,NORM_2,&norm2_inplace));
182     PetscCheck(PetscAbs(norm2 - norm2_inplace) <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ICC(0) %g and in-place ICC(0) %g give different residuals",(double)norm2,(double)norm2_inplace);
183     PetscCall(MatDestroy(&A));
184   }
185 
186   /* Free data structures */
187   PetscCall(ISDestroy(&row));
188   PetscCall(ISDestroy(&col));
189   PetscCall(MatDestroy(&C));
190   PetscCall(PetscViewerDestroy(&viewer1));
191   PetscCall(PetscViewerDestroy(&viewer2));
192   PetscCall(PetscRandomDestroy(&rdm));
193   PetscCall(VecDestroy(&x));
194   PetscCall(VecDestroy(&y));
195   PetscCall(VecDestroy(&ytmp));
196   PetscCall(VecDestroy(&b));
197   PetscCall(PetscFinalize());
198   return 0;
199 }
200 
201 /*TEST
202 
203    test:
204       args: -mat_ordering -display_matrices -nox
205       filter: grep -v " MPI process"
206 
207    test:
208       suffix: 2
209       args: -mat_ordering -display_matrices -nox -lu
210 
211    test:
212       suffix: 3
213       args: -mat_ordering -lu -triangular_solve
214 
215    test:
216       suffix: 4
217 
218    test:
219       suffix: 5
220       args: -lu
221 
222    test:
223       suffix: 6
224       args: -lu -triangular_solve
225       output_file: output/ex30_3.out
226 
227 TEST*/
228