Actual source code: ex7.c

petsc-3.13.6 2020-09-29
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,n,Ii,J;
 14:   PetscScalar    v,one = 1.0;
 15:   IS             perm,iperm;
 16:   Vec            x,u,b;
 17:   PetscReal      norm,fill;
 18:   MatFactorInfo  luinfo;

 20:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;

 22:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Mat test ex7 options","Mat");
 23:   m = 3; n = 3; fill = 2.0;
 24:   PetscOptionsInt("-m","Number of rows in grid",NULL,m,&m,NULL);
 25:   PetscOptionsInt("-n","Number of columns in grid",NULL,n,&n,NULL);
 26:   PetscOptionsReal("-fill","Expected fill ratio for factorization",NULL,fill,&fill,NULL);

 28:   PetscOptionsEnd();

 30:   /* Create the matrix for the five point stencil, YET AGAIN */
 31:   MatCreate(PETSC_COMM_SELF,&C);
 32:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 33:   MatSetFromOptions(C);
 34:   MatSetUp(C);
 35:   for (i=0; i<m; i++) {
 36:     for (j=0; j<n; j++) {
 37:       v = -1.0;  Ii = j + n*i;
 38:       if (i>0)   {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 39:       if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 40:       if (j>0)   {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 41:       if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 42:       v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 43:     }
 44:   }
 45:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 46:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 47:   MatGetOrdering(C,MATORDERINGRCM,&perm,&iperm);
 48:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);
 49:   ISView(perm,PETSC_VIEWER_STDOUT_SELF);

 51:   MatFactorInfoInitialize(&luinfo);

 53:   luinfo.fill          = fill;
 54:   luinfo.dtcol         = 0.0;
 55:   luinfo.zeropivot     = 1.e-14;
 56:   luinfo.pivotinblocks = 1.0;

 58:   MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&LU);
 59:   MatLUFactorSymbolic(LU,C,perm,iperm,&luinfo);
 60:   MatLUFactorNumeric(LU,C,&luinfo);

 62:   VecCreateSeq(PETSC_COMM_SELF,m*n,&u);
 63:   VecSet(u,one);
 64:   VecDuplicate(u,&x);
 65:   VecDuplicate(u,&b);

 67:   MatMult(C,u,b);
 68:   MatSolve(LU,b,x);

 70:   VecView(b,PETSC_VIEWER_STDOUT_SELF);
 71:   VecView(x,PETSC_VIEWER_STDOUT_SELF);

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

 77:   MatGetInfo(C,MAT_LOCAL,&info);
 78:   PetscPrintf(PETSC_COMM_SELF,"original matrix nonzeros = %D\n",(PetscInt)info.nz_used);
 79:   MatGetInfo(LU,MAT_LOCAL,&info);
 80:   PetscPrintf(PETSC_COMM_SELF,"factored matrix nonzeros = %D\n",(PetscInt)info.nz_used);

 82:   VecDestroy(&u);
 83:   VecDestroy(&b);
 84:   VecDestroy(&x);
 85:   ISDestroy(&perm);
 86:   ISDestroy(&iperm);
 87:   MatDestroy(&C);
 88:   MatDestroy(&LU);

 90:   PetscFinalize();
 91:   return ierr;
 92: }


 95: /*TEST

 97:    test:
 98:       suffix: 1
 99:       filter: grep -v "MPI processes"

101:    test:
102:       suffix: 2
103:       args: -m 1 -n 1 -fill 0.49
104:       filter: grep -v "MPI processes"

106: TEST*/