fe_utils package

Subpackages

Submodules

fe_utils.finite_elements module

class fe_utils.finite_elements.FiniteElement(cell, degree, nodes, entity_nodes=None)[source]

Bases: object

A finite element defined over cell.

Parameters:
  • cell – the ReferenceCell over which the element is defined.

  • degree – the polynomial degree of the element. We assume the element spans the complete polynomial space.

  • nodes – a list of coordinate tuples corresponding to point evaluation node locations on the element.

  • entity_nodes – a dictionary of dictionaries such that entity_nodes[d][i] is the list of nodes associated with entity (d, i) of dimension d and index i.

Most of the implementation of this class is left as exercises.

cell

The ReferenceCell over which the element is defined.

degree

The polynomial degree of the element. We assume the element spans the complete polynomial space.

entity_nodes

A dictionary of dictionaries such that entity_nodes[d][i] is the list of nodes associated with entity (d, i).

interpolate(fn)[source]

Interpolate fn onto this finite element by evaluating it at each of the nodes.

Parameters:

fn – A function fn(X) which takes a coordinate vector and returns a scalar value.

Returns:

A vector containing the value of fn at each node of this element.

The implementation of this method is left as an exercise.

node_count

The number of nodes in this element.

nodes

The list of coordinate tuples corresponding to the nodes of the element.

nodes_per_entity

nodes_per_entity[d] is the number of entities associated with an entity of dimension d.

tabulate(points, grad=False)[source]

Evaluate the basis functions of this finite element at the points provided.

Parameters:
  • points – a list of coordinate tuples at which to tabulate the basis.

  • grad – whether to return the tabulation of the basis or the tabulation of the gradient of the basis.

Result:

an array containing the value of each basis function at each point. If grad is True, the gradient vector of each basis vector at each point is returned as a rank 3 array. The shape of the array is (points, nodes) if grad is False and (points, nodes, dim) if grad is True.

The implementation of this method is left as an exercise.

class fe_utils.finite_elements.LagrangeElement(cell, degree)[source]

Bases: FiniteElement

An equispaced Lagrange finite element.

Parameters:
  • cell – the ReferenceCell over which the element is defined.

  • degree – the polynomial degree of the element. We assume the element spans the complete polynomial space.

The implementation of this class is left as an exercise.

fe_utils.finite_elements.lagrange_points(cell, degree)[source]

Construct the locations of the equispaced Lagrange nodes on cell.

Parameters:
  • cell – the ReferenceCell

  • degree – the degree of polynomials for which to construct nodes.

Returns:

a rank 2 array whose rows are the coordinates of the nodes.

The implementation of this function is left as an exercise.

fe_utils.finite_elements.vandermonde_matrix(cell, degree, points, grad=False)[source]

Construct the generalised Vandermonde matrix for polynomials of the specified degree on the cell provided.

Parameters:
  • cell – the ReferenceCell

  • degree – the degree of polynomials for which to construct the matrix.

  • points – a list of coordinate tuples corresponding to the points.

  • grad – whether to evaluate the Vandermonde matrix or its gradient.

Returns:

the generalised Vandermonde matrix

The implementation of this function is left as an exercise.

fe_utils.function_spaces module

class fe_utils.function_spaces.Function(function_space, name=None)[source]

Bases: object

A function in a finite element space. The main role of this object is to store the basis function coefficients associated with the nodes of the underlying function space.

Parameters:
  • function_space – The FunctionSpace in which this Function lives.

  • name – An optional label for this Function which will be used in output and is useful for debugging.

function_space

The FunctionSpace in which this Function lives.

integrate()[source]

Integrate this Function over the domain.

Result:

The integral (a scalar).

interpolate(fn)[source]

Interpolate a given Python function onto this finite element Function.

Parameters:

fn – A function fn(X) which takes a coordinate vector and returns a scalar value.

name

The (optional) name of this Function

plot(subdivisions=None)[source]

Plot the value of this Function. This is quite a low performance plotting routine so it will perform poorly on larger meshes, but it has the advantage of supporting higher order function spaces than many widely available libraries.

Parameters:

subdivisions – The number of points in each direction to use in representing each element. The default is \(2d+1\) where \(d\) is the degree of the FunctionSpace. Higher values produce prettier plots which render more slowly!

values

The basis function coefficient values for this Function

class fe_utils.function_spaces.FunctionSpace(mesh, element)[source]

Bases: object

A finite element space.

Parameters:
  • mesh – The Mesh on which this space is built.

  • element – The FiniteElement of this space.

Most of the implementation of this class is left as an exercise.

