[petsc-users] Element connectivity of a DMPlex
Matthew Knepley
knepley at gmail.com
Wed May 21 07:50:05 CDT 2025
On Mon, May 19, 2025 at 6:57 PM Noam T. via petsc-users <
petsc-users at mcs.anl.gov> wrote:
> Hello,
>
> I am trying to build the connectivity matrix for a mesh; i.e. the indices
> of the nodes that compose each cell.
>
> Example:
>
> 0-------1
> |.\.....|
> |...\...|
> |.....\.|
> 3-------2
>
> the matrix would look like [ [0, 3, 2], [2, 1, 0] ] (possibly different
> ordering).
>
> One option is using DMPlexGetTransitiveClosure, and accessing the last
> elements in the output "points", which contain the vertex indices.
> However, these indices are local per process (both for vertices and
> cells). Is it possible to get the global indices?
>
> I tried a mapping, as in example 14f (
> https://urldefense.us/v3/__https://petsc.org/release/src/ksp/ksp/tutorials/ex14f.F90.html__;!!G_uCfscf7eWS!cZBA-cU24g0I7WDWgJbyMqwJKxk50DrkIGAu36vgS0MNYlEAzNA354o5U-YVkcur8NsRsoSXABuuOsuln5el$
> <https://urldefense.us/v3/__https://petsc.org/release/src/ksp/ksp/tutorials/ex14f.F90.html__;!!G_uCfscf7eWS!ezMKlDED-8WxXpg_Elxus7WhYkZzCA5OmmLe6vHZGFoj4se5S5RD5ZgTrC3TYErx8znYKebSW-s0H0mt1_MEMAzpmFDwr0PH$>
> )
>
> ---
> DM :: dm
> ISLocalToGlobalMapping :: ltog_map
> PetscInt, pointer :: ltog_idx(:)
>
> ! ...
> ! Create a dm from a mesh file
> ! ...
> DMGetLocalToGlobalMapping(dm,ltog_map,ierr)
> ISLocalToGlobalMappingGetIndices(ltog_map, ltog_idx, ierr)
> ---
>
> but the returned array ltog is empty (ISLocalToGlobalMappingGetSize()
> returns zero). Are there other functions calls needed before being able
> to create this mapping? Or is this mapping simply not usable in this case?
>
> Is there perhaps better/simpler way to get this connectivity?
>
There are several ways you might do this. It helps to know what you are
aiming for.
1) If you just want this output, it might be easier to just DMView() with
the PETSC_VIEWER_HDF5_VIZ format, since that just outputs the cell-vertex
topology and coordinates
2) You can call DMPlexUninterpolate() to produce a mesh with just cells and
vertices, and output it in any format.
3) If you want it in memory, but still with global indices (I don't
understand this use case), then you can use DMPlexCreatePointNumbering()
for an overall global numbering, or DMPlexCreateCellNumbering() and
DMPlexCreateVertexNumbering() for separate global numberings.
Thanks,
Matt
> Thank you.
> Noam
>
>
--
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!cZBA-cU24g0I7WDWgJbyMqwJKxk50DrkIGAu36vgS0MNYlEAzNA354o5U-YVkcur8NsRsoSXABuuOvSIT50e$ <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!cZBA-cU24g0I7WDWgJbyMqwJKxk50DrkIGAu36vgS0MNYlEAzNA354o5U-YVkcur8NsRsoSXABuuOt6MQSJx$ >
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20250521/80798a69/attachment.html>
More information about the petsc-users
mailing list