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