Function spaces

Similarly to finite element codes function spaces form a central part of BEM++. To initialize a function space all we need is a grid object.

import bempp.api
grid = bempp.api.shapes.regular_sphere(3)

Let us now create a space of piecewise constant functions on the elements. This is the standard low-order space to represent Neumann data (that is normal derivatives) in boundary element methods.

space = bempp.api.function_space(grid, "DP", 0)

The first parameter of the function_space function is always the grid. The second gives the type of space, in this case “DP” for Discontinuous-Polynomial, and the third parameter is the order of the space (0 for piecewise constant).

We can now query the degrees of freedom of the space.

number_of_global_dofs = space.global_dof_count
print("Number of global dofs: {0}".format(number_of_global_dofs))
Number of global dofs: 512

For this space it is identical to the number of elements on the mesh. This is not necessarily always the case.

Let us now create a space of continuous, piecewise linear functions. This is the standard low-order space to represent Dirichlet data.

space = bempp.api.function_space(grid, "P", 1)
number_of_global_dofs = space.global_dof_count
print("Number of global dofs: {0}".format(number_of_global_dofs))
Number of global dofs: 258

The number of global dofs is now identical to the number of vertices in the grid.

Types of spaces

BEM++ supports the following types of spaces. The identifier for the function_space function is given in brackets.

  • Discontinuous polynomial spaces (DP). These are spaces of functions that are polynomial on each element but discontinuous across elements. The maximum order is 10.
  • Polynomial spaces (P). These are spaces of functions that are polynomial on each element and continuous across elements. The minimum order is zero. The maximum order is 10.
  • Polynomial spaces on barycentric grids (B-P). These are the same spaces as with the “P” identifier and the same number of degrees of freedom. But the underlying grid is a barycentric refinement of the original grid passed to the function (the barycentric refinement is created internally). This is needed in operator preconditioning. Currently, only order == 1 is supported.
  • Discontinuous polynomial spaces on barycentric grids (B-DP). As above, but corresponding to discontinuous polynomial spaces on the original grid. Currently, only order == 1 is supported.
  • Dual spaces of constant functions (DUAL). This is a space of constant functions defined on a dual grid (the dual grid is created internally from the grid object). These spaces form a stable dual pairing together with continuous, piecewise linear functions and are needed for certain opposite order preconditioners.
  • Raviart-Thomas Vector Space (“RT”). These are spaces of Raviart-Thomas basis functions. They are needed for integral operators in electromagnetic scattering. Currently, only low-order (order == 0) Raviart-Thomas spaces are supported.

For most scalar applications piecewise constant and continuous, piecewise linear spaces are sufficient. For the electric field and magnetic field operators spaces of Raviart-Thomas functions are required. The barycentric and dual spaces are for the assembly of certain types of preconditioned operators.

Local and global dofs

An important concept for spaces are global and local degrees of freedom. Global degrees of freedom are the actual dofs of the space associated with the basis functions, while local dofs are the coefficients of the restrictions of the basis functions to individual elements. Consider for example a space of continuous, piecewise linear basis functions. Each global dof is associated with a vertex. The corresponding basis function lives on all elements that are adjacent to this vertex. Sometimes it is useful to find out to what global dofs the basis functions on an element contribute. This is shown in the following example.

space = bempp.api.function_space(grid, "P", 1)
elements = list(grid.leaf_view.entity_iterator(0))
element = elements[0]
global_dofs, weights = space.get_global_dofs(element, dof_weights=True)
print("Map from local to global dofs on element: {0}".format(global_dofs))
print("Dof weights: {0} ".format(weights))
Local2GlobalDofMap on element: [2, 66, 68]
Dof weights: [1.0, 1.0, 1.0]

The map from local to global dofs on the element denotes which local basis function is assocated with which global basis functions. In this example, the element has three local basis functions (associated with the three vertices of the element), and the values of global_dofs are the indices of the global basis functions that they map to. The array weights returns scalar multipliers that are the weights with which the local basis function contributes to the global basis function. Hence, to obtain the local coefficient of a local basis function the corresponding global dof coefficient needs to be multiplied with the weight. For most basis types the weight is always 1. But for example, for Raviart-Thomas basis functions it can differ. Weights are only returned if the parameter dof_weights=True is set.

Function and class reference

Function spaces are instantiated by the function bempp.api.function_space. They return bempp.api.space.Space objects, which contain the necessary information that define a space over a grid.

bempp.api.function_space(grid, kind, order, domains=None, closed=True, strictly_on_segment=False, reference_point_on_segment=True, element_on_segment=False)

Return a space defined over a given grid.

