[petsc-users] [EXTERNAL] Re: Question about PCFieldSplit

Patrick Sanan patrick.sanan at gmail.com
Mon Jan 31 09:47:23 CST 2022


The current behavior is that a single IS is returned for each stratum, so
if you have 2 unknowns on vertices, it should still return a single IS
including both of those unknowns per vertex. Am I understanding that that's
working as expected but you need *two* ISs in that case?


Am Mo., 31. Jan. 2022 um 16:29 Uhr schrieb Jorti, Zakariae <zjorti at lanl.gov
>:

> Hi Patrick,
>
>
> Thanks for your recent updates on DMStag.
> After getting the Finite Difference Coloring to work with our solver, I
> was testing that DMCreateFieldDecomposition routine that you added last
> week.
> It seems to work fine when there is only one unknown per location (i.e.
> one unknown on vertices and/or one unknown on edges and/or one unknown on
> faces and/or one unknown on elements).
> That being said, when there is more than one unknown in some location
> (let's say 2 unknowns on vertices for instance), I could not get the ISs
> for those two unknowns with that routine.
> Should I still rely on PCFieldSplitSetDetectSaddlePoint in this case?
> Many thanks once again.
>
>
> Kind regards,
>
>
> Zakariae
> ------------------------------
> *From:* Patrick Sanan <patrick.sanan at gmail.com>
> *Sent:* Tuesday, January 25, 2022 5:36:17 AM
> *To:* Matthew Knepley
> *Cc:* Jorti, Zakariae; petsc-users at mcs.anl.gov; Tang, Xianzhu; Dave May
> *Subject:* [EXTERNAL] Re: [petsc-users] Question about PCFieldSplit
>
> Here is an MR which intends to introduce some logic to support
> DMCreateFieldDecomposition(). It doesn't use the PetscSection approach,
> which might be preferable, but nonetheless is a necessary component so It'd
> like to get it in, even if it has room for further optimization. Hopefully
> this can be followed fairly soon with some more examples and tests using
> PCFieldSplit itself.
>
> https://gitlab.com/petsc/petsc/-/merge_requests/4740
>
> Am Mi., 23. Juni 2021 um 12:15 Uhr schrieb Matthew Knepley <
> knepley at gmail.com>:
>
>> On Wed, Jun 23, 2021 at 12:51 AM Patrick Sanan <patrick.sanan at gmail.com>
>> wrote:
>>
>>> Hi Zakariae -
>>>
>>> The usual way to do this is to define an IS (index set) with the degrees
>>> of freedom of interest for the rows, and another one for the columns, and
>>> then use MatCreateSubmatrix [1] .
>>>
>>> There's not a particularly convenient way to create an IS with the
>>> degrees of freedom corresponding to a particular "stratum" (i.e. elements,
>>> faces, edges, or vertices) of a DMStag, but fortunately I believe we have
>>> some code to do exactly this in a development branch.
>>>
>>> I'll track it down and see if it can quickly be added to the main branch.
>>>
>>
>> Note that an easy way to keep track of this would be to create a section
>> with the different locations as fields. This Section could then
>> easily create the ISes, and could automatically interface with
>> PCFIELDSPLIT.
>>
>>   Thanks,
>>
>>      Matt
>>
>>
>>>
>>> [1]:
>>> https://petsc.org/release/docs/manualpages/Mat/MatCreateSubMatrix.html
>>>
>>> Am 22.06.2021 um 22:29 schrieb Jorti, Zakariae <zjorti at lanl.gov>:
>>>
>>> Hello,
>>>
>>> I am working on DMStag and I have one dof on vertices (let us call
>>> it V), one dof on edges (let us call it E), one dof on faces ((let us
>>> call it F)) and one dof on cells (let us call it C).
>>> I build a matrix on this DM, and I was wondering if there was a way to
>>> get blocks (or sub matrices) of this matrix corresponding to specific
>>> degrees of freedom, for example rows corresponding to V dofs and columns
>>> corresponding to E dofs.
>>> I already asked this question before and the answer I got was I could
>>> call PCFieldSplitSetDetectSaddlePoint with the diagonal entries being
>>> of the matrix being zero or nonzero.
>>> That worked well. Nonetheless, I am curious to know if there
>>> was another alternative that does not require creating a dummy matrix
>>> with appropriate diagonal entries and solving a dummy linear system with
>>> this matrix to define the splits.
>>>
>>>
>>> Many thanks.
>>>
>>> Best regards,
>>>
>>> Zakariae
>>> ------------------------------
>>> *From:* petsc-users <petsc-users-bounces at mcs.anl.gov> on behalf of
>>> Tang, Qi <tangqi at msu.edu>
>>> *Sent:* Sunday, April 18, 2021 11:51:59 PM
>>> *To:* Patrick Sanan
>>> *Cc:* petsc-users at mcs.anl.gov; Tang, Xianzhu
>>> *Subject:* [EXTERNAL] Re: [petsc-users] Question about PCFieldSplit
>>>
>>> Thanks a lot, Patrick. We appreciate your help.
>>>
>>> Qi
>>>
>>>
>>>
>>> On Apr 18, 2021, at 11:30 PM, Patrick Sanan <patrick.sanan at gmail.com>
>>> wrote:
>>>
>>> We have this functionality in a branch, which I'm working on cleaning up
>>> to get to master. It doesn't use PETScSection. Sorry about the delay!
>>>
>>> You can only use PCFieldSplitSetDetectSaddlePoint when your diagonal
>>> entries being zero or non-zero defines the splits correctly.
>>>
>>> Am 17.04.2021 um 21:09 schrieb Matthew Knepley <knepley at gmail.com>:
>>>
>>> On Fri, Apr 16, 2021 at 8:39 PM Jorti, Zakariae via petsc-users <
>>> petsc-users at mcs.anl.gov> wrote:
>>>
>>>> Hello,
>>>>
>>>> I have a DMStag grid with one dof on each edge and face center.
>>>> I want to use a PCFieldSplit preconditioner on a Jacobian matrix that I
>>>> assume is already split but I am not sure how to determine the fields.
>>>> In the DMStag examples (ex2.c and ex3.c), the
>>>> function PCFieldSplitSetDetectSaddlePoint is used to determine those fields
>>>> based on zero diagonal entries. In my case, I have a Jacobian matrix that
>>>> does not have zero diagonal entries.
>>>> Can I use that PCFieldSplitSetDetectSaddlePoint in this case?
>>>> If not, how should I do?
>>>> Should I do like this example (
>>>> https://www.mcs.anl.gov/petsc/petsc-master/src/ksp/ksp/tutorials/ex43.c.html
>>>> <https://urldefense.com/v3/__https://www.mcs.anl.gov/petsc/petsc-master/src/ksp/ksp/tutorials/ex43.c.html__;!!HXCxUKc!jbBwV2h9luOW4dtBcNh6n_W1ULQnSVeXpxl0Ef1752s4Hlef-nC2JcnksFSO3Q$>
>>>> ):
>>>> const PetscInt Bfields[1] = {0},Efields[1] = {1};
>>>> KSPGetPC(ksp,&pc);
>>>> PCFieldSplitSetBlockSize(pc,2);
>>>> PCFieldSplitSetFields(pc,"B",1,Bfields,Bfields);
>>>> PCFieldSplitSetFields(pc,"E",1,Efields,Efields);
>>>> where my B unknowns are defined on face centers and E unknowns are
>>>> defined on edge centers?
>>>>
>>> That will not work.That interface only works for colocated fields that
>>> you get from DMDA.
>>>
>>> Patrick, does DMSTAG use PetscSection? Then the field split would be
>>> automatically calculated. If not, does it maintain the
>>> field division so that it could be given to PCFIELDSPLIT as ISes?
>>>
>>>   Thanks,
>>>
>>>      Matt
>>>
>>>> One last thing, I do not know which field comes first. Is it the one
>>>> defined for face dofs or edge dofs.
>>>>
>>>> Thank you.
>>>> Best regards,
>>>>
>>>> Zakariae
>>>>
>>>>
>>>
>>> --
>>> 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/
>>> <https://urldefense.com/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!HXCxUKc!jbBwV2h9luOW4dtBcNh6n_W1ULQnSVeXpxl0Ef1752s4Hlef-nC2JcmGgSwfag$>
>>>
>>>
>>>
>>
>> --
>> 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/20220131/cc7f9232/attachment.html>


More information about the petsc-users mailing list