xref: /petsc/src/mat/tests/ex7.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscErrorCode ierr;
14c4762a1bSJed Brown   PetscScalar    v,one = 1.0;
15c4762a1bSJed Brown   IS             perm,iperm;
16c4762a1bSJed Brown   Vec            x,u,b;
171a6d72e3SSatish Balay   PetscReal      norm,fill;
18c4762a1bSJed Brown   MatFactorInfo  luinfo;
19c4762a1bSJed Brown 
20*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
211a6d72e3SSatish Balay 
221a6d72e3SSatish Balay   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Mat test ex7 options","Mat");CHKERRQ(ierr);
231a6d72e3SSatish Balay   m = 3; n = 3; fill = 2.0;
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-m","Number of rows in grid",NULL,m,&m,NULL));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-n","Number of columns in grid",NULL,n,&n,NULL));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-fill","Expected fill ratio for factorization",NULL,fill,&fill,NULL));
271a6d72e3SSatish Balay 
281a6d72e3SSatish Balay   ierr = PetscOptionsEnd();CHKERRQ(ierr);
291a6d72e3SSatish Balay 
30c4762a1bSJed Brown   /* Create the matrix for the five point stencil, YET AGAIN */
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&C));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
385f80ce2aSJacob Faibussowitsch       if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
395f80ce2aSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
405f80ce2aSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
415f80ce2aSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
425f80ce2aSJacob Faibussowitsch       v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
43c4762a1bSJed Brown     }
44c4762a1bSJed Brown   }
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(C,MATORDERINGRCM,&perm,&iperm));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(perm,PETSC_VIEWER_STDOUT_SELF));
50c4762a1bSJed Brown 
515f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&LU));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactorSymbolic(LU,C,perm,iperm,&luinfo));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactorNumeric(LU,C,&luinfo));
61c4762a1bSJed Brown 
625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m*n,&u));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(u,one));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(u,&x));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(u,&b));
66c4762a1bSJed Brown 
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(C,u,b));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSolve(LU,b,x));
69c4762a1bSJed Brown 
705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(b,PETSC_VIEWER_STDOUT_SELF));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_SELF));
72c4762a1bSJed Brown 
735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(x,-1.0,u));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(x,NORM_2,&norm));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",(double)norm));
76c4762a1bSJed Brown 
775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetInfo(C,MAT_LOCAL,&info));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"original matrix nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetInfo(LU,MAT_LOCAL,&info));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"factored matrix nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used));
81c4762a1bSJed Brown 
825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&perm));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iperm));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&LU));
89c4762a1bSJed Brown 
90*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
91*b122ec5aSJacob Faibussowitsch   return 0;
92c4762a1bSJed Brown }
93c4762a1bSJed Brown 
94c4762a1bSJed Brown /*TEST
95c4762a1bSJed Brown 
96c4762a1bSJed Brown    test:
971a6d72e3SSatish Balay       suffix: 1
981a6d72e3SSatish Balay       filter: grep -v "MPI processes"
991a6d72e3SSatish Balay 
1001a6d72e3SSatish Balay    test:
1011a6d72e3SSatish Balay       suffix: 2
1021a6d72e3SSatish Balay       args: -m 1 -n 1 -fill 0.49
103c4762a1bSJed Brown       filter: grep -v "MPI processes"
104c4762a1bSJed Brown 
105c4762a1bSJed Brown TEST*/
106