[petsc-users] Example for adaptive high-order FEM in PETSc

Matthew Knepley knepley at gmail.com
Sun Jan 13 14:01:59 CST 2019


On Sat, Jan 12, 2019 at 6:59 PM Manav Bhatia via petsc-users <
petsc-users at mcs.anl.gov> wrote:

> I have been studying the source code and the manual entries in the
> documentation and have a better understanding of setting up a high-order DG
> analysis.
>
> I have been able to identify the routines that set the order of the space
> and its continuity (CG/DG).
>
> The following is still unclear and I would appreciate some guidance:
>
> — A implicit DG solve of the inviscid Euler equations would require
> Jacobian contributions from the Riemann Solver at element interfaces. I do
> not see the ability to define Jacobian contribution in
> PetscDSSetRiemannSolver. What would be the best approach to do this?
>

You would set a BdResidual/BdJacobian function. You need my branch that
makes these integrals 2-sided. However, you would also probably need to
make an aux field with the geometric quantities you need. I can help you do
that. Then we could roll that up and put it in the library.


> — It appears that the only two choices for polynomials are Lagrange and
> “simple”. I am not clear what the “simple” basis implies. Are these just
> powers of x, y and z?
>

Correct. We would have to stick in another polynomial space. This is
actually pretty easy. We already have code for Legendre in the
PetscQuadrature code, and we should also probably add Chebyshev. In
addition, you probably want GLL quadrature. Some frustrating
people put GLL code in TS, which should be pulled out and made a
PetscQuadrature. You can see that this is some coding to do, but
none of it is design/hard, it just needs to be done.

  Thanks,

    Matt


> — Are Legendre polynomials available elsewhere in PETSc? I saw some
> mention of it in TS, but not in FE.
>
> Thanks,
> Manav
>
>
>
> > On Jan 12, 2019, at 2:46 AM, Manav Bhatia <bhatiamanav at gmail.com> wrote:
> >
> > While looking through some examples I see that both PetscDSSetJacobian
> and DMPlexSNESComputeJacobianFEM are being used to set the functions to
> evaluate the Jacobian contributions on each rank.
> >
> > This is my understanding after studying these (please correct if I am
> wrong):
> >
> > — The PetscDSSetJacobian (and its corresponding residual routine) does
> not seem to have the concept of an element. Instead, it provides the
> derivatives of solution at a point and requires that the coefficients for
> each form g_0, g_1, … be evaluated. Then, PETSc combines that in the
> background.
> >
> > — DMPlexSNESComputeJacobianFEM is a method that calls PetscDSSetJacobian
> to compute the Jacobian.
> >
> >
> > If this is correct, I am a bit confused about the following:
> >
> > — For a DG solve, I can see that PetscDSSetRiemannSolver can be used to
> compute the convective flux at the interface of two elements. However, I
> don’t see how the Jacobian contribution of this flux can be computed and
> added to the system Jacobian.
> >
> > — A more conventional FEM assembly iterates over elements, then over the
> quadrature points, where the shape functions and their derivatives are
> initialized and used for computation of residual and Jacobian. Is it
> possible to follow this procedure in PETSc with DMPlex/DMForest and still
> use multigrid? Is there an example that demonstrates this?
> >
> > — Is there some form of an index for elements that exist on the local
> processor that can be used to iterate on local elements?
> >
> >
> > I have tried to look, but could not find a document that outlines these
> concepts. Is there one that exists and I have missed it?
> >
> > I would greatly appreciate some direction.
> >
> > Regards,
> > Manav
> >
> >
> >
> >
> >> On Jan 11, 2019, at 10:03 PM, Manav Bhatia <bhatiamanav at gmail.com>
> wrote:
> >>
> >> Hi,
> >>
> >>  I am interested in setting up a high-order DG FEM solution with PETSc
> on an octree mesh with adaptivity and geometric multigrid preconditioning.
> Is there an example that demonstrates this with either DMForest or DMPlex.
> >>
> >> Thanks,
> >> Manav
> >>
> >>
> >
>
>

-- 
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/20190113/bda7cfc0/attachment-0001.html>


More information about the petsc-users mailing list