<div dir="ltr"><div dir="ltr">On Wed, Mar 18, 2020 at 12:58 PM Yann Jobic <<a href="mailto:yann.jobic@univ-amu.fr">yann.jobic@univ-amu.fr</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">Hi matt,<br>
<br>
Le 3/17/2020 à 4:00 PM, Matthew Knepley a écrit :<br>
> On Mon, Mar 16, 2020 at 5:20 PM Yann Jobic <<a href="mailto:yann.jobic@univ-amu.fr" target="_blank">yann.jobic@univ-amu.fr</a> <br>
> <mailto:<a href="mailto:yann.jobic@univ-amu.fr" target="_blank">yann.jobic@univ-amu.fr</a>>> wrote:<br>
> <br>
> Hi all,<br>
> <br>
> I would like to implement a nodal DG with the DMPlex interface.<br>
> Therefore, i must add the internal nodes to the DM (GLL nodes), with<br>
> the<br>
> constrains :<br>
> 1) Add them as solution points, with correct coordinates (and keep the<br>
> good rotational ordering)<br>
> 2) Find the shared nodes at faces in order to compute the fluxes<br>
> 3) For parallel use, so synchronize the ghost node at each time steps<br>
> <br>
> <br>
> Let me get the fundamentals straight before advising, since I have never <br>
> implemented nodal DG.<br>
> <br>
> 1) What is shared?<br>
I need to duplicate an edge in 2D, or a facet in 3D, and to sync it <br>
after a time step, in order to compute the numerical fluxes <br>
(Lax-Friedrichs at the beginning).<br></blockquote><div><br></div><div>I should have been more specific, but I think I see what you want. You do not "share" unknowns between cells,</div><div>so all unknowns should be associated with some cell in the Section.</div><div><br></div><div>You think of some cell unknowns as being "connected" to a face, so when you want to calculate a flux, you need</div><div>the unknowns from the adjacent cell in order to do it. In order to do this, I would partition with overlap=1, which</div><div>is what we do for finite volume, which has the same adjacency needs. You might also set</div><div><br></div><div> <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMSetAdjacency.html">https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/DMSetAdjacency.html</a></div><div> </div><div>to PETSC_TRUE, PETSC_FALSE, but you are probably doing everything matrix-free if you are using DG.</div><div>The above is optimal for FV, but not for DG because you communicate more than you absolutely have to.</div><div><br></div><div>A more complicated, but optimal, thing to do would be to assign interior dofs to the cell, and two sets of dofs</div><div>to each face, one for each cell. Then you only communicate the face dofs. Its just more bookkeeping for you,</div><div>but it will work in parallel just fine.</div><div><br></div><div>I don't think you need extra vertices, or coordinates, and for output I recommend using DMPlexProject() to get</div><div>the solution in some space that can be plotted like P1, or anything else supported by your visualization.</div><div><br></div><div> Thanks,</div><div><br></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">
> <br>
> We have an implementation of spectral element ordering <br>
> (<a href="https://gitlab.com/petsc/petsc/-/blob/master/src/dm/impls/plex/examples/tutorials/ex6.c" rel="noreferrer" target="_blank">https://gitlab.com/petsc/petsc/-/blob/master/src/dm/impls/plex/examples/tutorials/ex6.c</a>). <br>
> Those share<br>
> the whole element boundary.<br>
> <br>
> 2) What ghosts do you need?<br>
In order to compute the numerical fluxes of one element, i need the <br>
values of the surrounding nodes connected to the adjacent elements.<br>
> <br>
> 3) You want to store real space coordinates for a quadrature?<br>
It should be basically the same as PetscFE of higher order.<br>
I add some vertex needed to compute a polynomal solution of the desired <br>
order. That means that if i have a N, order of the local approximation, <br>
i need 0.5*(N+1)*(N+2) vertex to store in the DMPlex (in 2D), in order to :<br>
1) have the correct number of dof<br>
2) use ghost nodes to sync the values of the vertex/edge/facet for <br>
1D/2D/3D problem<br>
2) save correctly the solution<br>
<br>
Does it make sense to you ?<br>
<br>
Maybe like <br>
<a href="https://www.mcs.anl.gov/petsc/petsc-current/src/ts/examples/tutorials/ex11.c.html" rel="noreferrer" target="_blank">https://www.mcs.anl.gov/petsc/petsc-current/src/ts/examples/tutorials/ex11.c.html</a><br>
With the use of the function SplitFaces, which i didn't fully understood <br>
so far.<br>
<br>
Thanks,<br>
<br>
Yann<br>
<br>
> <br>
> We usually define a quadrature on the reference element once.<br>
> <br>
> Thanks,<br>
> <br>
> Matt<br>
> <br>
> I found elements of answers in those threads :<br>
> <a href="https://lists.mcs.anl.gov/pipermail/petsc-users/2016-August/029985.html" rel="noreferrer" target="_blank">https://lists.mcs.anl.gov/pipermail/petsc-users/2016-August/029985.html</a><br>
> <a href="https://lists.mcs.anl.gov/mailman/htdig/petsc-users/2019-October/039581.html" rel="noreferrer" target="_blank">https://lists.mcs.anl.gov/mailman/htdig/petsc-users/2019-October/039581.html</a><br>
> <br>
> However, it's not clear for me where to begin.<br>
> <br>
> Quoting Matt, i should :<br>
> " DMGetCoordinateDM(dm, &cdm);<br>
> <Set field information into cdm instead of dm><br>
> DMCreateLocalVector(cdm, &coordinatesLocal);<br>
> <Fill in higher order coordinate values><br>
> DMSetCoordinatesLocal(dm, coordinatesLocal);"<br>
> <br>
> However, i will not create ghost nodes this way. And i'm not sure to<br>
> keep the good ordering.<br>
> This part should be implemented in the PetscFE interface, for high<br>
> order<br>
> discrete solutions.<br>
> I did not succeed in finding the correct part of the source doing it.<br>
> <br>
> Could you please give me some hint to begin correctly thoses tasks ?<br>
> <br>
> Thanks,<br>
> <br>
> Yann<br>
> <br>
> <br>
> <br>
> -- <br>
> What most experimenters take for granted before they begin their <br>
> experiments is infinitely more interesting than any results to which <br>
> their experiments lead.<br>
> -- Norbert Wiener<br>
> <br>
> <a href="https://www.cse.buffalo.edu/~knepley/" rel="noreferrer" target="_blank">https://www.cse.buffalo.edu/~knepley/</a> <<a href="http://www.cse.buffalo.edu/~knepley/" rel="noreferrer" target="_blank">http://www.cse.buffalo.edu/~knepley/</a>><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>