[petsc-users] DMGetCoordinatesLocal and DMPlexGetCellCoordinates in PETSc > 3.18

Matthew Knepley knepley at gmail.com
Mon May 22 05:10:08 CDT 2023


On Mon, May 22, 2023 at 4:41 AM Berend van Wachem <berend.vanwachem at ovgu.de>
wrote:

> Dear Matt,
>
> I'm really sorry for this stupid bug!
>

No problem. You have really helped me get the bugs out of Plex.

  Thanks,

     Matt


> I can confirm that setting the coordinates with both
> CellCoordinatesLocal and CoordinatesLocal works.
>
> Best regards, Berend.
>
> Many thanks and best regards, Berend.
>
> On 5/17/23 23:04, Matthew Knepley wrote:
> > On Wed, May 17, 2023 at 2:01 PM Berend van Wachem
> > <berend.vanwachem at ovgu.de <mailto:berend.vanwachem at ovgu.de>> wrote:
> >
> >     Dear Matt,
> >
> >     I tried it, but it doesn't seem to work.
> >     Attached is a very small working example illustrating the problem.
> >     I create a DMPlexBox Mesh, periodic in the Y direction. I then scale
> >     the Y coordinates with a factor 10, and add 1.0 to it. Both
> >     DMGetCoordinatesLocal and DMGetCellCoordinatesLocal.
> >     Then I evaluate the coordinates with DMPlexGetCellCoordinates. Most
> >     of the Y coordinates are correct, but not all of them - for
> >     instance, the minimum Y coordinate is 0.0, and this should be 1.0.
> >
> >     Am I doing something wrong?
> >
> >
> > Quickly, I see that
> >
> >    a *= 10.0 + 1.0;
> >
> > is the same as
> >
> >    a *= 11.0;
> >
> > not multiply by 10 and add 1. I will send it back when I get everything
> > the way I want.
> >
> >    Thanks,
> >
> >      Matt
> >
> >     Thanks and best regards,
> >
> >     Berend.
> >
> >     On 5/17/23 17:58, Matthew Knepley wrote:
> >      > On Wed, May 17, 2023 at 11:20 AM Berend van Wachem
> >     <berend.vanwachem at ovgu.de <mailto:berend.vanwachem at ovgu.de>
> >     <mailto:berend.vanwachem at ovgu.de <mailto:berend.vanwachem at ovgu.de>>>
> >     wrote:
> >      >
> >      >     Dear Matt,
> >      >
> >      >     Is there a way to 'redo' the DMLocalizeCoordinates() ? Or to
> >     undo it?
> >      >     Alternatively, can we make the calling of
> >     DMLocalizeCoordinates() in the  DMPlexCreate...() routines optional?
> >      >
> >      >     Otherwise, we would have to copy all arrays of coordinates
> >     from DMGetCoordinatesLocal() and DMGetCellCoordinatesLocal() before
> >      >     scaling them.
> >      >
> >      >
> >      > I am likely not being clear. I think all you have to do is the
> >     following:
> >      >
> >      >    DMGetCoordinatesLocal(dm, &xl);
> >      >    VecScale(xl, scale);
> >      >    DMSetCoordinatesLocal(dm, xl);
> >      >    DMGetCellCoordinatesLocal(dm, &xl);
> >      >    VecScale(xl, scale);
> >      >    DMSetCellCoordinatesLocal(dm, xl);
> >      >
> >      > Does this not work?
> >      >
> >      >    Thanks,
> >      >
> >      >       Matt
> >      >
> >      >     Best regards, Berend.
> >      >
> >      >     On 5/17/23 16:35, Matthew Knepley wrote:
> >      >      > On Wed, May 17, 2023 at 10:21 AM Berend van Wachem
> >     <berend.vanwachem at ovgu.de <mailto:berend.vanwachem at ovgu.de>
> >     <mailto:berend.vanwachem at ovgu.de <mailto:berend.vanwachem at ovgu.de>>
> >      >     <mailto:berend.vanwachem at ovgu.de
> >     <mailto:berend.vanwachem at ovgu.de> <mailto:berend.vanwachem at ovgu.de
> >     <mailto:berend.vanwachem at ovgu.de>>>> wrote:
> >      >      >
> >      >      >     Dear Matt,
> >      >      >
> >      >      >     Thanks for getting back to me so quickly.
> >      >      >
> >      >      >     If I scale each of the coordinates of the mesh (say, I
> >     want to cube each
> >      >      >     co-ordinate), and I do this for both:
> >      >      >
> >      >      >     DMGetCoordinatesLocal();
> >      >      >     DMGetCellCoordinatesLocal();
> >      >      >
> >      >      >     How do I know I am not cubing one coordinate multiple
> >     times?
> >      >      >
> >      >      >
> >      >      > Good question. Right now, the only connection between the
> >     two sets of coordinates is DMLocalizeCoordinates(). Since
> >      >     sometimes
> >      >      > people want to do non-trivial things to
> >      >      > coordinates, I prefer not to push in an API for "just"
> >     scaling, but I could be convinced
> >      >      > the other way.
> >      >      >
> >      >      >    Thanks,
> >      >      >
> >      >      >       Matt
> >      >      >
> >      >      >     Thanks, Berend.
> >      >      >
> >      >      >     On 5/17/23 16:10, Matthew Knepley wrote:
> >      >      >      > On Wed, May 17, 2023 at 10:02 AM Berend van Wachem
> >      >      >      > <berend.vanwachem at ovgu.de
> >     <mailto:berend.vanwachem at ovgu.de> <mailto:berend.vanwachem at ovgu.de
> >     <mailto:berend.vanwachem at ovgu.de>> <mailto:berend.vanwachem at ovgu.de
> >     <mailto:berend.vanwachem at ovgu.de>
> >      >     <mailto:berend.vanwachem at ovgu.de
> >     <mailto:berend.vanwachem at ovgu.de>>> <mailto:berend.vanwachem at ovgu.de
> >     <mailto:berend.vanwachem at ovgu.de> <mailto:berend.vanwachem at ovgu.de
> >     <mailto:berend.vanwachem at ovgu.de>>
> >      >      >     <mailto:berend.vanwachem at ovgu.de
> >     <mailto:berend.vanwachem at ovgu.de> <mailto:berend.vanwachem at ovgu.de
> >     <mailto:berend.vanwachem at ovgu.de>>>>> wrote:
> >      >      >      >
> >      >      >      >     Dear PETSc Team,
> >      >      >      >
> >      >      >      >     We are using DMPlex, and we create a mesh using
> >      >      >      >
> >      >      >      >     DMPlexCreateBoxMesh (.... );
> >      >      >      >
> >      >      >      >     and get a uniform mesh. The mesh is periodic.
> >      >      >      >
> >      >      >      >     We typically want to "scale" the coordinates
> >     (vertices) of the mesh,
> >      >      >      >     and
> >      >      >      >     to achieve this, we call
> >      >      >      >
> >      >      >      >     DMGetCoordinatesLocal(dm, &coordinates);
> >      >      >      >
> >      >      >      >     and scale the entries in the Vector coordinates
> >     appropriately.
> >      >      >      >
> >      >      >      >     and then
> >      >      >      >
> >      >      >      >     DMSetCoordinatesLocal(dm, coordinates);
> >      >      >      >
> >      >      >      >
> >      >      >      >     After this, we localise the coordinates by
> calling
> >      >      >      >
> >      >      >      >     DMLocalizeCoordinates(dm);
> >      >      >      >
> >      >      >      >     This worked fine up to PETSc 3.18, but with
> >     versions after this, the
> >      >      >      >     coordinates we get from the call
> >      >      >      >
> >      >      >      >     DMPlexGetCellCoordinates(dm, CellID, &isDG,
> >     &CoordSize,
> >      >      >      >     &ArrayCoordinates, &Coordinates);
> >      >      >      >
> >      >      >      >     are no longer correct if the mesh is periodic.
> >     A number of the
> >      >      >      >     coordinates returned from calling
> >     DMPlexGetCellCoordinates are wrong.
> >      >      >      >
> >      >      >      >     I think, this is because DMLocalizeCoordinates
> >     is now automatically
> >      >      >      >     called within the routine DMPlexCreateBoxMesh.
> >      >      >      >
> >      >      >      >     So, my question is: How should we scale the
> >     coordinates from a periodic
> >      >      >      >     DMPlex mesh so that they are reflected
> >     correctly when calling both
> >      >      >      >     DMGetCoordinatesLocal and
> >     DMPlexGetCellCoordinates, with PETSc versions
> >      >      >      >       >= 3.18?
> >      >      >      >
> >      >      >      >
> >      >      >      > I think we might have to add an API function. For
> >     now, when you scale
> >      >      >      > the coordinates,
> >      >      >      > can you scale both copies?
> >      >      >      >
> >      >      >      >    DMGetCoordinatesLocal()
> >      >      >      >    DMGetCellCoordinatesLocal();
> >      >      >      >
> >      >      >      > and then set them back.
> >      >      >      >
> >      >      >      >    Thanks,
> >      >      >      >
> >      >      >      >       Matt
> >      >      >      >
> >      >      >      >     Many thanks, Berend.
> >      >      >      >
> >      >      >      > --
> >      >      >      > 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://www.cse.buffalo.edu/~knepley/>
> >     <https://www.cse.buffalo.edu/~knepley/
> >     <https://www.cse.buffalo.edu/~knepley/>>
> >      >     <https://www.cse.buffalo.edu/~knepley/
> >     <https://www.cse.buffalo.edu/~knepley/>
> >     <https://www.cse.buffalo.edu/~knepley/
> >     <https://www.cse.buffalo.edu/~knepley/>>>
> >     <http://www.cse.buffalo.edu/~knepley/
> >     <http://www.cse.buffalo.edu/~knepley/>
> >      >     <http://www.cse.buffalo.edu/~knepley/
> >     <http://www.cse.buffalo.edu/~knepley/>>
> >      >      >     <http://www.cse.buffalo.edu/~knepley/
> >     <http://www.cse.buffalo.edu/~knepley/>
> >     <http://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/
> >     <https://www.cse.buffalo.edu/~knepley/>
> >     <https://www.cse.buffalo.edu/~knepley/
> >     <https://www.cse.buffalo.edu/~knepley/>>
> >     <http://www.cse.buffalo.edu/~knepley/
> >     <http://www.cse.buffalo.edu/~knepley/>
> >      >     <http://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/
> >     <https://www.cse.buffalo.edu/~knepley/>
> >     <http://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/20230522/45a76793/attachment-0001.html>


More information about the petsc-users mailing list