Choosing and using a finite element

Contents

Choosing and using a finite element#

As seen in DOLFINx in general, we use the unified form language (UFL) to define variational forms. The power of this domain specific language is that it resembles mathematical syntax.

We will start with a standard problem, namely a projection:

(3)#\[\begin{align} u &= \frac{f(x,y,z)}{g(x,y,z)} \qquad \text{in } \Omega\subset \mathbb{R}^3. \end{align}\]

where \(\Omega\) is our computational domain, \(f\) and \(g\) are two known functions

Finite elements#

To solve this problem, we have to choose an appropriate finite element space to represent the function \(k\) and \(u\). There is a large variety of finite elements, for instance the Lagrange elements. The basis function for a first order Lagrange element is shown below


First order Lagrange basis functions

To symbolically represent this element, we use basix.ufl.element, that creates a UFL-compatible element using Basix.

import numpy as np
import basix.ufl

element = basix.ufl.element("Lagrange", "triangle", 1)

We note that we send in the finite element family, what cells we will use to represent the domain, and the degree of the function space. As this element is just a symbolic representation, it does not contain the information required to tabulate (evaluate) basis functions at an arbitrary point. UFL is written in this way, so that it can be used as a symbolic representation for a large range of finite element software, not just FENiCS.

Lets next tabulate the basis functions of our element at a given point p=(0, 0.5) in the reference element. the basix.finite_element.FiniteElement.tabulate function takes in two arguments, how many derivatives we want to compute, and at what points in the reference element we want to compute them.

points = np.array([[0., 0.5], [1, 0]], dtype=np.float64)
print(element.tabulate(0, points))
[[[ 5.00000000e-01 -2.77555756e-17  5.00000000e-01]
  [ 1.11022302e-16  1.00000000e+00 -2.15422807e-17]]]

We can also compute the derivatives at any points in the reference element

print(element.tabulate(1, points))
[[[ 5.00000000e-01 -2.77555756e-17  5.00000000e-01]
  [ 1.11022302e-16  1.00000000e+00 -2.15422807e-17]]

 [[-1.00000000e+00  1.00000000e+00 -4.30845614e-17]
  [-1.00000000e+00  1.00000000e+00 -4.30845614e-17]]

 [[-1.00000000e+00  0.00000000e+00  1.00000000e+00]
  [-1.00000000e+00  0.00000000e+00  1.00000000e+00]]]

Observe that the output we get from this command also includes the 0th order derivatives. Thus we note that the output has the shape (num_spatial_derivatives+1, num_points, num_basis_functions)

Not every function we want to represent is scalar valued. For instance, in fluid flow problems, the Taylor-Hood finite element pair is often used to represent the fluid velocity and pressure, where each component of the fluid velocity is in a Lagrange space. We represent this by adding a shape argument to the basix.ufl.element constructor.

vector_element = basix.ufl.element("Lagrange", "triangle", 2, shape=(2, ))

Basix allows for a large variety of extra options to tweak your finite elements, see for instance Variants of Lagrange elements for how to choose the node spacing in a Lagrange element.

To create the Taylor-Hood finite element pair, we use the basix.ufl.mixed_element

m_el = basix.ufl.mixed_element([vector_element, element])

There is a wide range of finite elements that are supported by ufl and basix. See for instance: Supported elements in basix/ufl.