[petsc-users] Correct way to set/track global numberings in DMPlex?
Matthew Knepley
knepley at gmail.com
Wed Apr 3 21:55:11 CDT 2024
On Wed, Apr 3, 2024 at 2:28 PM Ferrand, Jesus A. <FERRANJ2 at my.erau.edu>
wrote:
> Dear PETSc team: (Hoping to get a hold of Jed and/or Matt for this one)
> (Also, sorry for the mouthful below) I'm developing routines that will
> read/write CGNS files to DMPlex and vice versa. One of the recurring
> challenges is the bookkeeping
> ZjQcmQRYFpfptBannerStart
> This Message Is From an External Sender
> This message came from outside your organization.
>
> ZjQcmQRYFpfptBannerEnd
> Dear PETSc team:
>
> (Hoping to get a hold of Jed and/or Matt for this one)
> (Also, sorry for the mouthful below)
>
> I'm developing routines that will read/write CGNS files to DMPlex and vice
> versa.
> One of the recurring challenges is the bookkeeping of global numbering for
> vertices and cells.
> Currently, I am restricting my support to single Zone CGNS files, in which
> the file provides global numbers for vertices and cells.
>
I thought Jed had put in parallel CGNS loading. If so, maybe you can
transition to that. If not, we should get your stuff integrated.
> I used PetscHSetI as exemplified in DMPlexBuildFromCellListParallel() to
> obtain local DAG numbers from the global numbers provided by the file.
> I also used PetscSFCreateByMatchingIndices() to establish a basic DAG
> point distribution over the MPI processes.
> I use this PointSF to manually assemble a global PetscSection.
> For owned DAG points (per the PointSF) , I call
> PetscSectionSetOffset(section, point, file_offset);
> For ghost DAG points (per the PointSF) I call
> PetscSectionSetOffset(section, point, -(file_offset + 1));
>
This sounds alright to me, although I admit to not understanding exactly
what is being done.
All of what I have just described happens in my CGNS version of
> DMPlexTopologyLoad().
> My intention is to retain those numbers into the DMPlex, and reuse them in
> my CGNS analogues of DMPlexCoordinatesLoad(), DMPlexLabelsLoad(), and
> DMPlexGlobalVectorLoad().
> Anyhow, is this a good wait to track global numbers?
>
The way I think about it, for example in DMPlexBuildFromCellListParallel(),
you load all parallel things in blocks (cell adjacency, vertex coordinates,
etc). Then if you have to redistribute afterwards, you make a PetscSF to do
it. I first make one mapping points to points. With a PetscSection, you can
easily convert this into dofs to dofs. For example, we load sets of
vertices, but we want vertices distributed as they are attached to cells.
So we create a PetscSF mapping uniform blocks to the division attached to
cells. Then we use the PetscSection for coordinates to make a new PetscSF
and redistribute coordinates.
Thanks,
Matt
> Also, I need (for other applications) to eventually call
> DMPlexInterpolate() and DMPlexDistribute(), will the global PetscSection
> offsets be preserved after calling those two?
>
>
> Sincerely:
>
> *J.A. Ferrand*
>
> Embry-Riddle Aeronautical University - Daytona Beach - FL
> Ph.D. Candidate, Aerospace Engineering
>
> M.Sc. Aerospace Engineering
>
> B.Sc. Aerospace Engineering
>
> B.Sc. Computational Mathematics
>
>
> *Phone:* (386)-843-1829
>
> *Email(s):* ferranj2 at my.erau.edu
>
> jesus.ferrand at gmail.com
>
--
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!Z30NkeNY4hgcjCs5RtVH3AgiBI0E4BBkGDPYdLNB10LWOF050wW1AXJDMcOtZ0G3u9nPiKrc0MX9YIYzdUyK$ <https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!Z30NkeNY4hgcjCs5RtVH3AgiBI0E4BBkGDPYdLNB10LWOF050wW1AXJDMcOtZ0G3u9nPiKrc0MX9YIdWrfAJ$ >
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20240403/369284fd/attachment-0001.html>
More information about the petsc-users
mailing list