<div dir="ltr"><br><div class="gmail_extra"><br><div class="gmail_quote">On 25 July 2018 at 09:48, Matthew Knepley <span dir="ltr"><<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div class="gmail_quote"><span class="m_-3725184249599076478gmail-"><div dir="ltr">On Wed, Jul 25, 2018 at 4:24 AM Buesing, Henrik <<a href="mailto:hbuesing@eonerc.rwth-aachen.de" target="_blank">hbuesing@eonerc.rwth-aachen.d<wbr>e</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">Dear all, <br>
<br>
I would like to improve the iterative solver [1]. As I understand it I would need to improve the preconditioner for the Schur complement. <br>
<br>
How would I do that? <br></blockquote><div><br></div></span><div>1) I always start from the exact thing (full Schur factorization with exact solves and back off parts until I am happy)</div><div><br></div><div>2) Is the a11 block a good preconditioner for your Schur complement? If not, I would start by replacing that matrix</div><div>    with something better.</div></div></div></blockquote><div><br></div><div>Some additional info. If you want to pursue option 2, you need to do call <br></div><br></div><div class="gmail_quote">  PCFieldSplitSetSchurPre()<div><br></div><a href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCFieldSplitSetSchurPre.html#PCFieldSplitSetSchurPre" target="_blank">http://www.mcs.anl.gov/petsc/<wbr>petsc-current/docs/<wbr>manualpages/PC/<wbr>PCFieldSplitSetSchurPre.html#<wbr>PCFieldSplitSetSchurPre</a></div><div class="gmail_quote"><br></div><div class="gmail_quote">with <a href="http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCFieldSplitSchurPreType.html#PCFieldSplitSchurPreType" target="_blank">PC_FIELDSPLIT_SCHUR_PRE_USER</a>  (second arg) and your user defined schur complement preconditioner (last arg).</div>  <br><div class="gmail_quote">Thanks,</div><div class="gmail_quote">  Dave<br></div><div class="gmail_quote"><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 dir="ltr"><div class="gmail_quote"><div><br></div><div>  Thanks,</div><div><br></div><div>    Matt</div><span class="m_-3725184249599076478gmail-"><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">
Thank you for your help!<br>
Henrik<br>
<br>
<br>
<br>
[1]<br>
-ksp_max_it 100 -ksp_rtol 1e-6 -ksp_atol 1e-50 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_precondit<wbr>ion a11 -fieldsplit_p_w_ksp_type preonly -fieldsplit_S_n_ksp_type gmres -fieldsplit_p_w_pc_type hypre -fieldsplit_p_w_pc_hypre_type boomeramg -fieldsplit_S_n_pc_type hypre -fieldsplit_S_n_pc_hypre_type boomeramg -fieldsplit_S_n_ksp_max_it 100 fieldsplit_S_n_ksp_rtol 1e-2<br>
<br>
<br>
-- <br>
Dipl.-Math. Henrik Büsing<br>
Institute for Applied Geophysics and Geothermal Energy<br>
E.ON Energy Research Center<br>
RWTH Aachen University<br>
<br>
<a href="https://maps.google.com/?q=Mathieustr.+10&entry=gmail&source=g" target="_blank">Mathieustr. 10</a>        | Tel +49 (0)241 80 49907<br>
52074 Aachen, Germany | Fax +49 (0)241 80 49889<br>
<br>
<a href="http://www.eonerc.rwth-aachen.de/GGE" rel="noreferrer" target="_blank">http://www.eonerc.rwth-aachen.<wbr>de/GGE</a><br>
<a href="mailto:hbuesing@eonerc.rwth-aachen.de" target="_blank">hbuesing@eonerc.rwth-aachen.de</a><br>
<br>
</blockquote></span></div><span class="m_-3725184249599076478gmail-HOEnZb"><font color="#888888"><br clear="all"><div><br></div>-- <br><div dir="ltr" class="m_-3725184249599076478gmail-m_-2514012141401197864gmail_signature"><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.caam.rice.edu/~mk51/" target="_blank">https://www.cse.buffalo.edu/~k<wbr>nepley/</a><br></div></div></div></div></div></font></span></div>
</blockquote></div><br></div></div>