[petsc-dev] Periodic meshes with <3 elements per edge?
Matthew Knepley
knepley at gmail.com
Wed Aug 14 13:03:47 CDT 2019
On Wed, Aug 14, 2019 at 11:46 AM Jed Brown <jed at jedbrown.org> wrote:
> Matthew Knepley <knepley at gmail.com> writes:
>
> > On Tue, Aug 13, 2019 at 7:35 PM Stefano Zampini <
> stefano.zampini at gmail.com>
> > wrote:
> >
> >>
> >>
> >> On Aug 14, 2019, at 1:19 AM, Jed Brown via petsc-dev <
> >> petsc-dev at mcs.anl.gov> wrote:
> >>
> >> [Cc: petsc-dev]
> >>
> >> Also, why is our current mode of localized coordinates preferred over
> >> the coordinate DM being non-periodic? Is the intent to preserve that
> >> for every point in a DM, the point is also valid in the coordinate DM?
> >>
> >> Yes.
> >
> >> Can there be "gaps" in a chart?
> >>
> >> Yes.
> >
> >> I've been digging around in the implementation because there is no
> >> documentation of localized coordinates, but it feels more complicated
> >> than I'd have hoped.
> >>
> >>
> >> A while ago, “localization” of coordinates was supporting only very
> simple
> >> cases, where periodic points were identified though the ‘maxCell’
> parameter
> >> (used to compute the proper cell coordinates). I think this is the
> reason
> >> why you need at least 3 cells to support periodicity, since the BoxMesh
> >> constructor uses the maxCell trick.
> >>
> >
> > This was my original conception since it was the only fully automatic way
> > to apply periodicity on an unstructured mesh that
> > was read from a file or generator.
>
> I take it this is in reference to this logic.
>
> PetscErrorCode DMLocalizeCoordinateReal_Internal(DM dm, PetscInt dim,
> const PetscReal anchor[], const PetscReal in[], PetscReal out[])
> {
> PetscInt d;
>
> PetscFunctionBegin;
> if (!dm->maxCell) {
> for (d = 0; d < dim; ++d) out[d] = in[d];
> } else {
> for (d = 0; d < dim; ++d) {
> if ((dm->bdtype[d] != DM_BOUNDARY_NONE) && (PetscAbsReal(anchor[d] -
> in[d]) > dm->maxCell[d])) {
> out[d] = anchor[d] > in[d] ? dm->L[d] + in[d] : in[d] - dm->L[d];
> } else {
> out[d] = in[d];
> }
> }
> }
> PetscFunctionReturn(0);
> }
>
>
> This implies that the mesh be aligned to coordinates (at least in the
> periodic dimension)
No, it does not. That is the point, that you can have edges crossing the
periodic boundary.
> and that dm->L[d] be set (perhaps by the user?).
>
Yes. It something is on the periodic torus, you should know the length of
the torus.
> And if the mesh generator specifies periodicity, as with Gmsh below,
> this logic isn't needed.
>
True.
> >> Now, you can also inform Plex about periodicity without the maxCell
> trick
> >> , see e.g.
> >>
> https://bitbucket.org/petsc/petsc/src/6a494beb09767ff86fff34131928e076224d7569/src/dm/impls/plex/plexgmsh.c#lines-1468
> .
> >> In this case, it is user responsibility to populate the cell part of the
> >> coordinate section with the proper localized coordinates.
> >>
> >
> > This is a great addition from Stefano and Lisandro, but note that it is
> > nontrivial. The user has to identify the
> > periodic boundary.
>
> An (unstructured) mesh generator has to support periodicity anyway to
> generate a mesh with that topology.
>
Or you can fake it.
> >> The DMPlex code fully support coordinates localized only in those cells
> >> touching the periodic boundary. (I’m not a fan of this, since it
> requires a
> >> lot of ‘if’ ‘else’ switches )
> >>
> >
> > I believe that Lisandro wanted this to cut down on redundant storage of
> > coordinates.
>
> What was the rationale for having these special cells with
> localized/DG-style coordinates versus storing coordinates as a
> non-periodic continuous field?
I do not understand what you mean.
> The local points could be distinct for
> both fields and coordinates, with the global SF de-duplicating the
> periodic points for fields, versus leaving them distinct for
> coordinates.
Oh, no I would never do that.
> (Global coordinates are often not used.) It seems like a
> simpler model to me, and would eliminate a lot of if/else statements,
> but you've been thinking about this more than me.
>
> >> I think domain_box_size 1 is not possible, we can probably allow
> >> domain_box_size 2.
> >>
> >
> > Technically, now a single box is possible with higher order coordinate
> > spaces, but you have to do everything by hand
> > and it completely untested.
>
> Why would higher order coordinate spaces be needed? I could localize
> coordinates for all cells.
>
I mean if you wanted a single cell that was periodic, you could do that
with higher order coordinates.
Matt
I'm asking about this partly for spectral discretizations and partly as
> a hack to use code written for 3D to solve 2D problems.
>
--
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-dev/attachments/20190814/a43a885b/attachment.html>
More information about the petsc-dev
mailing list