<div dir="ltr"><div dir="ltr">On Fri, Mar 13, 2020 at 11:20 AM Pierre Jolivet <<a href="mailto:pierre.jolivet@enseeiht.fr">pierre.jolivet@enseeiht.fr</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"><br><div><br><blockquote type="cite"><div>On 13 Mar 2020, at 2:36 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:</div><br><div><div dir="ltr" style="font-family:Helvetica;font-size:12px;font-style:normal;font-variant-caps:normal;font-weight:normal;letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none"><div dir="ltr">On Fri, Mar 13, 2020 at 3:19 AM Pierre Jolivet <<a href="mailto:pierre.jolivet@enseeiht.fr" target="_blank">pierre.jolivet@enseeiht.fr</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"><br><div><br><blockquote type="cite"><div>On 12 Mar 2020, at 11:40 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>> wrote:</div><br><div><div dir="ltr"><div dir="ltr">On Thu, Mar 12, 2020 at 5:59 PM Jed Brown <<a href="mailto:jed@jedbrown.org" target="_blank">jed@jedbrown.org</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">Pierre Jolivet <<a href="mailto:pierre.jolivet@enseeiht.fr" target="_blank">pierre.jolivet@enseeiht.fr</a>> writes:<br><br>> Hello,<br>> Has there been any follow-up on this<span> </span><a href="https://lists.mcs.anl.gov/pipermail/petsc-users/2015-January/024020.html" rel="noreferrer" target="_blank">https://lists.mcs.anl.gov/pipermail/petsc-users/2015-January/024020.html</a><span> </span><<a href="https://lists.mcs.anl.gov/pipermail/petsc-users/2015-January/024020.html" rel="noreferrer" target="_blank">https://lists.mcs.anl.gov/pipermail/petsc-users/2015-January/024020.html</a>>?<br>> Given a 3x3 MatNest A = [A_00,0,0 ; 0,A_11,0 ; 0,0,A_22], I’d like to setup a two-way fieldsplit coupling [A_00,0 ; 0,A_11] and [A_22] but I can’t figure out the proper options.<br><br>Are you looking for a Schur split or additive/multiplicative?</blockquote></div></div></div></blockquote><div><br></div><div>Don’t know yet which will perform best, do you have a specific solution in mind for one scenario or the other?</div><div>I was mostly wondering if it was possible in a general context, not taking -pc_fieldsplit_type into account.</div><br><blockquote type="cite"><div><div dir="ltr"><div class="gmail_quote"><div>-pc_fieldsplit_field_0 0,1 -pc_fieldsplit_field_1 2 -pc_fieldsplit_type schur</div></div></div></div></blockquote><div><br></div><div>These flags, used with my .cpp, yield:</div><div><div>[0]PETSC ERROR: Arguments are incompatible</div><div>[0]PETSC ERROR: To use Schur complement preconditioner you must have exactly 2 fields</div><div>If I use -pc_fieldsplit_%d_fields <a,b,..> as advocated in the manual (instead of -pc_fieldsplit_field_%d as you suggested), I get the same error.</div></div></div></div></blockquote><div><br></div><div>Okay it is  -pc_fieldsplit_%d_fields. When you use this, how many fields does it think you have?</div></div></div></div></blockquote><div><br></div><div>Three.</div><div>In case it’s not clear, the MWE is at the bottom of my first email (nest.cpp).</div><div>$ mpicxx nest.cpp -I$PETSC_DIR/$PETSC_ARCH/include -I$PETSC_DIR/include -L$PETSC_DIR/$PETSC_ARCH/lib -lpetsc</div><div>$ ./a.out -pc_type fieldsplit -pc_fieldsplit_0_fields 0,1 -pc_fieldsplit_1_fields 2 -pc_fieldsplit_type schur</div><div>[…]</div><div><div>[0]PETSC ERROR: Arguments are incompatible</div><div>[0]PETSC ERROR: To use Schur complement preconditioner you must have exactly 2 fields</div><div>[…]</div><div><div>    Split info:</div><div>    Split number 0 Defined by IS</div><div>    Split number 1 Defined by IS</div><div>    Split number 2 Defined by IS</div><div>[…]</div></div></div></div></div></blockquote><div><br></div><div>I see now. If you call PCFieldsplitSetIS(), this overrides anything else. We do not even try to discover the split. If you want</div><div>a split determined by command line arguments, you have to defer the split. It looks like it can get splits from the MatNest,</div><div>so just do not call SetIS().</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><div><div>Thanks,</div><div>Pierre</div></div></div><br><blockquote type="cite"><div><div dir="ltr" style="font-family:Helvetica;font-size:12px;font-style:normal;font-variant-caps:normal;font-weight:normal;letter-spacing:normal;text-align:start;text-indent:0px;text-transform:none;white-space:normal;word-spacing:0px;text-decoration:none"><div class="gmail_quote"><div>   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 style="word-wrap:break-word"><div><blockquote type="cite"><div><div dir="ltr"><div class="gmail_quote"><div>I believe.</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">> Jed, in this answer<span> </span><a href="https://lists.mcs.anl.gov/pipermail/petsc-users/2015-January/023993.html" rel="noreferrer" target="_blank">https://lists.mcs.anl.gov/pipermail/petsc-users/2015-January/023993.html</a><<a href="https://lists.mcs.anl.gov/pipermail/petsc-users/2015-January/023993.html" rel="noreferrer" target="_blank">https://lists.mcs.anl.gov/pipermail/petsc-users/2015-January/023993.html</a>>, you recommend not to use MatNest. What would you recommend instead?<br><br>See src/snes/examples/tutorials/ex28.c for my preferred approach.</blockquote></div></div></div></blockquote><div><br></div><div>I’m not using DM, so in terms of Mat, I guess this example shows that I pretty much need to reimplement the field splitting myself?</div><div><br></div><div>Thanks,</div><div>Pierre</div><br><blockquote type="cite"><div><div dir="ltr">--<span> </span><br><div dir="ltr"><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></div></blockquote></div><br></div></blockquote></div><br clear="all"><div><br></div>--<span> </span><br><div dir="ltr"><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></div></div></div></div></div></div></div></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>