xref: /petsc/src/mat/tests/ex30.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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;
295f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
302c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-lf",&lf,NULL));
34c4762a1bSJed Brown 
355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,400,400,&viewer1));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,400,0,400,400,&viewer2));
37c4762a1bSJed Brown 
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&C));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,m*n,m*n,m*n,m*n));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
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;
475f80ce2aSJacob Faibussowitsch       J = Ii - n; if (J>=0)  CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
485f80ce2aSJacob Faibussowitsch       J = Ii + n; if (J<m*n) CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
495f80ce2aSJacob Faibussowitsch       J = Ii - 1; if (J>=0)  CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
505f80ce2aSJacob Faibussowitsch       J = Ii + 1; if (J<m*n) CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
515f80ce2aSJacob Faibussowitsch       v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
52c4762a1bSJed Brown     }
53c4762a1bSJed Brown   }
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
56c4762a1bSJed Brown 
575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatIsSymmetric(C,0.0,&flg));
58*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"C is non-symmetric");
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Create vectors for error checking */
615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(C,&x,&b));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&y));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&ytmp));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rdm));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(x,rdm));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(C,x,b));
68c4762a1bSJed Brown 
695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mat_ordering",&matordering));
70c4762a1bSJed Brown   if (matordering) {
715f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetOrdering(C,MATORDERINGRCM,&row,&col));
72c4762a1bSJed Brown   } else {
735f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetOrdering(C,MATORDERINGNATURAL,&row,&col));
74c4762a1bSJed Brown   }
75c4762a1bSJed Brown 
765f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-display_matrices",&MATDSPL));
77c4762a1bSJed Brown   if (MATDSPL) {
78c4762a1bSJed Brown     printf("original matrix:\n");
795f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_INFO));
805f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_SELF));
815f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF));
825f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_SELF));
835f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(C,viewer1));
84c4762a1bSJed Brown   }
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* Compute LU or ILU factor A */
875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatFactorInfoInitialize(&info));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   info.fill          = 1.0;
90c4762a1bSJed Brown   info.diagonal_fill = 0;
91c4762a1bSJed Brown   info.zeropivot     = 0.0;
92c4762a1bSJed Brown 
935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-lu",&LU));
94c4762a1bSJed Brown   if (LU) {
95c4762a1bSJed Brown     printf("Test LU...\n");
965f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&A));
975f80ce2aSJacob Faibussowitsch     CHKERRQ(MatLUFactorSymbolic(A,C,row,col,&info));
98c4762a1bSJed Brown   } else {
99c4762a1bSJed Brown     printf("Test ILU...\n");
100c4762a1bSJed Brown     info.levels = lf;
101c4762a1bSJed Brown 
1025f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_ILU,&A));
1035f80ce2aSJacob Faibussowitsch     CHKERRQ(MatILUFactorSymbolic(A,C,row,col,&info));
104c4762a1bSJed Brown   }
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactorNumeric(A,C,&info));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* Solve A*y = b, then check the error */
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSolve(A,b,y));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(y,-1.0,x));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(y,NORM_2,&norm2));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   /* Test in-place ILU(0) and compare it with the out-place ILU(0) */
114c4762a1bSJed Brown   if (!LU && lf==0) {
1155f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&A));
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(MatILUFactor(A,row,col,&info));
117c4762a1bSJed Brown     /*
118c4762a1bSJed Brown     printf("In-place factored matrix:\n");
1195f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_SELF));
120c4762a1bSJed Brown     */
1215f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolve(A,b,y));
1225f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(y,-1.0,x));
1235f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(y,NORM_2,&norm2_inplace));
1242c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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);
1255f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&A));
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;
1335f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&A));
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCholeskyFactorSymbolic(A,C,row,&info));
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 
1425f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_ICC,&A));
1435f80ce2aSJacob Faibussowitsch     CHKERRQ(MatICCFactorSymbolic(A,C,row,&info));
144c4762a1bSJed Brown   }
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCholeskyFactorNumeric(A,C,&info));
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering on aij matrix C */
148c4762a1bSJed Brown   if (lf == -1) {
1495f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsHasName(NULL,NULL,"-triangular_solve",&TRIANGULAR));
150c4762a1bSJed Brown     if (TRIANGULAR) {
151c4762a1bSJed Brown       printf("Test MatForwardSolve...\n");
1525f80ce2aSJacob Faibussowitsch       CHKERRQ(MatForwardSolve(A,b,ytmp));
153c4762a1bSJed Brown       printf("Test MatBackwardSolve...\n");
1545f80ce2aSJacob Faibussowitsch       CHKERRQ(MatBackwardSolve(A,ytmp,y));
1555f80ce2aSJacob Faibussowitsch       CHKERRQ(VecAXPY(y,-1.0,x));
1565f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(y,NORM_2,&norm2));
157c4762a1bSJed Brown       if (norm2 > tol) {
1585f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g\n",(double)norm2));
159c4762a1bSJed Brown       }
160c4762a1bSJed Brown     }
161c4762a1bSJed Brown   }
162c4762a1bSJed Brown 
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSolve(A,b,y));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(y,-1.0,x));
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(y,NORM_2,&norm2));
167c4762a1bSJed Brown   if (lf == -1 && norm2 > tol) {
1685f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ:   Cholesky/ICC levels %" PetscInt_FMT ", residual %g\n",lf,(double)norm2));
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) {
1735f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(C,MATSBAIJ,MAT_INITIAL_MATRIX,&A));
1745f80ce2aSJacob Faibussowitsch     CHKERRQ(MatICCFactor(A,row,&info));
175c4762a1bSJed Brown     /*
176c4762a1bSJed Brown     printf("In-place factored matrix:\n");
1775f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_SELF));
178c4762a1bSJed Brown     */
1795f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolve(A,b,y));
1805f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(y,-1.0,x));
1815f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(y,NORM_2,&norm2_inplace));
1822c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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);
1835f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&A));
184c4762a1bSJed Brown   }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* Free data structures */
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&row));
1885f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&col));
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer1));
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer2));
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
1935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ytmp));
1965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
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