[petsc-users] GMRES for outer solver

Matthew Knepley knepley at gmail.com
Mon May 2 07:12:03 CDT 2022


On Mon, May 2, 2022 at 12:23 AM Ramakrishnan Thirumalaisamy <
rthirumalaisam1857 at sdsu.edu> wrote:

> Thank you. I have a couple of questions. I am solving the low Mach
> Navier-Stokes system using a projection preconditioner (pc_shell type) with
> GMRES being the outer solver and Richardson being the Krylov
> preconditioner. The solver diverges when ksp_pc_type is "right”:
>
> Linear stokes_ solve did not converge due to DIVERGED_NANORINF iterations 0
>

NaN can always be tracked back. I recommend tracing it back to the first
NaN produced. My guess is that your equation of state if producing a NaN.

Also, we have an example of low Mach flow in TS ex76.

  Thanks,

    Matt


> and it converges when ksp_pc_type is "left":
>
> Residual norms for stokes_ solve.
>   0 KSP preconditioned resid norm 8.829128536017e+04 true resid norm
>     -nan ||r(i)||/||b||           -nan
>   1 KSP preconditioned resid norm 1.219313641627e+00 true resid norm
>     -nan ||r(i)||/||b||           -nan
>   2 KSP preconditioned resid norm 8.547033285706e-12 true resid norm
>     -nan ||r(i)||/||b||           -nan
> Linear stokes_ solve converged due to CONVERGED_RTOL iterations 2
>
>  I am curious to know why this is happening. The solver also diverges with
> "FGMRES" as the outer solver (which supports only right preconditioning).
>
> 2. Is it also possible to not get "-nan" when || b || = 0?
>
>
> Regards,
> Rama
>
> On Sun, May 1, 2022 at 12:12 AM Dave May <dave.mayhem23 at gmail.com> wrote:
>
>>
>>
>> On Sun 1. May 2022 at 07:03, Amneet Bhalla <mail2amneet at gmail.com> wrote:
>>
>>> How about using a fixed number of Richardson iterations as a Krylov
>>> preconditioner to a GMRES solver?
>>>
>>
>> That is fine.
>>
>> Would that lead to a linear operation?
>>>
>>
>> Yes.
>>
>>
>>
>>> On Sat, Apr 30, 2022 at 8:21 PM Jed Brown <jed at jedbrown.org> wrote:
>>>
>>>> In general, no. A fixed number of Krylov iterations (CG, GMRES, etc.)
>>>> is a nonlinear operation.
>>>>
>>>> A fixed number of iterations of a method with a fixed polynomial, such
>>>> as Chebyshev, is a linear operation so you don't need a flexible outer
>>>> method.
>>>>
>>>> Ramakrishnan Thirumalaisamy <rthirumalaisam1857 at sdsu.edu> writes:
>>>>
>>>> > Hi,
>>>> >
>>>> > I have a Krylov solver with a preconditioner that is also a Krylov
>>>> solver.
>>>> > I know I can use "fgmres" for the outer solver but can I use gmres
>>>> for the
>>>> > outer solver with a fixed number of iterations in the Krylov
>>>> > preconditioners?
>>>> >
>>>> >
>>>> > Thanks,
>>>> > Rama
>>>>
>>> --
>>> --Amneet
>>>
>>>
>>>
>>>

-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener

https://www.cse.buffalo.edu/~knepley/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20220502/d7e13188/attachment.html>


More information about the petsc-users mailing list