[petsc-users] Mat indices for DMPlex jacobian

Matthew Knepley knepley at gmail.com
Mon Sep 23 10:16:12 CDT 2024


On Mon, Sep 23, 2024 at 5:33 AM Matteo Semplice via petsc-users <
petsc-users at mcs.anl.gov> wrote:

> Dear petsc,
>
>      I need to hand-code a jacobian and I can't figure out how to
> translate the DMplex points/fields to matrix indices.
>
> My DMPlex has a section with m fields per cell (which makes for n>m dof
> per cell since some fields are vector). Say I want to insert an nxn
> block for the row corresponding to cell c and coloumn of its neighbour
> d. I guess that I should call either MatSetValues/MatSetValuesLocal or
> the blocked variants, but how do I find the row/col indices to pass in
> starting from the c/d dmplex points? And, while I am at it, which
> MatSetValues version standard/Local standard/blocked is best?
>
> I looked for a petsc example, but fails to find it: if there's one, can
> you just point me to it?
>

1. In FEM, the main unit of indices is usually the closure, so we provide

  https://urldefense.us/v3/__https://petsc.org/main/manualpages/DMPlex/DMPlexGetClosureIndices/__;!!G_uCfscf7eWS!fDvGoz5OloH7xryMIAIqyYsslt8U4HqegW7HeYDlwRPbkTloMU6k0kze8NlG7S2DxDuY2nUXRQcBIiuqLU1o$ 

which is what is used inside of

  https://urldefense.us/v3/__https://petsc.org/main/manualpages/DMPlex/DMPlexMatSetClosure/__;!!G_uCfscf7eWS!fDvGoz5OloH7xryMIAIqyYsslt8U4HqegW7HeYDlwRPbkTloMU6k0kze8NlG7S2DxDuY2nUXRQcBIi82N3yY$ 

2. These indices are calculated using the global section, so you can just

  DMGetGlobalSection(dm, &gs);

and then use PetscSectionGetDof() and PetscSectionGetOffset(), knowing that
off-process values are encoded as -(off + 1), so you need to convert those.

Does this make sense?

  Thanks,

     Matt


> Thanks in advance.
>
> Matteo Semplice
>
>

-- 
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://urldefense.us/v3/__https://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!fDvGoz5OloH7xryMIAIqyYsslt8U4HqegW7HeYDlwRPbkTloMU6k0kze8NlG7S2DxDuY2nUXRQcBIjjfDigI$  <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!fDvGoz5OloH7xryMIAIqyYsslt8U4HqegW7HeYDlwRPbkTloMU6k0kze8NlG7S2DxDuY2nUXRQcBIpAwlPEi$ >
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20240923/1b4fd704/attachment.html>


More information about the petsc-users mailing list