<div dir="ltr">Dear Matt and Barry,<div><br></div><div>Unfortunately DM is not used in our PETSc interface,which is essentially based on  the assembly matrix.<br>As Matt mentioned, the trick could be to define a DMShell in order to expose the splits we internally build based on IS. <br>If I understand well, I have to play around with DMShellSetCreateFieldDecomposition. Unfortunately, I do not find tests using this method. Am I missing something ?</div><div>Perhaps a good starting point would be to write a test like src/ksp/ksp/examples/tutorials/ex43.c in which one would consider the algebraic data  to be known (the matrix and the IS) and would set up a DMShell  (with an appropriate call to DMShellSetCreateFieldDecomposition).</div><div>Does that make sense?</div><div><br></div><div>Thanks,</div><div>Nicolas</div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">Le dim. 18 avr. 2021 à 15:54, Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>> a écrit :<br></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 dir="ltr">On Sat, Apr 17, 2021 at 6:13 PM Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank">bsmith@petsc.dev</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"><br>
  So you would like to be able to create three IS in your code and attach them with names to the PC.  Then have -pc_fieldsplit_XXX_fields be able to utilize the attached IS by name and use them to define the blocks.   <br>
<br>
  This is all doable and could be added to PCFIELDSPLIT without too much code, new code. The code would be largely like PCFieldSplitSetRuntimeSplits_Private.<br>
<br>
   The recursive part may also be doable but I think your syntax below is not exactly right. You would need something like <br>
<br>
-fieldsplit_0_pc_type fieldsplit   // split the first PC into a fieldsplit<br>
     -fieldsplit_0_pc_fieldsplit_0_fields xxx        -fieldsplit_0_fieldsplit_0_pc_type jacobi<br>
     -fieldsplit_0_pc_fieldsplit_1_fields yyy <br>
     etc<br>
<br>
this would split the first field into two fields and use jacobi on the first field.<br>
<br>
The problem is what to use for labelling the xxx and the yyy? <br>
<br>
I think one could achieve this by having the PCFIELDPLIT attach to each of its split PCs the appropriate modified IS with names attached to them.  There are two cases, <br>
<br>
  when building the split uses first all the entries from fieldsplit_v_, then from fieldsplit_p_ then the new ISs it needs to attach to the first split PC are two sets of integers the first from 0 to the len(v)-1 and the second from len(v) to len(v)+len(p)-1. <br>
<br>
  when building the split it interlaces the indices from v and p (interlacing only make sense if the size of v and p is the same). Then the new v would be {0,2,4,...} and p would be {1,3,...}. <br>
<br>
  If you are ambitious and would like to add this to fieldsplit.c we'd be very happy to receive an MR. It might even lead to allowing us to simply how the PCFIELDPLIT interacts with DMs. If all the split type, stride, named, etc are handle in a single consistent manner.<br></blockquote><div><br></div><div>Barry, this is already working with DMs, which I think is the right way to do this.</div><div><br></div><div>Here is the code:</div><div><br></div><div>   <a href="https://gitlab.com/petsc/petsc/-/blob/main/src/ksp/pc/impls/fieldsplit/fieldsplit.c#L420" target="_blank">https://gitlab.com/petsc/petsc/-/blob/main/src/ksp/pc/impls/fieldsplit/fieldsplit.c#L420</a><br></div><div><br></div><div>The DM must respond to DMCreateSubDM(). The interface is that the call provides a list of fields [f_0, f_1, ...]</div><div>and the DM returns an IS for that combination and a subdm for it. The subdm part is what allows you to do this</div><div>recursively, since you make the same call on the subdm.</div><div><br></div><div>If you use a PetscSection to keep track of the data layout, the code to do field selection is already done:</div><div><br></div><div>  <a href="https://gitlab.com/petsc/petsc/-/blob/main/src/dm/interface/dmi.c#L105" target="_blank">https://gitlab.com/petsc/petsc/-/blob/main/src/dm/interface/dmi.c#L105</a></div><div><br></div><div>which can just be called from the DMCreateSubDM() implementation.</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">
  Barry<br>
<br>
<br>
<br>
> On Apr 17, 2021, at 11:53 AM, Karin&NiKo <<a href="mailto:niko.karin@gmail.com" target="_blank">niko.karin@gmail.com</a>> wrote:<br>
> <br>
> Dear PETSc users,<br>
> <br>
> I use the fieldsplit PC in an application where the splits are programmatically defined by IS using PCFieldSplitSetIS. Then the user can specify its own PC at runtime using PETSc options.<br>
> My question : is it possible to define nested splits in this case as it can be done with strided splits (see snes/examples/tutorials/ex19.c with test suffix fieldsplit_4).<br>
> <br>
> In order to be perfectly clear : let's say I have a 3 fields problem : velocity (v split), pressure (p split) and temperature (t split).<br>
> What I would like to do is something like the following but it fails : <br>
> <br>
> -ksp_type fgmres<br>
> -pc_fieldsplit_type multiplicative <br>
> -pc_type fieldsplit    -pc_fieldsplit_0_fields fieldsplit_v_, fieldsplit_p_    -pc_fieldsplit_1_fields fieldsplit_t_ <br>
> <br>
> -prefix_push  fieldsplit_0_  <br>
> -ksp_type fgmres<br>
> -pc_fieldsplit_schur_factorization_type upper <br>
> -pc_type fieldsplit <br>
> -pc_fieldsplit_type schur <br>
> -pc_fieldsplit_schur_precondition a11 <br>
> -prefix_pop<br>
> <br>
> -prefix_push  fieldsplit_1_  <br>
> -ksp_type fgmres<br>
> -pc_type jacobi<br>
> -prefix_pop<br>
> <br>
> -prefix_push  fieldsplit_v_  <br>
> -ksp_type fgmres<br>
> -pc_type gamg<br>
> -prefix_pop<br>
> <br>
> -prefix_push  fieldsplit_p_  <br>
> -ksp_type fgmres<br>
> -pc_type jacobi<br>
> -prefix_pop<br>
> <br>
> I thank you for your help,<br>
> Nicolas<br>
>  <br>
> <br>
<br>
</blockquote></div><br clear="all"><div><br></div>-- <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>
</blockquote></div>