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