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