New API for various higher-order element kinds

Issue #92 resolved
Miklós Homolya created an issue

The immediate concern for us is spectral elements, although this proposal is also relevant for other element types, e.g., Bernstein, etc.

So FEniCS/Firedrake currently has equispaced Lagrange elements and their discontinuous counterparts. Gauss-Legendre (GL) and Gauss-Lobatto-Legendre (GLL) elements have recently been added to FIAT and UFL. They are also Lagrange elements (their basis functions are Lagrange polynomials), but with different node locations, providing a better condition number at higher order, and also functioning as quadrature rules (with appropriate weights).

So far so good, however, they are only for intervals, and 1-dimensional problems are not very exciting. In theory one could easily construct elements for hypercubes as tensor product elements, e.g., for quadrilaterals:

  • H1: GLL x GLL
  • H(curl): HCurl(GLL x GL) + HCurl(GL x GLL)
  • H(div): HDiv(GLL x GL) + HDiv(GL x GLL)
  • L2: GL x GL

Sounds straightforward, right? However, we support both extruded interval meshes (with interval x interval cells) and unstructured quadrilateral meshes. While in the extruded case such manual construction makes sense (though inconvenient), in the unstructured case it opens a whole can of worms. One may want to use, for example, P1 x P2 elements in the extruded case, i.e. P1 elements in the horizontal direction and P2 elements in the vertical. Such product elements do not make sense with unstructured quadrilaterals because there are no distinguished horizontal and vertical directions. This is why we had earlier decided not to accept product elements on quadrilaterals, only the names Q, RTCE, RTCF and DQ. Which was fine while equispaced Lagrange was the only available interval element to construct these from. Considering triangular prism cells and hexahedra, and also other possible (not yet supported) element types such as Chebyshev, Bernstein, Frechet, collapsed etc. we need an API that avoids combinatorial explosion, while sufficiently complete and safe.

Therefore we propose an extension of the API for constructing finite elements with a new optional argument node_set_hint (name subject to discussion) which:

  • defaults to None.
  • is of type string, e.g., "equispaced", "spectral".
  • can be a tuple of strings for tensor product elements (one string per subelement).

This way one could constructed a spectral RTCF elements, e.g., as follows:

FiniteElement('RTCF', quadrilateral, 6, node_set_hint='spectral')

which is in principle equivalent to HDiv(GLL6 x GL5) + HDiv(GL5 x GLL6). If one wants an equispaced higher-order elements, they could say, e.g.:

FiniteElement('Q', quadrilateral, 4, node_set_hint='equispaced')

Comments (10)

  1. Martin Sandve Alnæs

    I don't have time to think about this in any depth, @logg @meg @blechta @garth-wells opinions?

  2. Marie Elisabeth Rognes

    I understand the need. The name node_set_hint could (and should) be improved. Essentially, this is a way of specifying variants of degrees of freedom right?

  3. Joshua Chen

    Can I ask, how does one construct spectral elements with FEniCS/Firedrake? Is this possible currently?

  4. Miklós Homolya reporter

    @joshuawchen: node_set_hint has been renamed to variant, otherwise do as specified in the description.

    FiniteElement('Q', quadrilateral, 7, variant='spectral')
    FiniteElement('DQ', quadrilateral, 7, variant='spectral')
    FiniteElement('RTCF', quadrilateral, 7, variant='spectral')
    FiniteElement('RTCE', quadrilateral, 7, variant='spectral')
    

    You can definitely specify these in UFL. I am not sure if they work in FEniCS, but they definitely work in Firedrake.

  5. Miklós Homolya reporter

    Do you mean finite elements whose basis functions are frequency modes? No, they are not supported yet, and I don't know of anyone who is working on them. However, it should be fairly easy to implement them, one just needs to:

    • Implement 1D Fourier elements in FIAT.
    • Give it a variant= name, e.g., 'fourier', and wire in the elements in the FIAT/FInAT interface (which translates UFL elements to FIAT/FInAT elements).
  6. Joshua Chen

    My interest is using FEniCS for spectral and pseudo spectral Fourier Methods. Would this, assuming a uniform mesh, be amenable to this?

  7. Miklós Homolya reporter

    I would guess so, but I don't know enough of pseudo spectral Fourier methods to be able to answer this question with certainty.

  8. Log in to comment