<div dir="ltr"><div dir="ltr">On Fri, Aug 8, 2025 at 12:22 PM Frank Bramkamp <<a href="mailto:bramkamp@nsc.liu.se">bramkamp@nsc.liu.se</a>> wrote:</div><div class="gmail_quote gmail_quote_container"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">Dear PETSc Team,<br>
<br>
I have the following general questions.<br>
When using GMRES with classic gram schmidt, you have an optional refinement.<br>
<br>
The default is not to use any refinement, is that correct ? </blockquote><div><br></div><div>Yes. It has been this way at least since 2013 from checking the log.</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">So you use the standard classic gram schmidt <br>
without any additional refinement step. <br>
<br>
I originally thought that you use one refinement step if needed as default. Or was that somehow changed one day<br>
and I remember it from older settings ?!<br></blockquote><div><br></div><div>It is available, but we changed it since refinement was almost never happening.</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">
With refinement one also needs about twice the MPI communication, if refinement is done in each iteration,<br>
which can be a bit more costly due to the MPI_ALLREDUCE for the dot products.<br></blockquote><div><br></div><div>It is reduction that is expensive on large runs.</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">
I solve the compressible navier stokes equations with turbulence. In which cases do you have the experience one should<br>
better use a refinement step ?! So far I typically used the default options, which seemed to work well for me so far.<br></blockquote><div><br></div><div>If it is converging, everything is likely fine. You can always check the true residual to detect breakdown, or turn on refinement if necessary with the option.</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">
for classic gram schmidt one also has to compute multiple dot products which is VecMDot in PETSC if I see this correct.<br></blockquote><div><br></div><div>Yes.</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">
Do you by standard use blas GEMV to compute multiple dot products at once or do you typically prefer your<br>
custom implementation. I think you also have a version where you unroll manual loops by a factor of 4 or so.<br>
I am not quite sure what is the typical default for this, or in which situation one would not use GEMV.<br></blockquote><div><br></div><div>It turns out that Junchao did extensive tests on this last year (I believe). By default, we use GEMV for MDot, but not for MAXPY (since performance here was worse with vendor GEMV). It is difficult (or impossible) to tell which algorithm will be better without testing, so there are options to select each one.</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">
Thanks a lot. I wish you a nice weekend, Frank<br>
<br>
<br>
<br>
<br>
<br>
<br>
<br>
<br>
<br>
<br>
</blockquote></div><div><br clear="all"></div><div><br></div><span class="gmail_signature_prefix">-- </span><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="https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!fhbzgnVYBw6ICZM6N0XXNLIcj9qe9bcopUK7lxSsqxqVlv-SvpS8bdrcwYq0TvTbHmbgbSinCLXrvJcM-AWO$" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>