cell_nodes

The global cell node list. This is a two-dimensional array in which each row lists the global nodes incident to the corresponding cell. The implementation of this member is left as an exercise

element

The FiniteElement of this space.

mesh

The Mesh on which this space is built.

node_count

The total number of nodes in the function space.

fe_utils.mesh module

class fe_utils.mesh.Mesh(vertex_coords, cell_vertices)[source]

Bases: object

A one or two dimensional mesh composed of intervals or triangles respectively.

Parameters:
  • vertex_coords – a vertex_count x dim array of the coordinates of the vertices in the mesh.

  • cell_vertices – a cell_count x (dim+1) array of the indices of the vertices of which each cell is made up.

adjacency(dim1, dim2)[source]

Return the set of dim2 entities adjacent to each dim1 entity. For example if dim1==2 and dim2==1 then return the list of edges (1D entities) adjacent to each triangle (2D entity).

The return value is a rank 2 numpy.array such that adjacency(dim1, dim2)[e1, :] is the list of dim2 entities adjacent to entity (dim1, e1).

This operation is only defined where self.dim >= dim1 > dim2.

This method is simply a more systematic way of accessing edge_vertices, cell_edges and cell_vertices.

cell

The ReferenceCell of which this Mesh is composed.

cell_edges

The indices of the edges incident to each cell (only for 2D meshes).

cell_vertices

The indices of the vertices incident to cell.

dim

The geometric and topological dimension of the mesh. Immersed manifolds are not supported.

edge_vertices

The indices of the vertices incident to edge (only for 2D meshes).

entity_counts

The number of entities of each dimension in the mesh. So entity_counts[0] is the number of vertices in the mesh.

jacobian(c)[source]

Return the Jacobian matrix for the specified cell.

Parameters:

c – The number of the cell for which to return the Jacobian.

Result:

The Jacobian for cell c.

vertex_coords

The coordinates of all the vertices in the mesh.

class fe_utils.mesh.UnitIntervalMesh(nx)[source]

Bases: Mesh

A mesh of the unit interval.

Parameters:

nx – The number of cells.

class fe_utils.mesh.UnitSquareMesh(nx, ny)[source]

Bases: Mesh

A triangulated Mesh of the unit square.

Parameters:
  • nx – The number of cells in the x direction.

  • ny – The number of cells in the y direction.

fe_utils.quadrature module

class fe_utils.quadrature.QuadratureRule(cell, degree, points, weights)[source]

Bases: object

A quadrature rule implementing integration of the reference cell provided.

Parameters:
Points:

a list of the position vectors of the quadrature points.

Weights:

the corresponding vector of quadrature weights.

cell

The ReferenceCell over which this quadrature rule is defined.

degree

The degree of precision of the quadrature rule.

integrate(function)[source]

Integrate the function provided using this quadrature rule.

Parameters:

function – A Python function taking a position vector as its single argument and returning a scalar value.

The implementation of this method is left as an exercise.

points

Two dimensional array, the rows of which form the position vectors of the quadrature points.

weights

The corresponding array of quadrature weights.

fe_utils.quadrature.gauss_quadrature(cell, degree)[source]

Return a Gauss-Legendre QuadratureRule.

Parameters:

fe_utils.reference_elements module

class fe_utils.reference_elements.ReferenceCell(vertices, topology, name)[source]

Bases: object

An object storing the geometry and topology of the reference cell.

Parameters:
  • vertices – a list of coordinate vectors corresponding to the coordinates of the vertices of the cell.

  • topology – a dictionary of dictionaries such that topology[d][i] is the list of vertices incident to the i-th entity of dimension d.

dim

The geometric and topological dimension of the reference cell.

entity_counts

The number of entities of each dimension.

point_in_entity(x, e)[source]

Return true if the point x lies on the entity e.

Parameters:
  • x – The coordinate vector of the point.

  • e – The (d, i) pair describing the entity of dimension d and index i.

topology

The vertices making up each topological entity of the reference cell.

vertices

The list of coordinate veritices of the reference cell.

fe_utils.reference_elements.ReferenceInterval = ReferenceInterval

A ReferenceCell storing the geometry and topology of the interval [0, 1].

fe_utils.reference_elements.ReferenceTriangle = ReferenceTriangle

A ReferenceCell storing the geometry and topology of the triangle with vertices [[0., 0.], [1., 0.], [0., 1.]].

fe_utils.utils module

fe_utils.utils.errornorm(f1, f2)[source]

Calculate the L^2 norm of the difference between f1 and f2.

Module contents