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

Matthew Knepley knepley at gmail.com
Sun Jan 13 13:57:27 CST 2019


On Sat, Jan 12, 2019 at 2:47 AM Manav Bhatia via petsc-users <
petsc-users at mcs.anl.gov> 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.
>

Let me say first that we do not support DG. We have plans to do this, but
it does not currently work.


> 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.
>

Yes, this is true. Here is the idea: DS abstracts the problem itself,
without regard to discretization.The division is not exactly what you
expect currently, but I think it is the correct one. You are solving the
weak form of the problem, which is definitely different than the strong
form, so the DS assumes a weak form for the specification, but independent
of the bases involved.


> — DMPlexSNESComputeJacobianFEM is a method that calls PetscDSSetJacobian
> to compute the Jacobian.
>

Yes.


> 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.
>

That is really only designed for FV. You are right that we could repurpose
it to support DG. That was not how I was going to proceed. Instead I have a
branch that fixes boundary integration so that all interface integrals get
2-sided information. This would allow all the facets integrals that you
want with proper transforms, etc. However, I would also have to put in
support for geometric quantities needed for DG penalization. That is
simple, I know, but everything takes time.


> — 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?
>

Yes, that is how all the FEM examples work. For instance, you can look at
SNES ex12 or ex17 or ex62 or ex77. They iterates over batches
of cells, then projectbasis function to the quad points, evaluate user
point functions on the quad points, and then take the inner product with
test functions. After the element matrix is computed, it is set into the
global Jacobian.


> — Is there some form of an index for elements that exist on the local
> processor that can be used to iterate on local elements?
>

We just use local point number to iterate over cells.


> 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?
>

There ia a manual chapter on the low-level organization, but nothing on the
traversal because I am still experimenting.

  Thanks,

   Matt


> 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/a26c0106/attachment.html>


More information about the petsc-users mailing list