xref: /petsc/src/mat/tests/ex7.c (revision 4e278199b78715991f5c71ebbd945c1489263e6c)
1 
2 static char help[] = "Tests matrix factorization.  Note that most users should\n\
3 employ the KSP  interface to the linear solvers instead of using the factorization\n\
4 routines directly.\n\n";
5 
6 #include <petscmat.h>
7 
8 int main(int argc,char **args)
9 {
10   Mat            C,LU;
11   MatInfo        info;
12   PetscInt       i,j,m,n,Ii,J;
13   PetscErrorCode ierr;
14   PetscScalar    v,one = 1.0;
15   IS             perm,iperm;
16   Vec            x,u,b;
17   PetscReal      norm,fill;
18   MatFactorInfo  luinfo;
19 
20   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
21 
22   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Mat test ex7 options","Mat");CHKERRQ(ierr);
23   m = 3; n = 3; fill = 2.0;
24   ierr = PetscOptionsInt("-m","Number of rows in grid",NULL,m,&m,NULL);CHKERRQ(ierr);
25   ierr = PetscOptionsInt("-n","Number of columns in grid",NULL,n,&n,NULL);CHKERRQ(ierr);
26   ierr = PetscOptionsReal("-fill","Expected fill ratio for factorization",NULL,fill,&fill,NULL);CHKERRQ(ierr);
27 
28   ierr = PetscOptionsEnd();CHKERRQ(ierr);
29 
30   /* Create the matrix for the five point stencil, YET AGAIN */
31   ierr = MatCreate(PETSC_COMM_SELF,&C);CHKERRQ(ierr);
32   ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
33   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
34   ierr = MatSetUp(C);CHKERRQ(ierr);
35   for (i=0; i<m; i++) {
36     for (j=0; j<n; j++) {
37       v = -1.0;  Ii = j + n*i;
38       if (i>0)   {J = Ii - n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
39       if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
40       if (j>0)   {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
41       if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
42       v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
43     }
44   }
45   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
46   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
47   ierr = MatGetOrdering(C,MATORDERINGRCM,&perm,&iperm);CHKERRQ(ierr);
48   ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
49   ierr = ISView(perm,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
50 
51   ierr = MatFactorInfoInitialize(&luinfo);CHKERRQ(ierr);
52 
53   luinfo.fill          = fill;
54   luinfo.dtcol         = 0.0;
55   luinfo.zeropivot     = 1.e-14;
56   luinfo.pivotinblocks = 1.0;
57 
58   ierr = MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&LU);CHKERRQ(ierr);
59   ierr = MatLUFactorSymbolic(LU,C,perm,iperm,&luinfo);CHKERRQ(ierr);
60   ierr = MatLUFactorNumeric(LU,C,&luinfo);CHKERRQ(ierr);
61 
62   ierr = VecCreateSeq(PETSC_COMM_SELF,m*n,&u);CHKERRQ(ierr);
63   ierr = VecSet(u,one);CHKERRQ(ierr);
64   ierr = VecDuplicate(u,&x);CHKERRQ(ierr);
65   ierr = VecDuplicate(u,&b);CHKERRQ(ierr);
66 
67   ierr = MatMult(C,u,b);CHKERRQ(ierr);
68   ierr = MatSolve(LU,b,x);CHKERRQ(ierr);
69 
70   ierr = VecView(b,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
71   ierr = VecView(x,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
72 
73   ierr = VecAXPY(x,-1.0,u);CHKERRQ(ierr);
74   ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr);
75   ierr = PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",(double)norm);CHKERRQ(ierr);
76 
77   ierr = MatGetInfo(C,MAT_LOCAL,&info);CHKERRQ(ierr);
78   ierr = PetscPrintf(PETSC_COMM_SELF,"original matrix nonzeros = %D\n",(PetscInt)info.nz_used);CHKERRQ(ierr);
79   ierr = MatGetInfo(LU,MAT_LOCAL,&info);CHKERRQ(ierr);
80   ierr = PetscPrintf(PETSC_COMM_SELF,"factored matrix nonzeros = %D\n",(PetscInt)info.nz_used);CHKERRQ(ierr);
81 
82   ierr = VecDestroy(&u);CHKERRQ(ierr);
83   ierr = VecDestroy(&b);CHKERRQ(ierr);
84   ierr = VecDestroy(&x);CHKERRQ(ierr);
85   ierr = ISDestroy(&perm);CHKERRQ(ierr);
86   ierr = ISDestroy(&iperm);CHKERRQ(ierr);
87   ierr = MatDestroy(&C);CHKERRQ(ierr);
88   ierr = MatDestroy(&LU);CHKERRQ(ierr);
89 
90   ierr = PetscFinalize();
91   return ierr;
92 }
93 
94 /*TEST
95 
96    test:
97       suffix: 1
98       filter: grep -v "MPI processes"
99 
100    test:
101       suffix: 2
102       args: -m 1 -n 1 -fill 0.49
103       filter: grep -v "MPI processes"
104 
105 TEST*/
106