<div dir="ltr"><div dir="ltr">On Wed, Nov 6, 2019 at 1:32 PM Pierre Jolivet via petsc-dev <<a href="mailto:petsc-dev@mcs.anl.gov">petsc-dev@mcs.anl.gov</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div style="word-wrap:break-word"><div><blockquote type="cite"><div>On 6 Nov 2019, at 7:30 PM, Fande Kong <<a href="mailto:fdkong.jd@gmail.com" target="_blank">fdkong.jd@gmail.com</a>> wrote:</div><br><div><div dir="ltr"><div dir="ltr">Even if both the linear operator and the preconditioner  operator are symmetric, the preconditioned operator by RAS will  NOT be symmetric so that CG is not applicable. <div><br></div><div>RAS keeps overlapping elements during restriction but discards them when doing interpolation, which leads to an unsymmetric preconditioned operator.</div></div></div></div></blockquote><div><br></div><div>Of course.</div><div>What’s your point?</div></div></div></blockquote><div><br></div><div>I guess Fande's point is that we should consider switching to CG or MINRES, since we are being careful about the PC.</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div style="word-wrap:break-word"><div><div>Thanks,</div><div>Pierre</div><br><blockquote type="cite"><div><div dir="ltr"><div dir="ltr"><div>Fande,<br><div><br></div><div> </div></div></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Wed, Nov 6, 2019 at 11:13 AM Pierre Jolivet via petsc-dev <<a href="mailto:petsc-dev@mcs.anl.gov" target="_blank">petsc-dev@mcs.anl.gov</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">I need to figure this out myself first.<br>
In the meantime, here is another (PCASM related) question for you: why is PCASMType switched to BASIC when using a symmetric Pmat? (I guess I don’t have to tell you about the performance of RAS vs. ASM)<br>
To me, that would make sense if the default KSP was also switched to CG instead of GMRES if the Pmat and Amat were symmetric, but I don’t think this is the case.<br>
<br>
Thanks,<br>
Pierre<br>
<br>
> On 24 Oct 2019, at 5:40 PM, Smith, Barry F. <<a href="mailto:bsmith@mcs.anl.gov" target="_blank">bsmith@mcs.anl.gov</a>> wrote:<br>
> <br>
> <br>
>  Send the code and exact instructions to run a "good" and a "bad" ASM <br>
> <br>
>   Barry<br>
> <br>
> <br>
>> On Oct 14, 2019, at 10:44 AM, Pierre Jolivet <<a href="mailto:pierre.jolivet@enseeiht.fr" target="_blank">pierre.jolivet@enseeiht.fr</a>> wrote:<br>
>> <br>
>> Here are the three logs.<br>
>> FGMRES also gives a wrong first iterate.<br>
>> I think Mark was right in the sense that the problem is _most likely_ in my RHS.<br>
>> But I need to figure out why I only get this problem with right-preconditioned KSPs with restrict or none.<br>
>> <br>
>> Thanks,<br>
>> Pierre<br>
>> <br>
>> <br>
>> <br>
>>> On 13 Oct 2019, at 8:16 PM, Smith, Barry F. <<a href="mailto:bsmith@mcs.anl.gov" target="_blank">bsmith@mcs.anl.gov</a>> wrote:<br>
>>> <br>
>>> <br>
>>> Is this one process with one subdomain? (And hence no meaningful overlap since there is nothing to overlap?) And you expect to get the "exact" answer on one iteration? <br>
>>> <br>
>>> Please run the right preconditioned GMRES with -pc_asm_type [restrict and basic and none]  -ksp_monitor_true_solution  and send the output for the three cases.<br>
>>> <br>
>>> For kicks you can also try FGMRES (which always uses right preconditioning) to see if the same problem appears.<br>
>>> <br>
>>> Barry<br>
>>> <br>
>>> <br>
>>>> On Oct 13, 2019, at 2:41 AM, Pierre Jolivet via petsc-dev <<a href="mailto:petsc-dev@mcs.anl.gov" target="_blank">petsc-dev@mcs.anl.gov</a>> wrote:<br>
>>>> <br>
>>>> Hello,<br>
>>>> I’m struggling to understand the following weirdness with PCASM with exact subdomain solvers.<br>
>>>> I’m dealing with a very simple Poisson problem with Dirichlet + Neumann BCs.<br>
>>>> If I use PCBJACOBI + KSPPREONLY or 1 iteration of GMRES either preconditioned on the right or on the left, I get the expected result, cf. attached screenshot.<br>
>>>> If I use PCASM + KSPPREONLY or 1 iteration of GMRES preconditioned on the left, I get the expected result as well.<br>
>>>> However, with PCASM + 1 iteration of GMRES preconditioned on the right, I don’t get what I should (I believe).<br>
>>>> Furthermore, this problem is specific to -pc_asm_type restrict,none (I get the expected result with basic,interpolate).<br>
>>>> <br>
>>>> Any hint?<br>
>>>> <br>
>>>> Thanks,<br>
>>>> Pierre<br>
>>>> <br>
>>>> $ -sub_pc_type lu -ksp_max_it 1 -ksp_type gmres -pc_type bjacobi -ksp_pc_side right -> bjacobi_OK<br>
>>>> $ -sub_pc_type lu -ksp_max_it 1 -ksp_type gmres -pc_type asm -ksp_pc_side left -> asm_OK<br>
>>>> $ -sub_pc_type lu -ksp_max_it 1 -ksp_type gmres -pc_type asm -ksp_pc_side right -> asm_KO<br>
>>>> <br>
>>>> <bjacobi_OK.png><asm_OK.png><asm_KO.png><br>
>>> <br>
>> <br>
>> <dump-basic><dump-none><dump-restrict><br>
> <br>
<br>
</blockquote></div>
</div></blockquote></div><br></div></blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>