Actual source code: ex128.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: static char help[] = "Tests ILU and ICC factorization with and without matrix ordering on seqsbaij format. Modified from ex30.c\n\
  3:   Input parameters are:\n\
  4:   -lf <level> : level of fill for ILU (default is 0)\n\
  5:   -lu : use full LU or Cholesky factorization\n\
  6:   -m <value>,-n <value> : grid dimensions\n\
  7: Note that most users should employ the KSP interface to the\n\
  8: linear solvers instead of using the factorization routines\n\
  9: directly.\n\n";

 11: #include <petscmat.h>

 15: int main(int argc,char **args)
 16: {
 17:   Mat            C,sC,sA;
 18:   PetscInt       i,j,m = 5,n = 5,Ii,J,lf = 0;
 20:   PetscBool      CHOLESKY=PETSC_FALSE,TRIANGULAR=PETSC_FALSE,flg;
 21:   PetscScalar    v;
 22:   IS             row,col;
 23:   MatFactorInfo  info;
 24:   Vec            x,y,b,ytmp;
 25:   PetscReal      norm2,tol = 100*PETSC_MACHINE_EPSILON;
 26:   PetscRandom    rdm;
 27:   PetscMPIInt    size;

 29:   PetscInitialize(&argc,&args,(char*)0,help);
 30:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 31:   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
 32:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 33:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 34:   PetscOptionsGetInt(NULL,NULL,"-lf",&lf,NULL);

 36:   MatCreate(PETSC_COMM_SELF,&C);
 37:   MatSetSizes(C,m*n,m*n,m*n,m*n);
 38:   MatSetFromOptions(C);
 39:   MatSetUp(C);

 41:   /* Create matrix C in seqaij format and sC in seqsbaij. (This is five-point stencil with some extra elements) */
 42:   for (i=0; i<m; i++) {
 43:     for (j=0; j<n; j++) {
 44:       v = -1.0;  Ii = j + n*i;
 45:       J = Ii - n; if (J>=0)  {MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 46:       J = Ii + n; if (J<m*n) {MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 47:       J = Ii - 1; if (J>=0)  {MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 48:       J = Ii + 1; if (J<m*n) {MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 49:       v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 50:     }
 51:   }
 52:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 53:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 55:   MatIsSymmetric(C,0.0,&flg);
 56:   if (!flg) SETERRQ(PETSC_COMM_SELF,1,"C is non-symmetric");
 57:   MatConvert(C,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sC);

 59:   /* Create vectors for error checking */
 60:   MatCreateVecs(C,&x,&b);
 61:   VecDuplicate(x,&y);
 62:   VecDuplicate(x,&ytmp);
 63:   PetscRandomCreate(PETSC_COMM_SELF,&rdm);
 64:   PetscRandomSetFromOptions(rdm);
 65:   VecSetRandom(x,rdm);
 66:   MatMult(C,x,b);

 68:   MatGetOrdering(C,MATORDERINGNATURAL,&row,&col);

 70:   /* Compute CHOLESKY or ICC factor sA */
 71:   MatFactorInfoInitialize(&info);

 73:   info.fill          = 1.0;
 74:   info.diagonal_fill = 0;
 75:   info.zeropivot     = 0.0;

 77:   PetscOptionsHasName(NULL,NULL,"-cholesky",&CHOLESKY);
 78:   if (CHOLESKY) {
 79:     printf("Test CHOLESKY...\n");
 80:     MatGetFactor(sC,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sA);
 81:     MatCholeskyFactorSymbolic(sA,sC,row,&info);
 82:   } else {
 83:     printf("Test ICC...\n");
 84:     info.levels = lf;

 86:     MatGetFactor(sC,MATSOLVERPETSC,MAT_FACTOR_ICC,&sA);
 87:     MatICCFactorSymbolic(sA,sC,row,&info);
 88:   }
 89:   MatCholeskyFactorNumeric(sA,sC,&info);

 91:   /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering on aij matrix C */
 92:   if (CHOLESKY) {
 93:     PetscOptionsHasName(NULL,NULL,"-triangular_solve",&TRIANGULAR);
 94:     if (TRIANGULAR) {
 95:       printf("Test MatForwardSolve...\n");
 96:       MatForwardSolve(sA,b,ytmp);
 97:       printf("Test MatBackwardSolve...\n");
 98:       MatBackwardSolve(sA,ytmp,y);
 99:       VecAXPY(y,-1.0,x);
100:       VecNorm(y,NORM_2,&norm2);
101:       if (norm2 > tol) {
102:         PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g\n",(double)norm2);
103:       }
104:     }
105:   }

107:   MatSolve(sA,b,y);
108:   MatDestroy(&sC);
109:   MatDestroy(&sA);
110:   VecAXPY(y,-1.0,x);
111:   VecNorm(y,NORM_2,&norm2);
112:   if (lf == -1 && norm2 > tol) {
113:     PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ:   Cholesky/ICC levels %D, residual %g\n",lf,(double)norm2);
114:   }

116:   /* Free data structures */
117:   MatDestroy(&C);
118:   ISDestroy(&row);
119:   ISDestroy(&col);
120:   PetscRandomDestroy(&rdm);
121:   VecDestroy(&x);
122:   VecDestroy(&y);
123:   VecDestroy(&ytmp);
124:   VecDestroy(&b);
125:   PetscFinalize();
126:   return 0;
127: }