Actual source code: ex128.c


  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>

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

 26:   PetscInitialize(&argc,&args,(char*)0,help);
 27:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 29:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 30:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 31:   PetscOptionsGetInt(NULL,NULL,"-lf",&lf,NULL);

 33:   MatCreate(PETSC_COMM_SELF,&C);
 34:   MatSetSizes(C,m*n,m*n,m*n,m*n);
 35:   MatSetFromOptions(C);
 36:   MatSetUp(C);

 38:   /* Create matrix C in seqaij format and sC in seqsbaij. (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;  Ii = j + n*i;
 42:       J = Ii - n; if (J>=0)  MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);
 43:       J = Ii + n; if (J<m*n) MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);
 44:       J = Ii - 1; if (J>=0)  MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);
 45:       J = Ii + 1; if (J<m*n) MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);
 46:       v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 47:     }
 48:   }
 49:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 50:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 52:   MatIsSymmetric(C,0.0,&flg);
 54:   MatConvert(C,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sC);

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

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

 67:   /* Compute CHOLESKY or ICC factor sA */
 68:   MatFactorInfoInitialize(&info);

 70:   info.fill          = 1.0;
 71:   info.diagonal_fill = 0;
 72:   info.zeropivot     = 0.0;

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

 83:     MatGetFactor(sC,MATSOLVERPETSC,MAT_FACTOR_ICC,&sA);
 84:     MatICCFactorSymbolic(sA,sC,row,&info);
 85:   }
 86:   MatCholeskyFactorNumeric(sA,sC,&info);

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

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

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

126: /*TEST

128:    test:
129:       output_file: output/ex128.out

131:    test:
132:       suffix: 2
133:       args: -cholesky -triangular_solve

135: TEST*/