static char help[] = "Tests ILU and ICC factorization with and without matrix ordering on seqsbaij format. Modified from ex30.c\n\ Input parameters are:\n\ -lf : level of fill for ILU (default is 0)\n\ -lu : use full LU or Cholesky factorization\n\ -m ,-n : grid dimensions\n\ Note that most users should employ the KSP interface to the\n\ linear solvers instead of using the factorization routines\n\ directly.\n\n"; #include int main(int argc,char **args) { Mat C,sC,sA; PetscInt i,j,m = 5,n = 5,Ii,J,lf = 0; PetscErrorCode ierr; PetscBool CHOLESKY=PETSC_FALSE,TRIANGULAR=PETSC_FALSE,flg; PetscScalar v; IS row,col; MatFactorInfo info; Vec x,y,b,ytmp; PetscReal norm2,tol = 100*PETSC_MACHINE_EPSILON; PetscRandom rdm; PetscMPIInt size; ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-lf",&lf,NULL);CHKERRQ(ierr); ierr = MatCreate(PETSC_COMM_SELF,&C);CHKERRQ(ierr); ierr = MatSetSizes(C,m*n,m*n,m*n,m*n);CHKERRQ(ierr); ierr = MatSetFromOptions(C);CHKERRQ(ierr); ierr = MatSetUp(C);CHKERRQ(ierr); /* Create matrix C in seqaij format and sC in seqsbaij. (This is five-point stencil with some extra elements) */ for (i=0; i=0) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} J = Ii + n; if (J=0) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} J = Ii + 1; if (J tol) { ierr = PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g\n",(double)norm2);CHKERRQ(ierr); } } } ierr = MatSolve(sA,b,y);CHKERRQ(ierr); ierr = MatDestroy(&sC);CHKERRQ(ierr); ierr = MatDestroy(&sA);CHKERRQ(ierr); ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); if (lf == -1 && norm2 > tol) { ierr = PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ: Cholesky/ICC levels %D, residual %g\n",lf,(double)norm2);CHKERRQ(ierr); } /* Free data structures */ ierr = MatDestroy(&C);CHKERRQ(ierr); ierr = ISDestroy(&row);CHKERRQ(ierr); ierr = ISDestroy(&col);CHKERRQ(ierr); ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&y);CHKERRQ(ierr); ierr = VecDestroy(&ytmp);CHKERRQ(ierr); ierr = VecDestroy(&b);CHKERRQ(ierr); ierr = PetscFinalize(); return ierr; } /*TEST test: output_file: output/ex128.out test: suffix: 2 args: -cholesky -triangular_solve TEST*/