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