Actual source code: ex7.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2: static char help[] = "Tests matrix factorization.  Note that most users should\n\
  3: employ the KSP  interface to the linear solvers instead of using the factorization\n\
  4: routines directly.\n\n";

  6:  #include <petscmat.h>

  8: int main(int argc,char **args)
  9: {
 10:   Mat            C,LU;
 11:   MatInfo        info;
 12:   PetscInt       i,j,m = 3,n = 3,Ii,J;
 14:   PetscScalar    v,one = 1.0;
 15:   IS             perm,iperm;
 16:   Vec            x,u,b;
 17:   PetscReal      norm;
 18:   MatFactorInfo  luinfo;

 20:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 21:   /* Create the matrix for the five point stencil, YET AGAIN */
 22:   MatCreate(PETSC_COMM_SELF,&C);
 23:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 24:   MatSetFromOptions(C);
 25:   MatSetUp(C);
 26:   for (i=0; i<m; i++) {
 27:     for (j=0; j<n; j++) {
 28:       v = -1.0;  Ii = j + n*i;
 29:       if (i>0)   {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 30:       if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 31:       if (j>0)   {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 32:       if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 33:       v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 34:     }
 35:   }
 36:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 37:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 38:   MatGetOrdering(C,MATORDERINGRCM,&perm,&iperm);
 39:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);
 40:   ISView(perm,PETSC_VIEWER_STDOUT_SELF);

 42:   MatFactorInfoInitialize(&luinfo);

 44:   luinfo.fill          = 2.0;
 45:   luinfo.dtcol         = 0.0;
 46:   luinfo.zeropivot     = 1.e-14;
 47:   luinfo.pivotinblocks = 1.0;

 49:   MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&LU);
 50:   MatLUFactorSymbolic(LU,C,perm,iperm,&luinfo);
 51:   MatLUFactorNumeric(LU,C,&luinfo);

 53:   VecCreateSeq(PETSC_COMM_SELF,m*n,&u);
 54:   VecSet(u,one);
 55:   VecDuplicate(u,&x);
 56:   VecDuplicate(u,&b);

 58:   MatMult(C,u,b);
 59:   MatSolve(LU,b,x);

 61:   VecView(b,PETSC_VIEWER_STDOUT_SELF);
 62:   VecView(x,PETSC_VIEWER_STDOUT_SELF);

 64:   VecAXPY(x,-1.0,u);
 65:   VecNorm(x,NORM_2,&norm);
 66:   PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",(double)norm);

 68:   MatGetInfo(C,MAT_LOCAL,&info);
 69:   PetscPrintf(PETSC_COMM_SELF,"original matrix nonzeros = %D\n",(PetscInt)info.nz_used);
 70:   MatGetInfo(LU,MAT_LOCAL,&info);
 71:   PetscPrintf(PETSC_COMM_SELF,"factored matrix nonzeros = %D\n",(PetscInt)info.nz_used);

 73:   VecDestroy(&u);
 74:   VecDestroy(&b);
 75:   VecDestroy(&x);
 76:   ISDestroy(&perm);
 77:   ISDestroy(&iperm);
 78:   MatDestroy(&C);
 79:   MatDestroy(&LU);

 81:   PetscFinalize();
 82:   return ierr;
 83: }