I wasn't expecting this. Try

-pc_mg_type full -ksp_type richardson -mg_levels_pc_type bjacobi 
-mg_levels_ksp_type gmres -mg_levels_ksp_max_it 3 
-mg_coarse_pc_factor_mat_solver_package superlu_dist -mg_coarse_pc_type lu 
-pc_mg_galerkin -pc_mg_levels 5 -pc_type mg 
-log_summary  -pc_mg_log  -ksp_monitor_true_residual  -options_left   -ksp_view 
-ksp_max_it 30




On Sep 30, 2013, at 8:00 PM, Michele Rosso <[email protected]> wrote:

> Barry,
> 
> sorry again for the very late answer. I tried all the variations you 
> proposed: all of them converge very slow , except the last one (CG instead of 
> fgmres) diverges.
> I attached the diagnostics for the two options that convergence: each of the 
> attached files starts with a list of the options I used for the run.
> As you pointed out earlier, the residual I used to require was too small, 
> therefore I increased atol to e-9.
> After some tests, I noticed that any further increase of the absolute 
> tolerance changes significantly the solution.
> What would you suggest to try next?
> Thank you very much,
> 
> Michele
> 
> 
> 
> 
> 
> 
> On 09/24/2013 05:08 PM, Barry Smith wrote:
>>  Thanks. The balance of work on the different levels and across processes 
>> looks ok. So it is is a matter of improving the convergence rate.
>> 
>>  The initial residual norm is very small. Are you sure you need to decrease 
>> it to 10^-12 ????
>> 
>>  Start with a really robust multigrid smoother use
>> 
>> -pc_mg_type full -ksp_type fgmres -mg_levels_pc_type bjacobi 
>> -mg_levels_ksp_type gmres -mg_levels_ksp_max_it 3   PLUS 
>> -mg_coarse_pc_factor_mat_solver_package superlu_dist
>> -mg_coarse_pc_type lu -pc_mg_galerkin -pc_mg_levels 5 -pc_mg_log -pc_type mg
>> 
>>  run with the -log_summary and -pc_mg_log
>> 
>> Now back off a little on the smoother and use -mg_levels_pc_type sor instead 
>>  how does that change the convergence and time.
>> 
>> Back off even more an replace the -ksp_type fgmres with -ksp_type cg and the 
>> -mg_levels_ksp_type gmres with -mg_levels_ksp_type richardson   how does 
>> that change the convergence and the time?
>> 
>>  There are some additional variants we can try based on the results from 
>> above.
>> 
>>  Barry
>> 
>> 
>> 
>> On Sep 24, 2013, at 4:29 PM, Michele Rosso <[email protected]> wrote:
>> 
>>> Barry,
>>> 
>>> I re-rerun the test case with the option -pc_mg_log as you suggested.
>>> I attached the new output ("final_new.txt').
>>> Thanks for your help.
>>> 
>>> Michele
>>> 
>>> On 09/23/2013 09:35 AM, Barry Smith wrote:
>>>>   Run with the additional option -pc_mg_log and send us the log file.
>>>> 
>>>>   Barry
>>>> 
>>>> Maybe we should make this the default somehow.
>>>> 
>>>> 
>>>> On Sep 23, 2013, at 10:55 AM, Michele Rosso <[email protected]> wrote:
>>>> 
>>>>> Hi,
>>>>> 
>>>>> I am successfully using PETSc to solve a 3D Poisson's equation with CG + 
>>>>> MG .  Such equation arises from a projection algorithm for a multiphase 
>>>>> incompressible flow simulation.
>>>>> I set up the solver  as I was suggested to do in a previous thread 
>>>>> (title: "GAMG speed") and run a test case (liquid droplet with surface 
>>>>> tension falling under the effect of gravity in a quiescent fluid).
>>>>> The solution of the Poisson Equation via multigrid is correct but it 
>>>>> becomes progressively slower and slower as the simulation progresses (I 
>>>>> am performing successive solves) due to an increase in the number of 
>>>>> iterations.
>>>>> Since the solution of the Poisson equation is mission-critical, I need to 
>>>>> speed it up as much as I can.
>>>>> Could you please help me out with this?
>>>>> 
>>>>> I run the test case with the following options:
>>>>> 
>>>>> -pc_type mg  -pc_mg_galerkin  -pc_mg_levels 5   -mg_levels_ksp_type 
>>>>> richardson -mg_levels_ksp_max_it 1
>>>>> -mg_coarse_pc_type lu   -mg_coarse_pc_factor_mat_solver_package 
>>>>> superlu_dist
>>>>> -log_summary -ksp_view  -ksp_monitor_true_residual  -options_left
>>>>> 
>>>>> Please find the diagnostic for the final solve in the attached file 
>>>>> "final.txt'.
>>>>> Thank you,
>>>>> 
>>>>> Michele
>>>>> <final.txt>
>>> <final_new.txt>
>> 
> 
> <final1.txt><final2.txt>

Reply via email to