[petsc-users] Nesting splits based on IS

Matthew Knepley knepley at gmail.com
Tue Apr 27 09:55:23 CDT 2021


On Tue, Apr 27, 2021 at 9:58 AM Karin&NiKo <niko.karin at gmail.com> wrote:

> I will have a look at PetscSection (that I have never used before).
> The FE I use are P2-P1-P1-P1 Lagrange elements. The point is that the
> numbering of the unknowns is done by the application (I do not have the
> hand on it). At the moment, I build IS in order to define the different
> fields and I nest them programmatically in order to set some recursive
> preconditioner.
> To enable the use of the command line options in order to nest the splits
> is a matter of elegance and smoothness since the nesting of splits I have
> programmed are fully fixed.
>

Okay, PetscSection should be completely adequate for that. You just

1) Create a Section with 4 fields and give 1 dof/vertex for each field, and
1dof/edge for the first field (have to also set the cumulative dofs per
point)

2) SectionSetUp()

3) Add a permutation if necessary to the points to get the dof ordering you
want

Does this make sense?

  Thanks,

     Matt


> Nicolas
>
> Le mar. 27 avr. 2021 à 14:15, Matthew Knepley <knepley at gmail.com> a
> écrit :
>
>> On Tue, Apr 27, 2021 at 2:51 AM Karin&NiKo <niko.karin at gmail.com> wrote:
>>
>>> Dear PETSc Team,
>>>
>>> I am coming back to you to get some clue on the (missing?) Fortran
>>> interface to the DMShellSetCreateSubDM routine.
>>>
>>
>> Yes, we will have to put that in. Fortran bindings for callbacks are
>> somewhat involved.
>>
>> However, it is easier if you layout can be described by a
>> PetscSection (which has all the Fortran bindings). Then it will
>> work automatically from the DMShell (I think). What kind of layout were
>> you looking for?
>>
>>   Thanks,
>>
>>      Matt
>>
>>
>>> Thanks,
>>> Nicolas
>>>
>>> Le dim. 25 avr. 2021 à 18:20, Karin&NiKo <niko.karin at gmail.com> a
>>> écrit :
>>>
>>>> Dear Petsc Team,
>>>>
>>>> I have begun to write a small test to illustrate the use of a DMShell
>>>> in order to nest splits from the command line.
>>>> I am writing it in Fortran because the targeted application is also in
>>>> Fortran.
>>>> => I am afraid that DMShellSetCreateSubDM lacks a Fortran interface.
>>>> Am I wrong ? Can't the interface be generated by bfort ?
>>>>
>>>> Thanks,
>>>> Nicolas
>>>>
>>>> Le lun. 19 avr. 2021 à 12:28, Matthew Knepley <knepley at gmail.com> a
>>>> écrit :
>>>>
>>>>> On Sun, Apr 18, 2021 at 11:54 PM Barry Smith <bsmith at petsc.dev> wrote:
>>>>>
>>>>>>
>>>>>>   Matt,
>>>>>>
>>>>>>    I agree it can be done with DMs, but I am not convinced that is
>>>>>> necessarily a superior way when one has to make DMSHELLS and mess around
>>>>>> with their behavior.
>>>>>>
>>>>>>    One way to look at PCFIELDSPLIT is that its basic splitting tool
>>>>>> is IS and DM is a way of generating IS that can be used for PCFIELDSPLIT.
>>>>>> In this approach PCFIELDSPLIT needs a bit more support for IS and nesting.
>>>>>> To also help using multiple levels of DM.
>>>>>>
>>>>>>    The other way is to view DM and subDM as the basic splitting tool
>>>>>> for PCFIELDSPLIT but then one has to ask the question how does DM
>>>>>> communicate its splits to PCFIELDSPLIT?   I am too lazy to look at the code
>>>>>> but presumably with IS, hence the approach above.
>>>>>>
>>>>>>    In the traditional Unix software stack approach to code, each
>>>>>> layer uses the simpler layer below to communicate; so DM uses IS to
>>>>>> communicate to the layer below (the linear algebra and preconditioners).
>>>>>> DM is a hugely complicated layer and making it communicate directly with
>>>>>> the vector level is complex; I like having the IS level in between to
>>>>>> simplify the software stack and programming model.
>>>>>>
>>>>>>    PetscSection makes live more complicated since it is a bit
>>>>>> disjoint from IS. I've never resolved in my mind what role PetscSection
>>>>>> plays, is it a level above IS in the software stack that generates IS, does
>>>>>> it sometimes "skip over" IS to directly talk to linear algebra?
>>>>>>
>>>>>>    If you cannot make a cartoon picture of the software stack, with
>>>>>> all the objects, then I think the software stack is not well designed or
>>>>>> defined. I fear we cannot make such a cartoon currently.
>>>>>>
>>>>>
>>>>> DM _does_ communicate with PCFIELDSPLIT using IS. I agree with you
>>>>> that IS is a good object for communication. In PETSc, IS is just a nice way
>>>>> to pass a list of integers.
>>>>>
>>>>> I don't think DM is hugely complicated. It does a few simple jobs.
>>>>> Here it's job is to remember which field each dof belongs to. That is all
>>>>> we have to worry about.
>>>>>
>>>>> PetscSection semantically is a linear space with some structure. We
>>>>> already know we want some structure like this, since we break all linear
>>>>> spaces in PETSc into processes. Section allows you to break it down a
>>>>> little finer into the pieces of the space for each "point", where you can
>>>>> use a point to mark anything you want, like a process, cell, edge, another
>>>>> dof, etc. Sections can make an IS when asked a question, such as "which
>>>>> dofs lie in the closure of this cell", or "which dofs are in this field",
>>>>> or "which dofs are owned by this process". I have written this in the
>>>>> manual years ago.
>>>>>
>>>>>     Matt
>>>>>
>>>>>
>>>>>>  Barry
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>> On Apr 18, 2021, at 8:54 AM, Matthew Knepley <knepley at gmail.com>
>>>>>> wrote:
>>>>>>
>>>>>> On Sat, Apr 17, 2021 at 6:13 PM Barry Smith <bsmith at petsc.dev> wrote:
>>>>>>
>>>>>>>
>>>>>>>   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.
>>>>>>>
>>>>>>>
>>>>>>>   This is all doable and could be added to PCFIELDSPLIT without too
>>>>>>> much code, new code. The code would be largely like
>>>>>>> PCFieldSplitSetRuntimeSplits_Private.
>>>>>>>
>>>>>>>    The recursive part may also be doable but I think your syntax
>>>>>>> below is not exactly right. You would need something like
>>>>>>>
>>>>>>> -fieldsplit_0_pc_type fieldsplit   // split the first PC into a
>>>>>>> fieldsplit
>>>>>>>      -fieldsplit_0_pc_fieldsplit_0_fields xxx
>>>>>>> -fieldsplit_0_fieldsplit_0_pc_type jacobi
>>>>>>>      -fieldsplit_0_pc_fieldsplit_1_fields yyy
>>>>>>>      etc
>>>>>>>
>>>>>>> this would split the first field into two fields and use jacobi on
>>>>>>> the first field.
>>>>>>>
>>>>>>> The problem is what to use for labelling the xxx and the yyy?
>>>>>>>
>>>>>>> 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,
>>>>>>>
>>>>>>>   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.
>>>>>>>
>>>>>>>   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,...}.
>>>>>>>
>>>>>>>   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.
>>>>>>>
>>>>>>
>>>>>> Barry, this is already working with DMs, which I think is the right
>>>>>> way to do this.
>>>>>>
>>>>>> Here is the code:
>>>>>>
>>>>>>
>>>>>> https://gitlab.com/petsc/petsc/-/blob/main/src/ksp/pc/impls/fieldsplit/fieldsplit.c#L420
>>>>>>
>>>>>> The DM must respond to DMCreateSubDM(). The interface is that the
>>>>>> call provides a list of fields [f_0, f_1, ...]
>>>>>> and the DM returns an IS for that combination and a subdm for it. The
>>>>>> subdm part is what allows you to do this
>>>>>> recursively, since you make the same call on the subdm.
>>>>>>
>>>>>> If you use a PetscSection to keep track of the data layout, the code
>>>>>> to do field selection is already done:
>>>>>>
>>>>>>
>>>>>> https://gitlab.com/petsc/petsc/-/blob/main/src/dm/interface/dmi.c#L105
>>>>>>
>>>>>> which can just be called from the DMCreateSubDM() implementation.
>>>>>>
>>>>>>   Thanks,
>>>>>>
>>>>>>      Matt
>>>>>>
>>>>>>
>>>>>>>   Barry
>>>>>>>
>>>>>>>
>>>>>>>
>>>>>>> > On Apr 17, 2021, at 11:53 AM, Karin&NiKo <niko.karin at gmail.com>
>>>>>>> wrote:
>>>>>>> >
>>>>>>> > Dear PETSc users,
>>>>>>> >
>>>>>>> > 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.
>>>>>>> > 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).
>>>>>>> >
>>>>>>> > 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).
>>>>>>> > What I would like to do is something like the following but it
>>>>>>> fails :
>>>>>>> >
>>>>>>> > -ksp_type fgmres
>>>>>>> > -pc_fieldsplit_type multiplicative
>>>>>>> > -pc_type fieldsplit    -pc_fieldsplit_0_fields fieldsplit_v_,
>>>>>>> fieldsplit_p_    -pc_fieldsplit_1_fields fieldsplit_t_
>>>>>>> >
>>>>>>> > -prefix_push  fieldsplit_0_
>>>>>>> > -ksp_type fgmres
>>>>>>> > -pc_fieldsplit_schur_factorization_type upper
>>>>>>> > -pc_type fieldsplit
>>>>>>> > -pc_fieldsplit_type schur
>>>>>>> > -pc_fieldsplit_schur_precondition a11
>>>>>>> > -prefix_pop
>>>>>>> >
>>>>>>> > -prefix_push  fieldsplit_1_
>>>>>>> > -ksp_type fgmres
>>>>>>> > -pc_type jacobi
>>>>>>> > -prefix_pop
>>>>>>> >
>>>>>>> > -prefix_push  fieldsplit_v_
>>>>>>> > -ksp_type fgmres
>>>>>>> > -pc_type gamg
>>>>>>> > -prefix_pop
>>>>>>> >
>>>>>>> > -prefix_push  fieldsplit_p_
>>>>>>> > -ksp_type fgmres
>>>>>>> > -pc_type jacobi
>>>>>>> > -prefix_pop
>>>>>>> >
>>>>>>> > I thank you for your help,
>>>>>>> > Nicolas
>>>>>>> >
>>>>>>> >
>>>>>>>
>>>>>>>
>>>>>>
>>>>>> --
>>>>>> What most experimenters take for granted before they begin their
>>>>>> experiments is infinitely more interesting than any results to which their
>>>>>> experiments lead.
>>>>>> -- Norbert Wiener
>>>>>>
>>>>>> https://www.cse.buffalo.edu/~knepley/
>>>>>> <http://www.cse.buffalo.edu/~knepley/>
>>>>>>
>>>>>>
>>>>>>
>>>>>
>>>>> --
>>>>> What most experimenters take for granted before they begin their
>>>>> experiments is infinitely more interesting than any results to which their
>>>>> experiments lead.
>>>>> -- Norbert Wiener
>>>>>
>>>>> https://www.cse.buffalo.edu/~knepley/
>>>>> <http://www.cse.buffalo.edu/~knepley/>
>>>>>
>>>>
>>
>> --
>> What most experimenters take for granted before they begin their
>> experiments is infinitely more interesting than any results to which their
>> experiments lead.
>> -- Norbert Wiener
>>
>> https://www.cse.buffalo.edu/~knepley/
>> <http://www.cse.buffalo.edu/~knepley/>
>>
>

-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener

https://www.cse.buffalo.edu/~knepley/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20210427/c511e0e6/attachment-0001.html>


More information about the petsc-users mailing list