[petsc-dev] Periodic meshes with <3 elements per edge?
Patrick Sanan
patrick.sanan at gmail.com
Wed Aug 14 14:17:59 CDT 2019
Am Mi., 14. Aug. 2019 um 20:15 Uhr schrieb Jed Brown via petsc-dev <
petsc-dev at mcs.anl.gov>:
> Matthew Knepley <knepley at gmail.com> writes:
>
> > 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.
>
> I mean that u(x) = u(x+(L,0)) rather than u(x) = u(x+v) where ||v||=L is
> arbitrary.
>
> >> An (unstructured) mesh generator has to support periodicity anyway to
> >> generate a mesh with that topology.
> >
> > Or you can fake it.
>
> How would you handle a different number/distribution of points on
> surfaces that are meant to be joined periodically?
>
> >> >> 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.
>
> Can you help me understand why that model is bad?
>
I'm also interested in the answer to this question, because I am
considering something similar for DMStag; if DM has a periodic BC, the
corresponding coordinate DM has a "none" BC, so the boundary points are
duplicated - this would hopefully make it much easier to locate particles
in elements.
>
> >> (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.
>
> But it would also work with linear coordinates, no? (Assume the user
> localizes their own coordinates.)
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20190814/92526d84/attachment.html>
More information about the petsc-dev
mailing list