There is one particular integration method that may make the life easier — 
backward Euler. It has only one stage at each time step.

To Francesco: which TS method are you using? Is it backward Euler or the 
default theta method? The default one is not stiffly accurate, thus not good 
for stiff problems.

Hong 
 
> On Feb 10, 2016, at 2:06 PM, Barry Smith <[email protected]> wrote:
> 
> 
>   Francesco,
> 
>     Unfortunately the -snes_grid_sequence option only works with a single 
> SNES nonlinear solve and not directly within an outer time-integration scheme.
> 
>     We do not currently have library code to do grid sequence within a 
> time-step. I believe that providing this requires a great deal of 
> "replumbing" of the TS solvers to make them "grid sequence aware". This is 
> true I believe for both traditional grid sequencing and FAS (I would be happy 
> if Emil or Jed could prove me wrong).  There are seemingly multiple ways to 
> do the grid sequencing: 
> 
> 1) Restrict the current solution to a coarse mesh and do a single stage of 
> the ODE integrator (which means solve a single nonlinear solve) on that 
> coarse mesh, interpolate the solution to the next mesh and use it as an 
> initial guess for the single stage on that mesh and then continue to up the 
> finest mesh. Note that since the specific ODE integrator determines the 
> "recipe" for the SNES function evaluation this is not something that can be 
> done in pure user code, it must use the "plumbing" of the specific ODE 
> integrator.
> 
> an alternative is
> 2) Restrict the current solution to a coarse mesh and do the entire timestep 
> of the ODE integrator (which means multiple stage and multiple nonlinear 
> solvers) saving all the intermediate solutions (for each stage). Interpolate 
> all the intermediate solutions (as well as the final solution) to the finer 
> mesh and repeat the entire timestep computation but using the interpolated 
> intermediate solutions as the initial guesses for the nonlinear solves). 
> Continue to the finest mesh. Again this cannot be done in user code and 
> requires replumbing the ODE integrator.
> 
>   I think 2 will produce different "results" from 1 since for 1 the coarse 
> grid initial values on (for example) the second stage will be obtained from 
> restricting the finest grid solution from the first stage while with 2 the 
> initial values on the second stage will be directly the coarse solution from 
> the first stage.
> 
> To do the FAS is even more involved, I can only conceive of doing it for a 
> single stage. Again the FAS has to be plumbed directly into the specific ODE 
> solver to define the correct functions for the various levels.
> 
> With the above approaches the coarser nonlinear problems should be "easier" 
> than the finest because the time-step remains small but the space 
> discretization is larger. But the only way the approaches can work is if the 
> interpolated coarse grid solution used on the finer mesh is within the basin 
> of attraction on that mesh, otherwise you will still have the same problem 
> you have now of lack of convergence.
> 
> I don't know if anyone actually does this and how much it helps or could 
> help. I welcome input from others.
> 
> Based on your convergence experience it sounds like your solution is changing 
> rather rapidly in time? Is the solution somewhat smooth in space? If not then 
> I don't think grid sequencing will help much. 
> 
> 
>  Barry
> 
> 
> 
> 
> 
>> On Feb 10, 2016, at 9:38 AM, Francesco Magaletti 
>> <[email protected]> wrote:
>> 
>> Hi everyone,
>> 
>> I’m trying to solve a system of time dependent PDE’s, strongly non linear 
>> and quite stiff. The system is 1D, but I need a large number of grid points 
>> (1-10 million points). I’m facing some convergence problem of the Newton 
>> solver: if I use a small timestep (10^-5) the SNES converges correctly, but 
>> if I increase it, the Newton method diverges. Moreover, when I increase the 
>> number of grid points, the timestep for Newton convergence decreases. Such a 
>> restriction on the timestep is a big problem for my simulations.
>> 
>> I thought that the option -snes_grid_sequence could help the SNES to 
>> converge, but when I run the code the first grid converges and then the 
>> simulations abort with errors. Here is an example with 60000grid points:
>> 
>> mpiexec -n 4 ./code_petsc -ts_monitor -snes_monitor -snes_converged_reason 
>> -ksp_converged_reason -ftime 8.0 -ksp_rtol 1.0e-8 -snes_atol 1.0e-12 
>> -snes_rtol 1.0e-14 -snes_stol 1.0e-12 -snes_max_it 2000 -dt 2.0e-4 -nx 60000 
>> -snes_grid_sequence 1
>> 
>> 0 TS dt 0.0002 time 0
>>     0 SNES Function norm 4.597340391780e+04 
>>     Linear solve converged due to CONVERGED_RTOL iterations 9
>>     1 SNES Function norm 3.986169680241e+03 
>>     Linear solve converged due to CONVERGED_RTOL iterations 14
>>     2 SNES Function norm 4.292840674611e+01 
>>     Linear solve converged due to CONVERGED_RTOL iterations 15
>>     3 SNES Function norm 1.044935520416e-02 
>>     Linear solve converged due to CONVERGED_RTOL iterations 14
>>     4 SNES Function norm 2.936744161136e-06 
>>     Linear solve converged due to CONVERGED_RTOL iterations 10
>>     5 SNES Function norm 4.471992600451e-08 
>>   Nonlinear solve converged due to CONVERGED_SNORM_RELATIVE iterations 5
>> [0]PETSC ERROR: --------------------- Error Message 
>> --------------------------------------------------------------
>> [0]PETSC ERROR:   
>> [0]PETSC ERROR: Must call TSSetRHSFunction() and / or TSSetIFunction()
>> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for 
>> trouble shooting.
>> [0]PETSC ERROR: Petsc Release Version 3.6.3, Dec, 03, 2015 
>> [0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++ 
>> --with-fc=gfortran --download-fblaslapack --download-mpich 
>> --download-mpich-device=ch3:nemesis
>> [0]PETSC ERROR: #1 TSComputeIFunction() line 723 in 
>> /home/magaletto/Scaricati/petsc-3.6.3/src/ts/interface/ts.c
>> [0]PETSC ERROR: #2 SNESTSFormFunction_Theta() line 483 in 
>> /home/magaletto/Scaricati/petsc-3.6.3/src/ts/impls/implicit/theta/theta.c
>> [0]PETSC ERROR: #3 SNESTSFormFunction() line 4138 in 
>> /home/magaletto/Scaricati/petsc-3.6.3/src/ts/interface/ts.c
>> [0]PETSC ERROR: #4 SNESComputeFunction() line 2067 in 
>> /home/magaletto/Scaricati/petsc-3.6.3/src/snes/interface/snes.c
>> [0]PETSC ERROR: #5 SNESSolve_NEWTONLS() line 184 in 
>> /home/magaletto/Scaricati/petsc-3.6.3/src/snes/impls/ls/ls.c
>> [0]PETSC ERROR: #6 SNESSolve() line 3906 in 
>> /home/magaletto/Scaricati/petsc-3.6.3/src/snes/interface/snes.c
>> [0]PETSC ERROR: #7 TSStep_Theta() line 198 in 
>> /home/magaletto/Scaricati/petsc-3.6.3/src/ts/impls/implicit/theta/theta.c
>> [0]PETSC ERROR: #8 TSStep() line 3101 in 
>> /home/magaletto/Scaricati/petsc-3.6.3/src/ts/interface/ts.c
>> [0]PETSC ERROR: #9 TSSolve() line 3285 in 
>> /home/magaletto/Scaricati/petsc-3.6.3/src/ts/interface/ts.c
>> 
>> ===================================================================================
>> =   BAD TERMINATION OF ONE OF YOUR APPLICATION PROCESSES
>> =   EXIT CODE: 6
>> 
>> 
>> Where is the mistake? 
>> Could you suggest a "more convergent” SNES? Maybe a FAS as a nonlinear 
>> preconditioner? 
>> 
>> Hope someone can help here.
>> best regards
>> 
>> Francesco
> 

Reply via email to