<html><head><meta http-equiv="Content-Type" content="text/html; charset=utf-8"></head><body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class=""><br class=""></div><blockquote type="cite" class=""><div class=""><div class=""><div dir="ltr" class=""><div class="gmail_quote"><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div dir="ltr" class=""><div id="gmail-m_-9053906504632661544divtagdefaultwrapper" dir="ltr" class="" style="font-size: 12pt; font-family: Calibri, Helvetica, sans-serif;"><div style="margin-top: 0px; margin-bottom: 0px;" class="">P is a diffusion matrix, which itself is inverted by <i class="">KSPCG</i>)</div></div></div></blockquote></div></div></div></div></blockquote><div class=""><br class=""></div>This worries me. Unless solved to full precision the action of solving with CG is not a linear operator in the input variable b, this means that the action of your Schur complement is not a linear operator and so iterative eigenvalue algorithms may not work correctly (even if the linear system seems to be converging). <div class=""><br class=""></div><div class="">I suggest you try a KSP tolerance of 10^-14 for the inner KSP solve when you attempt the eigenvalue computation.</div><div class=""><br class=""></div><div class="">Barry</div><div class=""><br class=""><div><br class=""><blockquote type="cite" class=""><div class="">On Nov 10, 2021, at 9:09 AM, Vladislav Pimanov <<a href="mailto:Vladislav.Pimanov@skoltech.ru" class="">Vladislav.Pimanov@skoltech.ru</a>> wrote:</div><br class="Apple-interchange-newline"><div class="">
<meta http-equiv="Content-Type" content="text/html; charset=koi8-r" class="">
<div class="">
<style type="text/css" style="display:none;" class=""><!-- P {margin-top:0;margin-bottom:0;} --></style>
<div id="divtagdefaultwrapper" style="font-size: 12pt; font-family: Calibri, Helvetica, sans-serif;" dir="ltr" class=""><p class="">Great thanks, Matt!</p><p class="">The second option is what I was looking for.</p><p class=""><br class="">
</p><p class="">Best Regards,</p><p class="">Vlad</p>
</div>
<hr style="display:inline-block;width:98%" tabindex="-1" class="">
<div id="divRplyFwdMsg" dir="ltr" class=""><font face="Calibri, sans-serif" style="font-size:11pt" class=""><b class="">От:</b> Matthew Knepley <<a href="mailto:knepley@gmail.com" class="">knepley@gmail.com</a>><br class="">
<b class="">Отправлено:</b> 10 ноября 2021 г. 16:45:00<br class="">
<b class="">Кому:</b> Vladislav Pimanov<br class="">
<b class="">Копия:</b> <a href="mailto:petsc-users@mcs.anl.gov" class="">petsc-users@mcs.anl.gov</a><br class="">
<b class="">Тема:</b> Re: [petsc-users] How to compute the condition number of SchurComplementMat preconditioned with PCSHELL.</font>
<div class=""> </div>
</div>
<div class="">
<div dir="ltr" class="">
<div dir="ltr" class="">On Wed, Nov 10, 2021 at 8:42 AM Vladislav Pimanov <<a href="mailto:Vladislav.Pimanov@skoltech.ru" class="">Vladislav.Pimanov@skoltech.ru</a>> wrote:<br class="">
</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 dir="ltr" class="">
<div id="gmail-m_-9053906504632661544divtagdefaultwrapper" style="font-size: 12pt; font-family: Calibri, Helvetica, sans-serif;" dir="ltr" class=""><p class="">Dear PETSc community,</p><p class=""><br class="">
</p><div class=""><br class="webkit-block-placeholder"></div>
<span class="">I wonder if you could give me a hint on how to compute the condition number of a preconditioned matrix in a proper way.</span><p class="">I have a <i class="">MatSchurComplement</i> matrix S and a preconditioner P of the type
<i class="">PCSHELL</i> (P is a diffusion matrix, which itself is inverted by <i class="">KSPCG</i>).</p><p class="">I tried to compute the condition number of P^{-1}S "for free" during the outer PCG procedure using <span class=""><i class="">KSPComputeExtremeSingularValues()</i> routine.</span></p><p class="">Unfortunately, \sigma_min does not converge even if the solution is computed with very high precision.</p><p class=""><span class="">I also looked at SLEPc interface, but did not realised how PC should be included.</span><br class="">
</p><div class=""><br class="webkit-block-placeholder"></div>
</div>
</div>
</blockquote>
<div class="">You can do this at least two ways:</div>
<div class=""><br class="">
</div>
<div class=""> 1) Make a MatShell for P^{-1} S. This is easy, but you will not be able to use any factorization-type PC on that matrix.</div>
<div class=""><br class="">
</div>
<div class=""> 2) Solve instead the generalized EVP, S x = \lambda P x. Since you already have P^{-1}, this should work well.</div>
<div class=""><br class="">
</div>
<div class=""> Thanks,</div>
<div class=""><br class="">
</div>
<div class=""> Matt</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" class="">
<div id="gmail-m_-9053906504632661544divtagdefaultwrapper" style="font-size: 12pt; font-family: Calibri, Helvetica, sans-serif;" dir="ltr" class=""><p class="">Thanks!</p><p class=""><br class="">
</p><p class="">Sincerely,</p><p class="">Vladislav Pimanov</p><p class=""><br class="">
</p>
</div>
</div>
</blockquote>
</div>
<br clear="all" class="">
<div class=""><br class="">
</div>
-- <br class="">
<div dir="ltr" class="gmail_signature">
<div dir="ltr" class="">
<div class="">
<div dir="ltr" class="">
<div class="">
<div dir="ltr" class="">
<div class="">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br class="">
-- Norbert Wiener</div>
<div class=""><br class="">
</div>
<div class=""><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank" class="">https://www.cse.buffalo.edu/~knepley/</a><br class="">
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div>
</div></blockquote></div><br class=""></div></body></html>