xref: /petsc/src/mat/tests/ex7.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
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   CHKERRQ(PetscOptionsInt("-m","Number of rows in grid",NULL,m,&m,NULL));
25   CHKERRQ(PetscOptionsInt("-n","Number of columns in grid",NULL,n,&n,NULL));
26   CHKERRQ(PetscOptionsReal("-fill","Expected fill ratio for factorization",NULL,fill,&fill,NULL));
27 
28   ierr = PetscOptionsEnd();CHKERRQ(ierr);
29 
30   /* Create the matrix for the five point stencil, YET AGAIN */
31   CHKERRQ(MatCreate(PETSC_COMM_SELF,&C));
32   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
33   CHKERRQ(MatSetFromOptions(C));
34   CHKERRQ(MatSetUp(C));
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; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
39       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
40       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
41       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
42       v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
43     }
44   }
45   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
46   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
47   CHKERRQ(MatGetOrdering(C,MATORDERINGRCM,&perm,&iperm));
48   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
49   CHKERRQ(ISView(perm,PETSC_VIEWER_STDOUT_SELF));
50 
51   CHKERRQ(MatFactorInfoInitialize(&luinfo));
52 
53   luinfo.fill          = fill;
54   luinfo.dtcol         = 0.0;
55   luinfo.zeropivot     = 1.e-14;
56   luinfo.pivotinblocks = 1.0;
57 
58   CHKERRQ(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&LU));
59   CHKERRQ(MatLUFactorSymbolic(LU,C,perm,iperm,&luinfo));
60   CHKERRQ(MatLUFactorNumeric(LU,C,&luinfo));
61 
62   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m*n,&u));
63   CHKERRQ(VecSet(u,one));
64   CHKERRQ(VecDuplicate(u,&x));
65   CHKERRQ(VecDuplicate(u,&b));
66 
67   CHKERRQ(MatMult(C,u,b));
68   CHKERRQ(MatSolve(LU,b,x));
69 
70   CHKERRQ(VecView(b,PETSC_VIEWER_STDOUT_SELF));
71   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_SELF));
72 
73   CHKERRQ(VecAXPY(x,-1.0,u));
74   CHKERRQ(VecNorm(x,NORM_2,&norm));
75   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",(double)norm));
76 
77   CHKERRQ(MatGetInfo(C,MAT_LOCAL,&info));
78   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"original matrix nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used));
79   CHKERRQ(MatGetInfo(LU,MAT_LOCAL,&info));
80   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"factored matrix nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used));
81 
82   CHKERRQ(VecDestroy(&u));
83   CHKERRQ(VecDestroy(&b));
84   CHKERRQ(VecDestroy(&x));
85   CHKERRQ(ISDestroy(&perm));
86   CHKERRQ(ISDestroy(&iperm));
87   CHKERRQ(MatDestroy(&C));
88   CHKERRQ(MatDestroy(&LU));
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