Parameters:
  • grid (bempp.api.grid.Grid) – The grid object over which the space is defined.
  • kind (string) –

    The type of space. Currently, the following types are supported:

    “P” : Continuous and piecewise polynomial functions. “DP” : Discontinuous and elementwise polynomial functions. “RT”: Raviart-Thomas Vector spaces. “RWG”: RWG Vector spaces. “NC”: Nedelec Vector spaces. “SNC”: Scaled Nedelec Vector spaces. The Nedelec basis functions
    are scaled with the edge length so that they are identical to RWG functions crossed with the element normals.

    “B-P”: Polynomial spaces on barycentric grids. “B-DP”: Polynomial discontinuous spaces on barycentric grids. “B-RT”: Raviart-Thomas Vector spaces on barycentric grids. “B-RWG”: RWG Vector spaces on barycentric grids. “B-NC”: Nedelec Vector spaces on barycentric grids. “B-SNC”: Scaled Nedelec Vector spaces on barycentric grids.

    “DUAL”: Dual space on dual grid (only implemented for constants). “BC”: Buffa-Christian Vector space. “RBC”: Rotated Buffa-Christian Vector space of curl-conforming

    functions.
  • order (int) – The order of the space, e.g. 0 for piecewise const, 1 for piecewise linear functions.
  • domains (list) – List of integers specifying a list of physical entities of subdomains that should be included in the space.
  • closed (bool) – Specifies whether the space is defined on a closed or open subspace.
  • strictly_on_segment (bool) – Specifies whether local basis functions are truncated to the domains specified (True) or if they are allowed to extend past the domains (False). Default is False. This argument is only used for scalar continuous spaces.
  • reference_point_on_segment (bool) – If true only include a dof if its reference point (i.e. the dof position) is part of the segment. This argument is only used for discontinuous spaces (default is True).
  • element_on_segment (bool) – If true restrict the dofs to those whose support element is part of the segment (default is False).

Notes

The most frequent used types are the space of piecewise constant functions (kind=”DP”, order=0) and the space of continuous, piecewise linear functions (kind=”P”, order=1).

Either one of reference_point_on_segment or element_on_segment must be true for discontinuous spaces. For piecewise constant spaces neither of these two options has any effect.

This is a factory function that initializes a space object. To see a detailed help for space objects see the documentation of the instantiated object.

Examples

To initialize a space of piecewise constant functions use

>>> space = function_space(grid,"DP",0)

To initialize a space of continuous, piecewise linear functions, use

>>> space = function_space(grid,"P",1)
class bempp.api.space.Space(impl, comp_key)

Space of functions defined on a grid.

grid

bempp.api.grid.Grid – Grid over which to discretize the space.

dtype

numpy.dtype – Type of the basis functions in this space.

codomain_dimension

int – Number of components of values of functions in this space.

domain_dimension

int – Dimension of the domain on which the space is defined.

global_dof_count

int – Number of global degrees of freedom.

flat_local_dof_count

int – Total number of local degrees of freedom.

global_dof_interpolation_points

np.ndarray – (3xN) matrix of global interpolation points for the space, where each column is the coordinate of an interpolation point.

global_dof_interpolation_points

np.ndarray – (3xN) matrix of normal directions associated with the interpolation points.

codomain_dimension

Number of components of values of functions in this space.

For a scalar space this is one. For a vector space it is three.

discontinuous_space

Return the associated discontinuous scalar space.

domain_dimension

Dimension of the domain on which the space is defined.

dtype

Return the data type of the basis functions in the space.

evaluate_local_basis(element, local_coordinates, local_coefficients)

Evaluate a local basis on a given element.

evaluate_surface_gradient(element, local_coordinates, local_coefficients)

Evaluate the local surface gradient on a given element.

evaluation_functor

Functor transformation of shapesets from the reference element.

flat_local_dof_count

Return the total number of local degrees of freedom.

get_global_dofs(element, dof_weights=False)

Global dofs associated with the local dofs on the given element.

If dof_weights=True also return the associated dof_weights

global_dof_count

Return the number of global degrees of freedom.

global_dof_interpolation_points

(3xN) matrix of the N global interpolation points for the space.

global_dof_normals

(3xN) matrix of normal directions at interp. points.

global_to_local_dofs(global_dofs)

Local dofs and weights for the given list of global dofs.

grid

Return the underlying grid of the space.

has_local_support

True if support of each basis fct. is restricted to one element.

has_non_barycentric_space

True if space has an equivalent non-barycentric space.

inverse_mass_matrix()

Return the inverse mass matrix for this space.

is_barycentric

True if the space is defined over a barycentric refinement.

is_compatible(other)

Return true if spaces have the same number of global dofs.

is_identical(other)

Return true of spaces are identical.

mass_matrix()

Return the mass matrix associated with this space.

non_barycentric_space

The associated non-barycentric space if it exists.

order

Return the order of the space.

shapeset(element)

Return the Shapeset associated with a given element.

super_space

A super space that contains the space as subspace.

This attribute is useful for the compressed assembly of operators. The super space is typically a more local space (often discontinuous, elementwise defined) that is easier to compress than the actual space. The actual space is then obtained via projection.