<html><head><meta http-equiv="Content-Type" content="text/html; charset=us-ascii"></head><body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><div class=""><br class=""></div>  Schur complements are generally dense, unless your original matrix has specific structure that you know of, it is unlikely hypre AMG will work well on the Schur complement.<div class=""><br class=""></div><div class="">  Barry</div><div class=""><br class=""><div><br class=""><blockquote type="cite" class=""><div class="">On Apr 20, 2021, at 5:26 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com" class="">knepley@gmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><div class=""><div dir="ltr" class=""><div dir="ltr" class="">On Tue, Apr 20, 2021 at 4:33 PM Jorti, Zakariae via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" class="">petsc-users@mcs.anl.gov</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_-1889263841708740911divtagdefaultwrapper" style="font-size: 12pt; font-family: Calibri, Helvetica, sans-serif;" dir="ltr" class=""><p class="">Hello,</p><p class=""><br class="">
</p><p class="">I have a block matrix A,</p><p class=""><br class="">
</p><p class="">[A_{00}   A_{01}]</p><p class="">[A_{10}   A_{11}] , </p><p class=""><br class="">
</p><p class="">and <span style="font-size:12pt" class="">a</span><span style="font-size:12pt" class="">m trying to use PCFIELDSPLIT as a preconditioner on</span><span style="font-size:12pt" class=""> it. The preconditioner has this shape: </span></p><p class=""><span style="font-size:12pt" class=""><br class="">
</span></p><p class="">[inv(A_{00})       0     ]</p><p class="">[  0                -ksp(S)] ;</p><p class=""><br class="">
</p><p class="">where S is the Schur Complement matrix: </p><p class="">S = <span style="font-size:12pt" class="">A_{11} - </span><span style="font-size:12pt" class="">A_{10} * inv(</span><span style="font-size:12pt" class="">A_{00}</span><span style="font-size:12pt" class="">) * </span><span style="font-size:12pt" class="">A_{01},</span></p><p class=""><span style="font-size:12pt" class="">that I get by calling: </span></p><p class=""><br class="">
</p><p class="">MatGetSchurComplement(A, isrow0, iscol0, isrow1, iscol1, MAT_INITIAL_MATRIX, &S, MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_IGNORE_MATRIX, NULL);</p>
<div class=""><br class="">
</div><p class=""><span style="font-size:12pt" class="">I would like to set BoomerAMG preconditioner for </span>ksp(S). So, I am calling:</p><div class=""> <br class="webkit-block-placeholder"></div><p class="">PCFieldSplitGetSubKSP(pc, &n, &subksp);</p><p class="">KSPSetOperators(subksp[1], S, S);</p><p class="">KSPGetPC(<span style="font-size:12pt" class="">subksp[1],</span><span style="font-size:12pt" class=""> subpc</span><span style="font-size:12pt" class="">);</span></p><p class=""><span style="font-size:12pt" class="">PCSetType(subpc,"boomeramg");</span></p><p class=""><span style="font-size:12pt" class=""><br class="">
</span></p><p class=""><span style="font-size:12pt" class="">but I end up getting the following error:</span></p><p class=""><span style="font-size:12pt" class=""><br class="">
</span></p><p class="">[0]PETSC ERROR: No support for this operation for this object type</p><p class="">[0]PETSC ERROR: Mat type schurcomplement</p><p class="">[0]PETSC ERROR: See <a href="https://www.mcs.anl.gov/petsc/documentation/faq.html" target="_blank" class="">https://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.</p><p class="">[0]PETSC ERROR: Petsc Release Version 3.14.2, Dec 03, 2020 </p><p class="">[0]PETSC ERROR: ./mimeticcurleuler2 on a arch-linux-c-debug named nid00140 by zjorti Tue Apr 20 12:23:34 2021</p><p class="">[0]PETSC ERROR: Configure options CC=cc CXX=CC FC=ftn COPTFLAGS=-g CXXOPTFLAGS=-g FOPTFLAGS=-g --download-hypre --with-debugging=1 --with-cxx-dialect=C++11</p><p class="">[0]PETSC ERROR: #1 MatGetRow() line 561 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/mat/interface/matrix.c</p><p class="">[0]PETSC ERROR: #2 MatConvert_Basic() line 53 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/mat/utils/convert.c</p><p class="">[0]PETSC ERROR: #3 MatConvert() line 4254 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/mat/interface/matrix.c</p><p class="">[0]PETSC ERROR: #4 PCSetUp_HYPRE() line 231 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ksp/pc/impls/hypre/hypre.c</p><p class="">[0]PETSC ERROR: #5 PCSetUp() line 1009 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ksp/pc/interface/precon.c</p><p class="">[0]PETSC ERROR: #6 KSPSetUp() line 406 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ksp/ksp/interface/itfunc.c</p><p class="">[0]PETSC ERROR: #7 KSPSolve_Private() line 658 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ksp/ksp/interface/itfunc.c</p><p class="">[0]PETSC ERROR: #8 KSPSolve() line 889 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ksp/ksp/interface/itfunc.c</p><p class="">[0]PETSC ERROR: #9 PCApply_FieldSplit_Schur() line 1133 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ksp/pc/impls/fieldsplit/fieldsplit.c</p><p class="">[0]PETSC ERROR: #10 PCApply() line 444 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ksp/pc/interface/precon.c</p><p class="">[0]PETSC ERROR: #11 KSP_PCApply() line 283 in /global/u1/z/zjorti/Software/petsc-3.14.2/include/petsc/private/kspimpl.h</p><p class="">[0]PETSC ERROR: #12 KSPInitialResidual() line 65 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ksp/ksp/interface/itres.c</p><p class="">[0]PETSC ERROR: #13 KSPSolve_GMRES() line 245 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ksp/ksp/impls/gmres/gmres.c</p><p class="">[0]PETSC ERROR: #14 KSPSolve_Private() line 719 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ksp/ksp/interface/itfunc.c</p><p class="">[0]PETSC ERROR: #15 KSPSolve() line 889 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ksp/ksp/interface/itfunc.c</p><p class="">[0]PETSC ERROR: #16 SNESSolve_KSPONLY() line 51 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/snes/impls/ksponly/ksponly.c</p><p class="">[0]PETSC ERROR: #17 SNESSolve() line 4569 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/snes/interface/snes.c</p><p class="">[0]PETSC ERROR: #18 TSTheta_SNESSolve() line 185 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ts/impls/implicit/theta/theta.c</p><p class="">[0]PETSC ERROR: #19 TSStep_Theta() line 223 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ts/impls/implicit/theta/theta.c</p><p class="">[0]PETSC ERROR: #20 TSStep() line 3757 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ts/interface/ts.c</p><p class="">[0]PETSC ERROR: #21 TSSolve() line 4154 in /global/u1/z/zjorti/Software/petsc-3.14.2/src/ts/interface/ts.c</p><p class=""><span style="font-size:12pt" class=""><br class="">
</span></p><p class=""><span style="font-size:12pt" class="">Am I doing something </span>wrong here? How can I fix this?</p></div></div></blockquote><div class="">It really is best to do this with options. It will take a lot of programming to replicate what you can do.</div><div class=""><br class=""></div><div class="">You are trying to precondition an implicit matrix (S) with Boomeramg, which needs explicit entries. You can tell PCFIELDSPLIT to build an</div><div class="">explicit matrix for the Schur complement:</div><div class=""><br class=""></div><div class="">  <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCFieldSplitSetSchurPre.html" class="">https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCFieldSplitSetSchurPre.html</a></div><div class=""><br class=""></div><div class="">  Thanks,</div><div class=""><br class=""></div><div class="">     Matt <br class=""></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_-1889263841708740911divtagdefaultwrapper" style="font-size: 12pt; font-family: Calibri, Helvetica, sans-serif;" dir="ltr" class=""><p class="">Many thanks.</p><p class=""><br class="">
</p><p class="">Best regards,</p><p class=""><br class="">
</p><p class="">Zakariae Jorti</p><p class=""><span style="font-size:12pt" class=""><br class="">
</span></p><div class=""><span style="font-size:12pt" class=""></span><span style="font-size:12pt" class=""> </span><br class="webkit-block-placeholder"></div>
</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></blockquote></div><br class=""></div></body></html>