[petsc-users] incomplete LU preconditioner

Hong Zhang hzhang at mcs.anl.gov
Tue Jun 7 12:30:39 CDT 2011


Danesh:

>> Thanks, I changed my code to use same preconditioner. Some more
>> questions:
>>
>> 1. I changed the level from 1730 to 28 and I expected "no convergence"
>> or at convergence at more iterations or higher residual, but nothing
>> changed. What can be the reason?
>
>   28 levels of fill is equivalent/results in a comple LU factorization

Note, 28 levels of fill in petsc ilu() is different from what your
mentioned "M" as "maxrowcount"
in iludt().
What Barry means here is that ilu(28) is mathematically equivalent to
a complete LU,
but not as efficient as complete LU. Try to use a direct solver, e.g., mumps,
(use runtime option '-pc_type lu -pc_factor_mat_solver_package mumps),
likely you'll
get much faster time than ilu(28).

Few years ago, I wrote iludt for petsc and tested it on a collection
of matrices.
Comparing it with ilu(p) in petsc, I did not see  any advantage of
iludt over ilu.
Do any of you have reference that demonstrates the advantage of iludt over ilu?

Computationally, iludt requires full LU computation for each row,
then discard entries that either smaller than droptol or take more
spaces than "maxrowcount"
with row sort.
Analytically, dropping small entries not necessarily maintains the
quality of the LU
factor. Without row/column permutation, L factor gives larger
multipliers resulting
bad convergences.
We realized that much more are needed to make our iludt useful.
Due to time constraint and lack of good idea, I stopped this work :-(

Hong

>> 2. Using same preconditioner I see no change in the solution time.
>
>   Run with -log_summary to see where the code is spending the time. If the ILU numeric factorization is not taking much of the run time then doing it less often will not be noticable.
>
>   Barry
>
>>
>> In all cases I get the convergence at max iteration 2 with very low
>> residual, but it is a bit slow.
>>
>> Thanks,
>>
>> D.
>>
>>
>> On Tue, 2011-06-07 at 10:49 -0500, Barry Smith wrote:
>>> On Jun 7, 2011, at 10:44 AM, Danesh Daroui wrote:
>>>
>>>>
>>>> So, it is not implemented? That's too bad! ILU with drop threshold is
>>>> very useful so I am wondering why it is not implemented yet? Do you have
>>>> any plan to implement it?
>>>
>>>   We've found that providing a high-quality drop tolerance ILU is far from trivial and haven't had the time/person with expertise to devote the time needed to implement it.
>>>
>>>>
>>>> Setting levels did work OK for me. I have another question. If I want to
>>>> calculate the preconditioner once, and use it later to solve other
>>>> equations, since I know that the preconditioner would be OK because I
>>>> solve several equations in several iterations and I know that my
>>>> coefficient matrix will not change that much so calculating the
>>>> preconditioner at once would be OK. What routines should I call once and
>>>> how can I apply the pre-calculated preconditioner each time I want to
>>>> solve the equation?
>>>
>>>   Since you are using the linear solvers directly. If you change the right hand side but NOT the matrix you simply call KSPSolve() again with the new right hand side.
>>>
>>>   If you are changing the matrix but do not want to update the preconditioner call KSPSetOperators() with the new matrix but use a final argument of SAME_PRECONDITIONER
>>>
>>>   Barry
>>>
>>>>
>>>> Thanks,
>>>>
>>>> D.
>>>>
>>>>
>>>>
>>>> On Tue, 2011-06-07 at 10:28 -0500, Barry Smith wrote:
>>>>> 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