xref: /petsc/src/mat/tests/ex7.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
19*327415f7SBarry Smith   PetscFunctionBeginUser;
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
211a6d72e3SSatish Balay 
22d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Mat test ex7 options","Mat");
231a6d72e3SSatish Balay   m = 3; n = 3; fill = 2.0;
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-m","Number of rows in grid",NULL,m,&m,NULL));
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-n","Number of columns in grid",NULL,n,&n,NULL));
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-fill","Expected fill ratio for factorization",NULL,fill,&fill,NULL));
271a6d72e3SSatish Balay 
28d0609cedSBarry Smith   PetscOptionsEnd();
291a6d72e3SSatish Balay 
30c4762a1bSJed Brown   /* Create the matrix for the five point stencil, YET AGAIN */
319566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&C));
329566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
339566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
349566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
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;
389566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
399566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
409566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
419566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
429566063dSJacob Faibussowitsch       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
43c4762a1bSJed Brown     }
44c4762a1bSJed Brown   }
459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
479566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(C,MATORDERINGRCM,&perm,&iperm));
489566063dSJacob Faibussowitsch   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
499566063dSJacob Faibussowitsch   PetscCall(ISView(perm,PETSC_VIEWER_STDOUT_SELF));
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&luinfo));
52c4762a1bSJed Brown 
531a6d72e3SSatish 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 
589566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&LU));
599566063dSJacob Faibussowitsch   PetscCall(MatLUFactorSymbolic(LU,C,perm,iperm,&luinfo));
609566063dSJacob Faibussowitsch   PetscCall(MatLUFactorNumeric(LU,C,&luinfo));
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,m*n,&u));
639566063dSJacob Faibussowitsch   PetscCall(VecSet(u,one));
649566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&x));
659566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&b));
66c4762a1bSJed Brown 
679566063dSJacob Faibussowitsch   PetscCall(MatMult(C,u,b));
689566063dSJacob Faibussowitsch   PetscCall(MatSolve(LU,b,x));
69c4762a1bSJed Brown 
709566063dSJacob Faibussowitsch   PetscCall(VecView(b,PETSC_VIEWER_STDOUT_SELF));
719566063dSJacob Faibussowitsch   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_SELF));
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x,-1.0,u));
749566063dSJacob Faibussowitsch   PetscCall(VecNorm(x,NORM_2,&norm));
759566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",(double)norm));
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(C,MAT_LOCAL,&info));
789566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"original matrix nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used));
799566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(LU,MAT_LOCAL,&info));
809566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"factored matrix nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used));
81c4762a1bSJed Brown 
829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
859566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
869566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iperm));
879566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
889566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&LU));
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
91b122ec5aSJacob Faibussowitsch   return 0;
92c4762a1bSJed Brown }
93c4762a1bSJed Brown 
94c4762a1bSJed Brown /*TEST
95c4762a1bSJed Brown 
96c4762a1bSJed Brown    test:
971a6d72e3SSatish Balay       suffix: 1
988cc725e6SPierre Jolivet       filter: grep -v " MPI process"
991a6d72e3SSatish Balay 
1001a6d72e3SSatish Balay    test:
1011a6d72e3SSatish Balay       suffix: 2
1021a6d72e3SSatish Balay       args: -m 1 -n 1 -fill 0.49
1038cc725e6SPierre Jolivet       filter: grep -v " MPI process"
104c4762a1bSJed Brown 
105c4762a1bSJed Brown TEST*/
106