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