Actual source code: ex30.c
1: /*$Id: ex30.c,v 1.27 2001/08/07 03:03:07 balay Exp $*/
3: static char help[] = "Tests ILU factorization and illustrates drawing of matrix sparsity structure with MatView().n
4: Input parameters are:n
5: -lf <level> : level of fill for ILU (default is 0)n
6: -lu : use full LU factorizationn
7: -m <value>,-n <value> : grid dimensionsn
8: Note that most users should employ the SLES interface to then
9: linear solvers instead of using the factorization routinesn
10: directly.nn";
12: #include petscmat.h
14: int main(int argc,char **args)
15: {
16: Mat C,A;
17: int i,j,m = 5,n = 5,I,J,ierr,lf = 0;
18: PetscTruth flg1;
19: PetscScalar v;
20: IS row,col;
21: PetscViewer viewer1,viewer2;
23: PetscInitialize(&argc,&args,(char *)0,help);
24: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
25: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
26: PetscOptionsGetInt(PETSC_NULL,"-lf",&lf,PETSC_NULL);
28: PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,400,400,&viewer1);
29: PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,400,0,400,400,&viewer2);
31: MatCreate(PETSC_COMM_SELF,m*n,m*n,m*n,m*n,&C);
32: MatSetFromOptions(C);
33: MatSeqBDiagSetPreallocation(C,0,1,PETSC_NULL,PETSC_NULL);
34: MatSeqAIJSetPreallocation(C,5,PETSC_NULL);
36: /* Create the matrix. (This is five-point stencil with some extra elements) */
37: for (i=0; i<m; i++) {
38: for (j=0; j<n; j++) {
39: v = -1.0; I = j + n*i;
40: J = I - n; if (J>=0) {MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
41: J = I + n; if (J<m*n) {MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
42: J = I - 1; if (J>=0) {MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
43: J = I + 1; if (J<m*n) {MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
44: v = 4.0; MatSetValues(C,1,&I,1,&I,&v,INSERT_VALUES);
45: }
46: }
47: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
48: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
50: MatGetOrdering(C,MATORDERING_RCM,&row,&col);
51: printf("original matrix:n");
52: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_INFO);
53: MatView(C,PETSC_VIEWER_STDOUT_SELF);
54: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);
55: MatView(C,PETSC_VIEWER_STDOUT_SELF);
56: MatView(C,viewer1);
58: /* Compute factorization */
59: PetscOptionsHasName(PETSC_NULL,"-lu",&flg1);
60: if (flg1){
61: MatLUFactorSymbolic(C,row,col,PETSC_NULL,&A);
62: } else {
63: MatILUInfo info;
64: info.levels = lf;
65: info.fill = 1.0;
66: info.diagonal_fill = 0;
67: info.damping = 0;
68: info.damp = 0;
69: info.zeropivot = 0.0;
70: MatILUFactorSymbolic(C,row,col,&info,&A);
71: }
72: MatLUFactorNumeric(C,&A);
74: printf("factored matrix:n");
75: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_INFO);
76: MatView(A,PETSC_VIEWER_STDOUT_SELF);
77: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);
78: MatView(A,PETSC_VIEWER_STDOUT_SELF);
79: MatView(A,viewer2);
81: /* Free data structures */
82: MatDestroy(C);
83: MatDestroy(A);
84: ISDestroy(row);
85: ISDestroy(col);
86: PetscViewerDestroy(viewer1);
87: PetscViewerDestroy(viewer2);
88: PetscFinalize();
89: return 0;
90: }