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