GMRES left-preconditioned with ILU1 versus ILU0

Barry Smith bsmith at mcs.anl.gov
Fri May 16 11:42:47 CDT 2008


    When monkeying with ILU it is always useful to run tests with - 
ksp_monitor_true_residual, because bad pivots in ILU
can produce wildly scaled preconditioners so the preconditioned  
residual can be at a different scale then the actual residual.

   The block versions of the factorizations don't support the  
pc_factor_shiftXXXX options. You could try to add it.

   That is one nasty matrix, I don't know if you'll ever get  
reasonable performance with ILU. You may want to stick to
direct solvers on each process (it may be slow and use lots of memory  
but at least it might work).

   Barry


On May 16, 2008, at 11:22 AM, Stephane Aubert wrote:

> Hi,
> I tried to improve the convergence of a rather badly conditioned  
> linear system by increasing the ILU(k) level from k=0 to k=1.
> And, whereas I got a convergence for k=0, I ended up with an  
> explosive divergence for k=1!!!
> The PETCS version is 2.3.2-p8.
> The matrix type is MPIBAIJ (block size=7,mesh points=23010,non-empty  
> blocks=297454)
> The common command line options are:
>
>   * KSP="-ksp_type gmres -ksp_max_it 800 -ksp_gmres_restart 800
>     -ksp_rtol 1.0e-12 -ksp_left_pc -ksp_gmres_modifiedgramschmidt
>     -ksp_gmres_cgs_refinement_type REFINE_NEVER -ksp_singmonitor
>     -ksp_compute_eigenvalues": I'm forcing 800 krylovs without restart
>     to get the condition number of the pre-conditioned system (if I
>     understand correctly the man page of -ksp_singmonitor)
>   * PC="-pc_type asm -pc_asm_overlap 2": I'm planning to run with more
>     than 1 partition, but for the time being, only one partition is  
> used.
>   * BLK_KSP="-sub_ksp_type preonly": Because of GMRES+ILU
>
> For ILU(0), I'm using:
>
>   * BLK_PC="-sub_pc_type ilu -sub_pc_factor_levels 0
>     -sub_pc_factor_fill 1.00 -sub_pc_factor_shift_nonzero
>     -sub_pc_factor_mat_ordering_type rcm - 
> sub_pc_factor_pivot_in_blocks"
>
> and I got for convergence:
> 0 KSP Residual norm 1.687258996558e+00 % max 1 min 1 max/min 1
> 1 KSP Residual norm 1.687132576728e+00 % max 67.8829 min 67.8829 max/ 
> min 1
> 2 KSP Residual norm 1.685760733293e+00 % max 3496.78 min 19.5582 max/ 
> min 178.789
> 3 KSP Residual norm 1.668552043995e+00 % max 3604.26 min 12.0073 max/ 
> min 300.174
> 4 KSP Residual norm 1.578511835381e+00 % max 3639.92 min 7.35118 max/ 
> min 495.148
> ....
> 795 KSP Residual norm 1.465607932165e-09 % max 18209.9 min  
> 0.00612973 max/min 2.97075e+06
> 796 KSP Residual norm 1.390602265424e-09 % max 18227.3 min  
> 0.00612913 max/min 2.97388e+06
> 797 KSP Residual norm 1.320529491862e-09 % max 18231.9 min 0.0061286  
> max/min 2.97489e+06
> 798 KSP Residual norm 1.253371917713e-09 % max 18234.8 min  
> 0.00612856 max/min 2.97538e+06
> 799 KSP Residual norm 1.188955299647e-09 % max 18278.5 min  
> 0.00612594 max/min 2.98378e+06
> 800 KSP Residual norm 1.118294486519e-09 % max 18278.5 min  
> 0.00612475 max/min 2.98437e+06
>
> and the iterative solution compares very well with the one computed  
> using complete LU factorization.
>
> For ILU(1), I'm using:
>
>   * BLK_PC="-sub_pc_type ilu -sub_pc_factor_levels 1
>     -sub_pc_factor_fill 3.81 -sub_pc_factor_shift_nonzero
>     -sub_pc_factor_mat_ordering_type rcm
>     -sub_pc_factor_pivot_in_blocks": RCM gives the smallest fill  
> value.
>
> and I got for "convergence":
> 0 KSP Residual norm 7.095990612421e+126 % max 1 min 1 max/min 1
> 1 KSP Residual norm 3.313547979190e+123 % max 1.68012e+135 min  
> 1.68012e+135 max/min 1
> 2 KSP Residual norm 1.257750994639e+119 % max 6.34518e+135 min  
> 3.55953e+131 max/min 17825.9
> 3 KSP Residual norm 5.233083258710e+118 % max 1.25538e+136 min  
> 1.42732e+127 max/min 8.79538e+08
> 4 KSP Residual norm 1.938981257595e+118 % max 1.44472e+136 min  
> 4.82369e+125 max/min 2.99506e+10
> 5 KSP Residual norm 3.371270926617e+116 % max 1.45841e+136 min  
> 1.79839e+125 max/min 8.1095e+10
> 6 KSP Residual norm 2.179293254483e+115 % max 1.45842e+136 min  
> 5.16422e+122 max/min 2.82408e+13
> 7 KSP Residual norm 2.120598006522e+115 % max 1.46024e+136 min  
> 3.67337e+122 max/min 3.97521e+13
> 8 KSP Residual norm 1.486820601733e+115 % max 1.461e+136 min 1.02249e 
> +122 max/min 1.42886e+14
> 9 KSP Residual norm 7.653834441859e+114 % max 1.46138e+136 min  
> 4.93314e+121 max/min 2.96237e+14
> 10 KSP Residual norm 5.001204920218e+114 % max 1.47243e+136 min  
> 4.89675e+121 max/min 3.00694e+14
>
> My guess is that the ILU(1) is singular (zero as diagonal  
> elements?), but I thought that the options "- 
> sub_pc_factor_shift_nonzero -sub_pc_factor_pivot_in_blocks" were  
> taking care of that... I got lost in the source files trying to find  
> out who at the end is computing and applying the ILU for MPIBAIJ  
> format (replaced by SEQBAIJ with only one partition, I'm guessing).
>
> The question is: What am I doing wrong? I never heard that ILU(1)  
> was worst than ILU(0)!
> Stef.
>
> -- 
> ___________________________________________________________
> Dr. Stephane AUBERT, CEO & CTO
> FLUOREM s.a.s
> Centre Scientifique Auguste MOIROUX
> 64 chemin des MOUILLES
> F-69130 ECULLY, FRANCE
> International:      fax: +33 4.78.33.99.39      tel: +33 4.78.33.99.35
> France:             fax: 04.78.33.99.39         tel: 04.78.33.99.35
> email: stephane.aubert at fluorem.com
> web: www.fluorem.com
>
>
>




More information about the petsc-users mailing list