Actual source code: ex30.c

  1: /*$Id: ex30.c,v 1.27 2001/08/07 03:03:07 balay Exp $*/

  3: static char help[] = "Tests ILU factorization and illustrates drawing of matrix sparsity structure with MatView().\n\
  4:   Input parameters are:\n\
  5:   -lf <level> : level of fill for ILU (default is 0)\n\
  6:   -lu : use full LU factorization\n\
  7:   -m <value>,-n <value> : grid dimensions\n\
  8: Note that most users should employ the SLES interface to the\n\
  9: linear solvers instead of using the factorization routines\n\
 10: directly.\n\n";

 12:  #include petscmat.h

 16: int main(int argc,char **args)
 17: {
 18:   Mat         C,A;
 19:   int         i,j,m = 5,n = 5,I,J,ierr,lf = 0;
 20:   PetscTruth  flg1;
 21:   PetscScalar v;
 22:   IS          row,col;
 23:   PetscViewer viewer1,viewer2;

 25:   PetscInitialize(&argc,&args,(char *)0,help);
 26:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 27:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 28:   PetscOptionsGetInt(PETSC_NULL,"-lf",&lf,PETSC_NULL);

 30:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,400,400,&viewer1);
 31:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,400,0,400,400,&viewer2);

 33:   MatCreate(PETSC_COMM_SELF,m*n,m*n,m*n,m*n,&C);
 34:   MatSetFromOptions(C);
 35:   MatSeqBDiagSetPreallocation(C,0,1,PETSC_NULL,PETSC_NULL);
 36:   MatSeqAIJSetPreallocation(C,5,PETSC_NULL);

 38:   /* Create the matrix. (This is five-point stencil with some extra elements) */
 39:   for (i=0; i<m; i++) {
 40:     for (j=0; j<n; j++) {
 41:       v = -1.0;  I = j + n*i;
 42:       J = I - n; if (J>=0)  {MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 43:       J = I + n; if (J<m*n) {MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 44:       J = I - 1; if (J>=0)  {MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 45:       J = I + 1; if (J<m*n) {MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 46:       v = 4.0; MatSetValues(C,1,&I,1,&I,&v,INSERT_VALUES);
 47:     }
 48:   }
 49:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 50:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 52:   MatGetOrdering(C,MATORDERING_RCM,&row,&col);
 53:   printf("original matrix:\n");
 54:   PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_INFO);
 55:   MatView(C,PETSC_VIEWER_STDOUT_SELF);
 56:   PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);
 57:   MatView(C,PETSC_VIEWER_STDOUT_SELF);
 58:   MatView(C,viewer1);

 60:   /* Compute factorization */
 61:   PetscOptionsHasName(PETSC_NULL,"-lu",&flg1);
 62:   if (flg1){
 63:     MatLUFactorSymbolic(C,row,col,PETSC_NULL,&A);
 64:   } else {
 65:     MatFactorInfo info;
 66:     info.levels        = lf;
 67:     info.fill          = 1.0;
 68:     info.diagonal_fill = 0;
 69:     info.damping       = 0;
 70:     info.zeropivot     = 0.0;
 71:     MatILUFactorSymbolic(C,row,col,&info,&A);
 72:   }
 73:   MatLUFactorNumeric(C,&A);

 75:   printf("factored matrix:\n");
 76:   PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_INFO);
 77:   MatView(A,PETSC_VIEWER_STDOUT_SELF);
 78:   PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);
 79:   MatView(A,PETSC_VIEWER_STDOUT_SELF);
 80:   MatView(A,viewer2);

 82:   /* Free data structures */
 83:   MatDestroy(C);
 84:   MatDestroy(A);
 85:   ISDestroy(row);
 86:   ISDestroy(col);
 87:   PetscViewerDestroy(viewer1);
 88:   PetscViewerDestroy(viewer2);
 89:   PetscFinalize();
 90:   return 0;
 91: }