[petsc-users] incomplete LU preconditioner

Barry Smith bsmith at mcs.anl.gov
Tue Jun 7 10:28:54 CDT 2011


  Sorry for the confusion, it is our fault. We don't have a drop tolerance ILU hence PCFactorSetDropTolerance() isn't doing anything.

  You can run with -pc_factor_levels  M   or call PCFactorSetLevels(pc,M) in the code  where M is the size of the matrix.

   Barry

On Jun 7, 2011, at 10:20 AM, Danesh Daroui wrote:

> 
> Hi,
> 
> I want to use incomplete LU preconditioner. My code works OK when I use
> no preconditioner and direct solver for test:
> 
> 
> ierr=KSPCreate(PETSC_COMM_WORLD, &ksp);
> ierr=KSPSetOperators(ksp, Mp, Mp, DIFFERENT_NONZERO_PATTERN);
> ierr=KSPSetTolerances(ksp, 1.e-2/Msize, 1.e-50, PETSC_DEFAULT,
> PETSC_DEFAULT);
> 
> // GMRES iterative solver is used
> ierr=KSPSetType(ksp, KSPGMRES);
> 
> // no PETSc preconditioner is used
> ierr=KSPGetPC(ksp, &prec);
> ierr=PCSetType(prec, PCLU);
> 
> // set up the solver according to the options
> ierr=KSPSetFromOptions(ksp);
> 
> ierr=KSPSetUp(ksp);
> 
> // solve the equation using an iterative solver
> ierr=KSPSolve(ksp, bp, xp);
> 
> 
> Again for some tests, at the first step, I want to use ILU but create
> the exact LU and then change the parameters to tune the preconditioner
> and I almost know what values to set according to my tests in MATLAB. I
> use the code below to create the get exact LU, using ILU preconditioner:
> 
> 
> 
> 
> 
> ierr=KSPCreate(PETSC_COMM_WORLD, &ksp);
> ierr=KSPSetOperators(ksp, Mp, Mp, DIFFERENT_NONZERO_PATTERN);
> ierr=KSPSetTolerances(ksp, 1.e-2/Msize, 1.e-50, PETSC_DEFAULT,
> PETSC_DEFAULT);
> 
> // GMRES iterative solver is used
> ierr=KSPSetType(ksp, KSPGMRES);
> 
> // no PETSc preconditioner is used
> ierr=KSPGetPC(ksp, &prec);
> ierr=PCSetType(prec, PCILU);
> 
> PetscReal dt=1e-10;
> PetscReal dtcol=0.1;
> PetscInt maxrowcount=M;
> PCFactorSetDropTolerance(prec, dt, dtcol, maxrowcount);
> 
> // set up the solver according to the options
> ierr=KSPSetFromOptions(ksp);
> 
> ierr=KSPSetUp(ksp);
> 
> // solve the equation using an iterative solver
> ierr=KSPSolve(ksp, bp, xp);
> 
> 
> and my coefficient matrix is a square "MxM" matrix and I pass "M" as
> "maxrowcount". Using this code, the solver never converges while it
> solved the equation using previous code (direct solver). Can anybody let
> me know how can I get exact LU using ILU and what parameters should I
> change in the above code?
> 
> Regards,
> 
> D.
> 
> 



More information about the petsc-users mailing list