New nodal Hdiv finite element spaces

62 views
Skip to first unread message

Eldar Khattatov

unread,
Jun 12, 2017, 8:28:12 PM6/12/17
to deal.II developers
Dear all,

As a part of my phd thesis I've been working on mixed finite elements for Darcy problem that allow for localization of velocities leading to cell-centered pressure system. In the lowest order case, it is what's called a Multipoint Flux Mixed Finite Element method (it was published by Wheeler and Yotov in 2006, 2010), which works in 2d with BDM_1 spaces and in 3d with enhanced BDM_1 spaces and trapezoid quadrature rule. This idea extends to higher order cases as well, although it requires new spaces that are basically V_k = RT_(k-1) + B_k, where B_k is the part, that does not ruin the Hdiv conformity but adds extra DoFs so that the total dimension of V_k = d*dim(Q_k), hence it is possible to do local elimination of velocities around nodes provided the use of appropriate order Gauss-Lobatto quadrature rule. Unfortunately, the paper on these is not published yet (although I can send you a preprint which has all the details and proofs, as well as the numerical results), but all the theoretical results are proven already, and show that the elements V_k form stable pair with DGQ_(k-1) and have correct approximation properties. Numerical tests with mixed Darcy problem also show the correct convergence rates for arbitrary degree elements. 

As these spaces have much more degrees of freedom than RT, or BDM providing only the same order of accuracy as the latter, their only use case is in MFMFE type methods which can be considered as CCFD in the lowest order case or mixed finite volumes in higher order case. On the other hand, due to localization properties and the possibility to invert mass matrix locally, the overall time for obtaining a solution is significantly smaller compared to classical Schur complement approach for mixed problem, despite the spaces having many more DoFs. These spaces are also "nodal" as the Raviart-Thomas Nodal space in Deal.II and guarantee that the mass matrix (even with full tensor coefficient) is block diagonal, where each block is associated to a vertex (lowest order case) of "support point" in higher order case.

I've implemented these elements in Deal.II and would be happy to contribute them to the library if you find it valuable. The implementation consists of the new polynomials_rt_bubbles.cc and fe_rt_bubbles.cc (and corresponding .h, .inst.in) files. I also implemented tests similar to fe_rt_1.cc - fe_rt_16.cc, and polynomials_rt_1.cc. I could do more tests if needed. I can share my implementation of Darcy problem using the approach I described above with you if you are interested. I can also add it to the code gallery, but as soon as the paper gets submitted (should be pretty soon).

Thank you,
Eldar
2d_RTBubble_3.png
3d_RTBubble_2.png

Wolfgang Bangerth

unread,
Jun 14, 2017, 9:27:50 PM6/14/17
to dealii-d...@googlegroups.com, Eldar Khattatov

Eldar,

> As a part of my phd thesis I've been working on mixed finite elements for
> Darcy problem that allow for localization of velocities leading to
> cell-centered pressure system. In the lowest order case, it is what's called a
> Multipoint Flux Mixed Finite Element method (it was published by Wheeler and
> Yotov in 2006, 2010), which works in 2d with BDM_1 spaces and in 3d with
> enhanced BDM_1 spaces and trapezoid quadrature rule. This idea extends to
> higher order cases as well, although it requires new spaces that are basically
> V_k = RT_(k-1) + B_k, where B_k is the part, that does not ruin the Hdiv
> conformity but adds extra DoFs so that the total dimension of V_k =
> d*dim(Q_k), hence it is possible to do local elimination of velocities around
> nodes provided the use of appropriate order Gauss-Lobatto quadrature rule.
> Unfortunately, the paper on these is not published yet (although I can send
> you a preprint which has all the details and proofs, as well as the numerical
> results), but all the theoretical results are proven already, and show that
> the elements V_k form stable pair with DGQ_(k-1) and have correct
> approximation properties. Numerical tests with mixed Darcy problem also show
> the correct convergence rates for arbitrary degree elements.

I think I have seen Mary Wheeler talk about this combination of elements. It
looked really exciting to me!


> As these spaces have much more degrees of freedom than RT, or BDM providing
> only the same order of accuracy as the latter, their only use case is in MFMFE
> type methods which can be considered as CCFD in the lowest order case or mixed
> finite volumes in higher order case. On the other hand, due to localization
> properties and the possibility to invert mass matrix locally, the overall time
> for obtaining a solution is significantly smaller compared to classical Schur
> complement approach for mixed problem, despite the spaces having many more
> DoFs.

Yes, I think that's the way to see it. The number of DoFs doesn't matter. What
matters is the number of DoFs that you can't locally eliminate and that lead
to globally coupled linear systems.


> I've implemented these elements in Deal.II and would be happy to contribute
> them to the library if you find it valuable. The implementation consists of
> the new polynomials_rt_bubbles.cc and fe_rt_bubbles.cc (and corresponding .h,
> .inst.in) files. I also implemented tests similar to fe_rt_1.cc - fe_rt_16.cc,
> and polynomials_rt_1.cc.

Yes, we would be very interested in seeing these elements in deal.II.


> I could do more tests if needed. I can share my
> implementation of Darcy problem using the approach I described above with you
> if you are interested. I can also add it to the code gallery, but as soon as
> the paper gets submitted (should be pretty soon).

I think that, too, would be quite exciting. I'm going to add that we've always
had reviewers comment positively if we made the code for a paper available. I
imagine that would be the case also for your paper, assuming it's still at a
stage where you could make the promise to make your code available.

Best
Wolfgang


--
------------------------------------------------------------------------
Wolfgang Bangerth email: bang...@colostate.edu
www: http://www.math.colostate.edu/~bangerth/

Eldar Khattatov

unread,
Jun 15, 2017, 10:09:24 AM6/15/17
to deal.II developers, eld.kh...@gmail.com, bang...@colostate.edu
Dear Wolfgang,

Thank you for the reply. You can find the preprint here and my current implementation of the MFMFE method with deal.ii and these elements in my github repository here (darcy_mfmfe.cc and .h). 
I also pushed the elements to my branch of deal.ii which you can find here.

I can create a pull request if it is appropriate and/or more convenient to see changes that way.

Best,
Eldar
Reply all
Reply to author
Forward
0 new messages