xref: /petsc/src/mat/tests/ex7.c (revision d0609ced746bc51b019815ca91d747429db24893)
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;
121a6d72e3SSatish Balay   PetscInt       i,j,m,n,Ii,J;
13c4762a1bSJed Brown   PetscScalar    v,one = 1.0;
14c4762a1bSJed Brown   IS             perm,iperm;
15c4762a1bSJed Brown   Vec            x,u,b;
161a6d72e3SSatish Balay   PetscReal      norm,fill;
17c4762a1bSJed Brown   MatFactorInfo  luinfo;
18c4762a1bSJed Brown 
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
201a6d72e3SSatish Balay 
21*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Mat test ex7 options","Mat");
221a6d72e3SSatish Balay   m = 3; n = 3; fill = 2.0;
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-m","Number of rows in grid",NULL,m,&m,NULL));
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-n","Number of columns in grid",NULL,n,&n,NULL));
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-fill","Expected fill ratio for factorization",NULL,fill,&fill,NULL));
261a6d72e3SSatish Balay 
27*d0609cedSBarry Smith   PetscOptionsEnd();
281a6d72e3SSatish Balay 
29c4762a1bSJed Brown   /* Create the matrix for the five point stencil, YET AGAIN */
309566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&C));
319566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
329566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
339566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
34c4762a1bSJed Brown   for (i=0; i<m; i++) {
35c4762a1bSJed Brown     for (j=0; j<n; j++) {
36c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
379566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
389566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
399566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
409566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
419566063dSJacob Faibussowitsch       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
42c4762a1bSJed Brown     }
43c4762a1bSJed Brown   }
449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
469566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(C,MATORDERINGRCM,&perm,&iperm));
479566063dSJacob Faibussowitsch   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
489566063dSJacob Faibussowitsch   PetscCall(ISView(perm,PETSC_VIEWER_STDOUT_SELF));
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&luinfo));
51c4762a1bSJed Brown 
521a6d72e3SSatish Balay   luinfo.fill          = fill;
53c4762a1bSJed Brown   luinfo.dtcol         = 0.0;
54c4762a1bSJed Brown   luinfo.zeropivot     = 1.e-14;
55c4762a1bSJed Brown   luinfo.pivotinblocks = 1.0;
56c4762a1bSJed Brown 
579566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&LU));
589566063dSJacob Faibussowitsch   PetscCall(MatLUFactorSymbolic(LU,C,perm,iperm,&luinfo));
599566063dSJacob Faibussowitsch   PetscCall(MatLUFactorNumeric(LU,C,&luinfo));
60c4762a1bSJed Brown 
619566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,m*n,&u));
629566063dSJacob Faibussowitsch   PetscCall(VecSet(u,one));
639566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&x));
649566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&b));
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch   PetscCall(MatMult(C,u,b));
679566063dSJacob Faibussowitsch   PetscCall(MatSolve(LU,b,x));
68c4762a1bSJed Brown 
699566063dSJacob Faibussowitsch   PetscCall(VecView(b,PETSC_VIEWER_STDOUT_SELF));
709566063dSJacob Faibussowitsch   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_SELF));
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x,-1.0,u));
739566063dSJacob Faibussowitsch   PetscCall(VecNorm(x,NORM_2,&norm));
749566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",(double)norm));
75c4762a1bSJed Brown 
769566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(C,MAT_LOCAL,&info));
779566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"original matrix nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used));
789566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(LU,MAT_LOCAL,&info));
799566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"factored matrix nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used));
80c4762a1bSJed Brown 
819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
849566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
859566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iperm));
869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
879566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&LU));
88c4762a1bSJed Brown 
899566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
90b122ec5aSJacob Faibussowitsch   return 0;
91c4762a1bSJed Brown }
92c4762a1bSJed Brown 
93c4762a1bSJed Brown /*TEST
94c4762a1bSJed Brown 
95c4762a1bSJed Brown    test:
961a6d72e3SSatish Balay       suffix: 1
971a6d72e3SSatish Balay       filter: grep -v "MPI processes"
981a6d72e3SSatish Balay 
991a6d72e3SSatish Balay    test:
1001a6d72e3SSatish Balay       suffix: 2
1011a6d72e3SSatish Balay       args: -m 1 -n 1 -fill 0.49
102c4762a1bSJed Brown       filter: grep -v "MPI processes"
103c4762a1bSJed Brown 
104c4762a1bSJed Brown TEST*/
105