Finite Element Methods

The superscreen.fem module contains functions useful for performing computations on a triangular mesh. Many of the functions take as arguments two arrays, points and triangles, which together define the mesh. points, which has shape (n, 2), defines the coordinates of the triangle vertices. Each row in triangles, which has shape (m, 3) gives the indices of the three vertices in a single triangle, such that points[triangles] is a shape (m, 3, 2) array where each row gives the coordinates of the three vertices of a triangle in the mesh.

The weight matrix, mass matrix, and Laplacian operator are all sparse, meaning that most of their entries are zero. In order to save memory, functions that generate these arrays give the option to return a scipy.sparse matrix. All matrices are converted to dense numpy arrays when simulating a superscreen.Device, however.

superscreen.fem.triangle_areas(points, triangles)[source]

Calculates the area of each triangle.

Parameters:
  • points (ndarray) – Shape (n, 2) array of x, y coordinates of vertices

  • triangles (ndarray) – Shape (m, 3) array of triangle indices

Return type:

ndarray

Returns:

Shape (m, ) array of triangle areas

superscreen.fem.in_polygon(poly_points, query_points, radius=0)[source]

Returns a boolean array indicating which points in query_points lie inside the polygon defined by poly_points.

Parameters:
Return type:

Union[bool, ndarray]

Returns:

A shape (n, ) boolean array indicating which query_points lie inside the polygon. If only a single query point is given, then a single boolean value is returned.

superscreen.fem.centroids(points, triangles)[source]

Returns x, y coordinates for triangle centroids (centers of mass).

Parameters:
  • points (ndarray) – Shape (n, 2) array of x, y coordinates of vertices.

  • triangles (ndarray) – Shape (m, 3) array of triangle indices.

Return type:

ndarray

Returns:

Shape (m, 2) array of triangle centroid (center of mass) coordinates

superscreen.fem.adjacency_matrix(triangles, sparse=True)[source]

Computes the adjacency matrix for a given set of triangles.

Parameters:
  • triangles (ndarray) – Shape (m, 3) array of triangle indices

  • sparse (bool) – Whether to return a sparse array or numpy ndarray.

Return type:

Union[ndarray, csr_array]

Returns:

Shape (n, n) adjacency matrix, where n = triangles.max() + 1

superscreen.fem.adj_directed_tri_indices(triangles, num_sites)[source]

Construct the directed adjacency matrix.

Each element (i, j) represents an edge in the mesh, and the value at (i, j) is 1 + the index of a triangle containing that edge.

Parameters:
  • triangles (ndarray) – The triangle indices, shape (m, 3)

  • num_sites (int) – The number of sites in the mesh

Return type:

csc_array

Returns:

A directed adjacency matrix containing triangle indices + 1

superscreen.fem.weights_inv_euclidean(points, triangles, sparse=True)[source]

Weights edges by the inverse Euclidean distance of the edge lengths.

Parameters:
  • points (ndarray) – Shape (n, 2) array of x, y coordinates of vertices.

  • triangles (ndarray) – Shape (m, 3) array of triangle indices.

  • sparse (bool) – Whether to return a sparse matrix or numpy ndarray.

Return type:

Union[ndarray, lil_array]

Returns:

Shape (n, n) matrix of vertex weights

superscreen.fem.weights_half_cotangent(points, triangles, sparse=True)[source]

Weights edges by half of the cotangent of the two angles opposite the edge.

Parameters:
  • points (ndarray) – Shape (n, 2) array of x, y coordinates of vertices.

  • triangles (ndarray) – Shape (m, 3) array of triangle indices.

  • sparse (bool) – Whether to return a sparse matrix or numpy ndarray.

Return type:

Union[ndarray, lil_array]

Returns:

Shape (n, n) matrix of vertex weights

superscreen.fem.calculate_weights(points, triangles, method, sparse=True)[source]

Returns the weight matrix, calculated using the specified method.

Parameters:
  • points (ndarray) – Shape (n, 2) array of x, y coordinates of vertices.

  • triangles (ndarray) – Shape (m, 3) array of triangle indices.

  • method (str) – Method for calculating the weights. One of: “uniform”, “inv_euclidean”, or “half_cotangent”.

  • sparse (bool) – Whether to return a sparse matrix or numpy ndarray.

Return type:

Union[ndarray, csr_array]

Returns:

Shape (n, n) matrix of vertex weights

superscreen.fem.laplace_operator(points, triangles, masses=None, weight_method='half_cotangent')[source]

Laplacian operator for the mesh (sometimes called Laplace-Beltrami operator).

The Laplacian operator is defined as inv(M) @ L, where M is the mass matrix and L is the Laplacian matrix.

Parameters:
  • points (ndarray) – Shape (n, 2) array of x, y coordinates of vertices.

  • triangles (ndarray) – Shape (m, 3) array of triangle indices.

  • masses (Optional[ndarray]) – Pre-computed mass matrix, i.e., the vertex areas.

  • weight_method (Literal['uniform', 'half_cotangent', 'inv_euclidean']) – Method for calculating the weights. One of: “uniform”, “inv_euclidean”, or “half_cotangent”.

Return type:

csr_array

Returns:

Shape (n, n) Laplacian operator

superscreen.fem.gradient_triangles(points, triangles, areas=None)[source]

Returns the triangle gradient operators Gx and Gy.

Given a mesh with n vertices and m triangles and a scalar field f defined at the mesh vertices, Gx and Gy are shape (m, n) matrices such that Gx @ f and Gy @ f compute the gradients of f along x and y, evaluated at the triangle centroids.

Parameters:
  • points (ndarray) – Shape (n, 2) array of x, y coordinates of vertices

  • triangles (ndarray) – Shape (m, 3) array of triangle indices

  • _areas – Shape (m, ) array of triangle areas

  • areas (ndarray | None) –

Return type:

Tuple[csr_array, csr_array]

Returns:

x and y gradient matrices, both of which have shape (m, n)

superscreen.fem.gradient_vertices(points, triangles, areas=None)[source]

Returns the vertex gradient operators gx and gy.

Given a mesh with n vertices and m triangles and a scalar field f defined at the mesh vertices, gx and gy are shape (n, n) matrices such that gx @ f and gy @ f compute the gradients of f along x and y, evaluated at the vertices.

The vertex gradient operators are calculated by averaging the the triangle gradient operators over all triangles adjacent to each vertex.

Parameters:
  • points (ndarray) – Shape (n, 2) array of x, y coordinates of vertices

  • triangles (ndarray) – Shape (m, 3) array of triangle indices

  • areas (Optional[ndarray]) – Shape (m, ) array of triangle areas

Return type:

Tuple[csr_array, csr_array]

Returns:

x and y gradient matrices, both of which have shape (n, n)