xref: /petsc/src/mat/tests/ex128.c (revision daa037dfd3c3bec8dc8659548d2b20b07c1dc6de)
1 
2 static char help[] = "Tests ILU and ICC factorization with and without matrix ordering on seqsbaij format. Modified from ex30.c\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,sC,sA;
16   PetscInt       i,j,m = 5,n = 5,Ii,J,lf = 0;
17   PetscBool      CHOLESKY=PETSC_FALSE,TRIANGULAR=PETSC_FALSE,flg;
18   PetscScalar    v;
19   IS             row,col;
20   MatFactorInfo  info;
21   Vec            x,y,b,ytmp;
22   PetscReal      norm2,tol = 100*PETSC_MACHINE_EPSILON;
23   PetscRandom    rdm;
24   PetscMPIInt    size;
25 
26   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
27   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
28   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
29   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
30   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
31   PetscCall(PetscOptionsGetInt(NULL,NULL,"-lf",&lf,NULL));
32 
33   PetscCall(MatCreate(PETSC_COMM_SELF,&C));
34   PetscCall(MatSetSizes(C,m*n,m*n,m*n,m*n));
35   PetscCall(MatSetFromOptions(C));
36   PetscCall(MatSetUp(C));
37 
38   /* Create matrix C in seqaij format and sC in seqsbaij. (This is five-point stencil with some extra elements) */
39   for (i=0; i<m; i++) {
40     for (j=0; j<n; j++) {
41       v = -1.0;  Ii = j + n*i;
42       J = Ii - n; if (J>=0)  PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
43       J = Ii + n; if (J<m*n) PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
44       J = Ii - 1; if (J>=0)  PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
45       J = Ii + 1; if (J<m*n) PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
46       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
47     }
48   }
49   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
50   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
51 
52   PetscCall(MatIsSymmetric(C,0.0,&flg));
53   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"C is non-symmetric");
54   PetscCall(MatConvert(C,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sC));
55 
56   /* Create vectors for error checking */
57   PetscCall(MatCreateVecs(C,&x,&b));
58   PetscCall(VecDuplicate(x,&y));
59   PetscCall(VecDuplicate(x,&ytmp));
60   PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm));
61   PetscCall(PetscRandomSetFromOptions(rdm));
62   PetscCall(VecSetRandom(x,rdm));
63   PetscCall(MatMult(C,x,b));
64 
65   PetscCall(MatGetOrdering(C,MATORDERINGNATURAL,&row,&col));
66 
67   /* Compute CHOLESKY or ICC factor sA */
68   PetscCall(MatFactorInfoInitialize(&info));
69 
70   info.fill          = 1.0;
71   info.diagonal_fill = 0;
72   info.zeropivot     = 0.0;
73 
74   PetscCall(PetscOptionsHasName(NULL,NULL,"-cholesky",&CHOLESKY));
75   if (CHOLESKY) {
76     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Test CHOLESKY...\n"));
77     PetscCall(MatGetFactor(sC,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sA));
78     PetscCall(MatCholeskyFactorSymbolic(sA,sC,row,&info));
79   } else {
80     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Test ICC...\n"));
81     info.levels = lf;
82 
83     PetscCall(MatGetFactor(sC,MATSOLVERPETSC,MAT_FACTOR_ICC,&sA));
84     PetscCall(MatICCFactorSymbolic(sA,sC,row,&info));
85   }
86   PetscCall(MatCholeskyFactorNumeric(sA,sC,&info));
87 
88   /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering on aij matrix C */
89   if (CHOLESKY) {
90     PetscCall(PetscOptionsHasName(NULL,NULL,"-triangular_solve",&TRIANGULAR));
91     if (TRIANGULAR) {
92       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Test MatForwardSolve...\n"));
93       PetscCall(MatForwardSolve(sA,b,ytmp));
94       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Test MatBackwardSolve...\n"));
95       PetscCall(MatBackwardSolve(sA,ytmp,y));
96       PetscCall(VecAXPY(y,-1.0,x));
97       PetscCall(VecNorm(y,NORM_2,&norm2));
98       if (norm2 > tol) {
99         PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g\n",(double)norm2));
100       }
101     }
102   }
103 
104   PetscCall(MatSolve(sA,b,y));
105   PetscCall(MatDestroy(&sC));
106   PetscCall(MatDestroy(&sA));
107   PetscCall(VecAXPY(y,-1.0,x));
108   PetscCall(VecNorm(y,NORM_2,&norm2));
109   if (lf == -1 && norm2 > tol) {
110     PetscCall(PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ:   Cholesky/ICC levels %" PetscInt_FMT ", residual %g\n",lf,(double)norm2));
111   }
112 
113   /* Free data structures */
114   PetscCall(MatDestroy(&C));
115   PetscCall(ISDestroy(&row));
116   PetscCall(ISDestroy(&col));
117   PetscCall(PetscRandomDestroy(&rdm));
118   PetscCall(VecDestroy(&x));
119   PetscCall(VecDestroy(&y));
120   PetscCall(VecDestroy(&ytmp));
121   PetscCall(VecDestroy(&b));
122   PetscCall(PetscFinalize());
123   return 0;
124 }
125 
126 /*TEST
127 
128    test:
129       output_file: output/ex128.out
130 
131    test:
132       suffix: 2
133       args: -cholesky -triangular_solve
134 
135 TEST*/
136