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

import numpy as np
import basix.ufl
import ufl

element = ufl.FiniteElement("Lagrange", ufl.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 abitrary 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.

In FENiCSx, we use Basix to tabulate finite elements. Basix can convert any ufl-element into a basix finite element, which we in turn can evaluate.

basix_element = basix.ufl.convert_ufl_element(element)
/tmp/ipykernel_1314/1840541581.py:1: DeprecationWarning: Converting elements created in UFL to Basix elements is deprecated. You should create the elements directly using basix.ufl.element instead
  basix_element = basix.ufl.convert_ufl_element(element)

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(basix_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(basix_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 using a ufl.VectorElement, which you can give the number of components as an input argument.

vector_element = ufl.VectorElement("Lagrange", ufl.triangle, 2, dim=2)

We can also use basix for this

el = basix.ufl.element("Lagrange", "triangle", 2)
v_el = basix.ufl.blocked_element(el, shape=(2,))

Both these definitions of elements are compatible with DOLFINx. The reason for having both of them is that basix offers more tweaking of the finite elements that the UFL definition does. For more information regarding this see: Variants of Lagrange elements

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

m_el = ufl.MixedElement([vector_element, element])
m_el_basix = basix.ufl.mixed_element([v_el, basix_element])

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