Available operators

In the following we describe all operators available in BEM++.

Scalar non-local operators (Laplace, modified Helmholtz and Helmholtz)

Let \(g(x,y)\) be a given Green’s function. We define the single-layer potential operator and double layer-potential operator as follows:

\[\begin{split}[\mathcal{S}\phi](x) &= \int_{\Gamma}g(x,y)\phi(y)ds(y),~x\in\mathbb{R}^3\backslash\Gamma\nonumber\\ [\mathcal{K}\phi](x) &= \int_{\Gamma}\frac{\partial g(x,y)}{\partial\nu(y)}\phi(y)ds(y),~x\in\mathbb{R}^3\backslash\Gamma\end{split}\]

From this we derive the following boundary operators:

  • Single-Layer boundary operator: \([S\phi](x) = \int_{\Gamma}g(x,y)\phi(y)ds(y),~x\in\Gamma\)
  • Double-Layer boundary operator: \([K\phi](x) = \int_{\Gamma}\frac{\partial g(x,y)}{\partial\nu(y)}\phi(y)ds(y),~x\in\Gamma\)
  • Adjoint double layer boundary operator: \([K'\phi](x) = \int_{\Gamma}\frac{\partial g(x,y)}{\partial\nu(x)}\phi(y)ds(y),~x\in\Gamma\)
  • Hypersingular boundary operator: \([D\phi](x) = -\frac{\partial}{\partial \nu(x)}\int_{\Gamma}\frac{\partial g(x,y)}{\partial\nu(y)}\phi(y)ds(y),~x\in\Gamma\)

The actual implementation of boundary operators in BEM++ are based on variational formulations. Hence, the implementation for example of the single-layer boundary operator is given by

\[s(\phi, \psi) = \int_{\Gamma}\overline{\psi(x)}\int_{\Gamma}g(x,y)\phi(y)ds(y)ds(x).\]

For the hypersingular operator integration a slightly different formulation based on integration by parts of this variational form is used. Details can be found for example in the book Numerical Approximation Methods for Elliptic Boundary Value Problems by O. Steinbach.

The different types of boundary operators are dependent on the choice of the Green’s function \(g\). In BEM++ the definitions are as follows.

  • Laplace (\(-\Delta u = 0\)):

    \[g(x,y) = \frac{1}{4\pi |x-y|}\]
  • Modified Helmholtz (\(-\Delta u+\omega^2 u = 0\)):

    \[g(x,y) = \frac{e^{-\omega |x-y|}}{4\pi|x-y|}\]
  • Helmholtz (\(-\Delta u - k^2 u = 0\)):

    \[g(x,y) = \frac{e^{ik |x-y|}}{4\pi|x-y|}\]

The boundary operators are defined in the following modules:

  • Laplace: bempp.api.operators.boundary.laplace
  • Helmholtz: bempp.api.operators.boundary.helmholtz
  • Modified Helmholtz: bempp.api.operators.boundary.modified_helmholtz

In each of these modules the names of the corrresponding functions are single_layer, double_layer, adjoint_double_layer and hypersingular.

The packages for the corresponding potential operators are:

  • Laplace: bempp.api.operators.potential.laplace
  • Helmholtz: bempp.api.operators.potential.helmholtz
  • Modified Helmholtz: bempp.api.operators.potential.modified_helmholtz

The names of the corresponding potentials in each module are single_layer and double_layer.

For example, a Helmholtz hypersingular boundary operator is instantiated as

hyp = bempp.api.operators.boundary.helmholtz.hypersingular(domain, range, dual_to_range, k)

Here, \(k\) is the wavenumber as defined above.

In several applications, in particular transmission problems it is useful to define the multitrace operator

\[\begin{split}A = \begin{bmatrix} -K & S\\ D & K'\end{bmatrix}.\end{split}\]

This operator has the property that

\[\begin{split}\left[\frac{1}{2}I + A\right]\begin{bmatrix}u \\ u_n \end{bmatrix} = \begin{bmatrix}u \\ u_n\end{bmatrix}\end{split}\]

if \(u\) and \(u_n\) are the Dirichlet, respectively Neumann data, of an interior solution to the corresponding PDE. The operator \(\frac{1}{2}I + A\) is also called a Calderon projector. From this property it follows directly that \((2A)^2 = I\) and therefore \(A\) is self-regularizing.

Currently, the operator \(A\) is implemented using linear, piecewise continuous basis functions for the Dirichlet data and piecewise constant basis functions on a dual grid for the Neumann data. This choice of basis functions allows to perform the product \(A * A\) in a stable way.

The following defines a multitrace operator for a Helmholtz problem with wavenumber \(k=1\) and forms the interior Calderon projector \(\frac{1}{2}I+A\).

grid = bempp.api.shapes.regular_sphere(3)
A = bempp.api.operators.helmholtz.multitrace_operator(grid, 1)
ident = bempp.api.operators.boundary.sparse.multitrace_identity(grid)
calderon = .5 * ident + A

The assembly of the multitrace operator \(A\) requires the assembly of the single-layer and double-layer boundary operators using discontinuous piecewise linear basis functions on a barycentric refinement of the grid. This barycentric refinement has six times the number of elements as the original grid. This is necessary to assemble a basis for piecewise constant basis functions on the dual grid.

Scalar, sparse operators

BEM++ implements two sparse operators acting on scalar spaces, the identity operator \(I\) and the Laplace-Beltrami operator \(-\Delta_{\Gamma}\). They are given in their variational form as follows.

  • Identity Operator:

    \[m(\phi, \psi) = \int_{\Gamma} \overline{\Psi(y)} \phi(y)dy\]
  • Laplace-Beltrami Operator:

    \[k(\phi, \psi) = \int_{\Gamma} \overline{\nabla_{\Gamma}\Psi(y)}\cdot \nabla_{\Gamma}\phi(y)dy\]

The corresponding functions are bempp.api.boundary.sparse.identity and bempp.api.boundary.sparse.laplace_beltrami.

Maxwell operators

BEM++ supports the solution of the time-harmonic Maxwell equation of the form

\[\nabla\times\nabla\times u - k^2 u = 0.\]

We define the following two potentials:

  • Maxwell electric field potential operator:

    \[[\mathcal{E}\phi](x) = ik\int_{\Gamma}g(x,y)\phi(y)ds(y)-\frac{1}{ik}\nabla_x\int_{\Gamma}g(x,y)(\nabla_{\Gamma}\cdot\phi)(y)ds(y)\]
  • Maxwell magnetic field potential operator:

    \[[\mathcal{M}\phi](x) = \nabla_x\times\int_{\Gamma}g(x,y)\phi(y)ds(y)\]

    The corresponding functions are bempp.api.operators.potential.maxwell.electric_field and bempp.api.operators.potential.maxwell.magnetic_field.

The definition of the electric field operator given above includes a factor \(i\) in the nominator. This is less common, but has implementational advantages in BEM++.

Based on the electric and magnetic field potential operators we can derive the corresponding boundary operators. We will not give details but just state the variational formulations of the operators.

  • Maxwell electric field boundary operator:

    \[s(\phi, \psi) = \int_{\Gamma}\int_{\Gamma}g(x,y)\left[-ik\overline{\psi(x)}\cdot\phi(y)-\frac{1}{ik}\left(\nabla_{\Gamma}\cdot\overline{\psi}\right)(x)\left(\nabla_{\Gamma}\cdot\phi\right)(y)\right]ds(x)ds(y)\]
  • Maxwell magnetic field boundary operator:

    \[c(\phi, \psi) = \int_{\Gamma}\int_{\Gamma}\nabla_xg(x,y)\left[\overline{\psi(x)}\times \phi(y)\right]ds(x)ds(y)\]

    The corresponding BEM++ functions are bempp.operators.boundary.maxwell.electric_field and bempp.operators.boundary.maxwell.magnetic_field.

The discretisation of Maxwell operators uses an antisymmetric pseudo-inner product defined as

\[i(\phi, \psi) = \int_{\Gamma}\overline{\psi(y)}\left(\phi(y)\times \nu(y)\right)ds(y).\]

For Maxwell problems this should be used instead of the standard identity operator. The corresponding function in BEM++ is bempp.operators.boundary.sparse.maxwell_identity.

Far Field Operators

BEM++ implements far field operators for Maxwell and Helmholtz problems.

Approximate DtN and NtD operators

BEM++ implements the OSRC approximations to the Dirichlet-to-Neumann (DtN) and Neumann-to-Dirichlet (NtD) maps (see e.g. X. Antoine, M. Darbas, Generalized Combined Field Integral Equations for the Iterative Solution of the Three-Dimensional Helmholtz Equation, Mathematical Modelling and Numerical Analysis 41 (1) (2007), pp. 147-167). The corresponding functions are bempp.api.boundary.helmholtz.osrc_ntd and bempp.api.boundary.helmholtz.osrc_dtn.

Function and class reference (Laplace)

bempp.api.operators.boundary.laplace.single_layer(domain, range_, dual_to_range, label='SLP', symmetry='no_symmetry', parameters=None, use_projection_spaces=True, assemble_only_singular_part=False)

Return the Laplace single-layer boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false).
bempp.api.operators.boundary.laplace.double_layer(domain, range_, dual_to_range, label='DLP', symmetry='no_symmetry', parameters=None, use_projection_spaces=True, assemble_only_singular_part=False)

Return the Laplace double-layer boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false).
bempp.api.operators.boundary.laplace.adjoint_double_layer(domain, range_, dual_to_range, label='ADJ_DLP', symmetry='no_symmetry', parameters=None, use_projection_spaces=True, assemble_only_singular_part=False)

Return the Laplace adjoint double-layer boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false).
bempp.api.operators.boundary.laplace.hypersingular(domain, range_, dual_to_range, label='HYP', symmetry='no_symmetry', parameters=None, use_slp=False, use_projection_spaces=True, assemble_only_singular_part=False)

Return the Laplace hypersingular boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • use_slp (True/False or boundary operator object) – The hypersingular operator can be represented as a sparse transformation of a single-layer operator. If use_slp=True this representation is used. If use_slp=op for a single-layer boundary operator assembled on a suitable space this operator is used to assemble the hypersingular operator. Note that if use_slp=op is used no checks are performed if the slp operator is correctly defined for representing the hypersingular operator. Hence, if no care is taken this option can lead to a wrong operator. Also, use_slp=True or use_slp=op is only valid if the domain and dual_to_range spaces are identical.
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false). Note. This option is only used if use_slp is not specified.
bempp.api.operators.boundary.laplace.multitrace_operator(grid, parameters=None, spaces='linear', target=None)

Return the Laplace multitrace operator.

Parameters:
  • grid (bempp.api.grid.Grid) – The underlying grid for the multitrace operator
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • spaces (string) – Choose ‘linear’ to assemble the operator with continuous linear function spaces for the Dirichlet and Neumann component (default). For a dual pairing of a linear space for the Dirichlet data and piecewise constant space for the Neumann data choose ‘dual’.
  • target (bempp.api.grid.Grid) – Specifies a target grid. If it is different from ‘grid’ then the operator maps from ‘grid’ to ‘target’.
bempp.api.operators.potential.laplace.single_layer(*args, **kwargs)

Wrapper for potentials to emit log messages.

bempp.api.operators.potential.laplace.double_layer(*args, **kwargs)

Wrapper for potentials to emit log messages.

Function and class reference (Modified Helmholtz)

bempp.api.operators.boundary.modified_helmholtz.single_layer(domain, range_, dual_to_range, wave_number, label='SLP', symmetry='no_symmetry', parameters=None, use_projection_spaces=True, assemble_only_singular_part=False)

Return the modified Helmholtz single-layer boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • wave_number (complex) – Wavenumber for the Helmholtz problem.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false).
bempp.api.operators.boundary.modified_helmholtz.double_layer(domain, range_, dual_to_range, wave_number, label='DLP', symmetry='no_symmetry', parameters=None, use_projection_spaces=True, assemble_only_singular_part=False)

Return the modified Helmholtz double-layer boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • wave_number (complex) – Wavenumber for the Helmholtz problem.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false).
bempp.api.operators.boundary.modified_helmholtz.adjoint_double_layer(domain, range_, dual_to_range, wave_number, label='ADJ_DLP', symmetry='no_symmetry', parameters=None, use_projection_spaces=True, assemble_only_singular_part=False)

Return the modified Helmholtz adjoint double-layer boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • wave_number (complex) – Wavenumber for the Helmholtz problem.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false).
bempp.api.operators.boundary.modified_helmholtz.hypersingular(domain, range_, dual_to_range, wave_number, label='HYP', symmetry='no_symmetry', parameters=None, use_slp=False, use_projection_spaces=True, assemble_only_singular_part=False)

Return the modified Helmholtz hypersingular boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • wave_number (complex) – Wavenumber for the Helmholtz problem.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • use_slp (True/False or boundary operator object) – The hypersingular operator can be represented as a sparse transformation of a single-layer operator. If use_slp=True this representation is used. If use_slp=op for a single-layer boundary operator assembled on a suitable space this operator is used to assemble the hypersingular operator. Note that if use_slp=op is used no checks are performed if the slp operator is correctly defined for representing the hypersingular operator. Hence, if no care is taken this option can lead to a wrong operator. Also, use_slp=True or use_slp=op is only valid if the domain and dual_to_range spaces are identical.
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false). Note. This option is only used if use_slp is not specified.
bempp.api.operators.boundary.modified_helmholtz.multitrace_operator(grid, wave_number, parameters=None, spaces='linear', target=None)

Return the modified Helmholtz multitrace operator.

Parameters:
  • grid (bempp.api.grid.Grid) – The underlying grid for the multitrace operator
  • wave_number (complex) – Wavenumber of the operator.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • spaces (string) – Choose ‘linear’ to assemble the operator with continuous linear function spaces for the Dirichlet and Neumann component (default). For a dual pairing of a linear space for the Dirichlet data and piecewise constant space for the Neumann data choose ‘dual’.
  • target (bempp.api.grid.Grid) – Specifies a target grid. If it is different from ‘grid’ then the operator maps from ‘grid’ to ‘target’.
bempp.api.operators.potential.modified_helmholtz.single_layer(*args, **kwargs)

Wrapper for potentials to emit log messages.

bempp.api.operators.potential.modified_helmholtz.double_layer(*args, **kwargs)

Wrapper for potentials to emit log messages.

Function and class reference (Helmholtz)

bempp.api.operators.boundary.helmholtz.single_layer(domain, range_, dual_to_range, wave_number, label='SLP', symmetry='no_symmetry', parameters=None, use_projection_spaces=True, assemble_only_singular_part=False)

Return the Helmholtz single-layer boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • wave_number (complex) – Wavenumber for the Helmholtz problem.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false).
bempp.api.operators.boundary.helmholtz.double_layer(domain, range_, dual_to_range, wave_number, label='DLP', symmetry='no_symmetry', parameters=None, use_projection_spaces=True, assemble_only_singular_part=False)

Return the Helmholtz double-layer boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • wave_number (complex) – Wavenumber for the Helmholtz problem.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false).
bempp.api.operators.boundary.helmholtz.adjoint_double_layer(domain, range_, dual_to_range, wave_number, label='ADJ_DLP', symmetry='no_symmetry', parameters=None, use_projection_spaces=True, assemble_only_singular_part=False)

Return the Helmholtz adjoint double-layer boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • wave_number (complex) – Wavenumber for the Helmholtz problem.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false).
bempp.api.operators.boundary.helmholtz.hypersingular(domain, range_, dual_to_range, wave_number, label='HYP', symmetry='no_symmetry', parameters=None, use_slp=False, use_projection_spaces=True, assemble_only_singular_part=False)

Return the Helmholtz hypersingular boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • wave_number (complex) – Wavenumber for the Helmholtz problem.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • use_slp (True/False or boundary operator object) – The hypersingular operator can be represented as a sparse transformation of a single-layer operator. If use_slp=True this representation is used. If use_slp=op for a single-layer boundary operator assembled on a suitable space this operator is used to assemble the hypersingular operator. Note that if use_slp=op is used no checks are performed if the slp operator is correctly defined for representing the hypersingular operator. Hence, if no care is taken this option can lead to a wrong operator. Also, use_slp=True or use_slp=op is only valid if the domain and dual_to_range spaces are identical.
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false). Note. This option is only used if use_slp is not specified.
bempp.api.operators.boundary.helmholtz.multitrace_operator(grid, wave_number, parameters=None, spaces='linear', target=None)

Return the Helmholtz multitrace operator.

Parameters:
  • grid (bempp.api.grid.Grid) – The underlying grid for the multitrace operator
  • wave_number (complex) – Wavenumber of the operator.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • spaces (string) – Choose ‘linear’ to assemble the operator with continuous linear function spaces for the Dirichlet and Neumann component (default). For a dual pairing of a linear space for the Dirichlet data and piecewise constant space for the Neumann data choose ‘dual’.
  • target (bempp.api.grid.Grid) – Specifies a target grid. If it is different from ‘grid’ then the operator maps from ‘grid’ to ‘target’.
bempp.api.operators.boundary.helmholtz.osrc_dtn(space, wave_number, npade=2, theta=1.0471975511965976, parameters=None, label='osrc_dtn', damped_wavenumber=None)

Return the OSRC approximation to the DtN map.

Parameters:
  • space (bempp.api.space.Space) – Underlying function space. Needs to be a continuous space.
  • wave_number (complex) – Wavenumber of the operator.
  • npade (int) – Number of Pade terms to be used in the approximation
  • theta (float64) – Angle of the branch cut of the square root operator.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • label (string) – Label for the operator.
  • damped_wavenumber (complex) – Explicitly specify the damped wavenumber of the OSRC operator.
bempp.api.operators.boundary.helmholtz.osrc_ntd(space, wave_number, npade=2, theta=1.0471975511965976, parameters=None, label='osrc_ntd', damped_wavenumber=None)

Return the OSRC approximation to the NtD map.

Parameters:
  • space (bempp.api.space.Space) – Underlying function space. Needs to be a continuous space.
  • wave_number (complex) – Wavenumber of the operator.
  • npade (int) – Number of Pade terms to be used in the approximation
  • theta (float64) – Angle of the branch cut of the square root operator.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • label (string) – Label for the operator.
  • damped_wavenumber (complex) – Explicitly specify the damped wavenumber of the OSRC operator.
bempp.api.operators.potential.helmholtz.single_layer(space, evaluation_points, wave_number, parameters=None)

Return the Helmholtz single-layer potential operator

Parameters:
  • space (bempp.api.space.Space) – The function space over which to assemble the potential.
  • evaluation_points (numpy.ndarray) – A (3 x N) array of N evaluation points, where each column corresponds to the coordinates of one evaluation point.
  • wave_number (complex) – Wavenumber of the operator.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
bempp.api.operators.potential.helmholtz.double_layer(space, evaluation_points, wave_number, parameters=None)

Return the Helmholtz double-layer potential operator

Parameters:
  • space (bempp.api.space.Space) – The function space over which to assemble the potential.
  • evaluation_points (numpy.ndarray) – A (3 x N) array of N evaluation points, where each column corresponds to the coordinates of one evaluation point.
  • wave_number (complex) – Wavenumber of the operator.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
bempp.api.operators.far_field.helmholtz.single_layer(*args, **kwargs)

Wrapper for potentials to emit log messages.

bempp.api.operators.far_field.helmholtz.double_layer(*args, **kwargs)

Wrapper for potentials to emit log messages.

Function and class reference (Maxwell)

bempp.api.operators.boundary.maxwell.electric_field(domain, range_, dual_to_range, wave_number, label='EFIE', symmetry='no_symmetry', parameters=None, use_projection_spaces=True, assemble_only_singular_part=False)

Return the Maxwell electric field boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • wave_number (complex) – Wavenumber for the Helmholtz problem.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false).
bempp.api.operators.boundary.maxwell.magnetic_field(domain, range_, dual_to_range, wave_number, label='MFIE', symmetry='no_symmetry', parameters=None, use_projection_spaces=True, assemble_only_singular_part=False)

Return the Maxwell magnetic field boundary operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • wave_number (complex) – Wavenumber for the Helmholtz problem.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • use_projection_spaces (bool) – Represent operator by projection from higher dimensional space if available. This parameter can speed up fast assembly routines, such as H-Matrices or FMM (default true).
  • assemble_only_singular_part (bool) – When assembled the operator will only contain components for adjacent or overlapping test and trial functions (default false).
bempp.api.operators.potential.maxwell.electric_field(*args, **kwargs)

Wrapper for potentials to emit log messages.

bempp.api.operators.potential.maxwell.magnetic_field(*args, **kwargs)

Wrapper for potentials to emit log messages.

bempp.api.operators.far_field.maxwell.electric_field(*args, **kwargs)

Wrapper for potentials to emit log messages.

bempp.api.operators.far_field.maxwell.magnetic_field(*args, **kwargs)

Wrapper for potentials to emit log messages.

Function and class reference (sparse operators)

bempp.api.operators.boundary.sparse.identity(domain, range_, dual_to_range, label='IDENTITY', symmetry='no_symmetry', parameters=None)

Return the identity operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
bempp.api.operators.boundary.sparse.laplace_beltrami(domain, range_, dual_to_range, label='LAPLACE_BELTRAMI', symmetry='no_symmetry', parameters=None)

Return the Laplace Beltrami operator.

Parameters:
  • domain (bempp.api.space.Space) – Domain space.
  • range (bempp.api.space.Space) – Range space.
  • dual_to_range (bempp.api.space.Space) – Dual space to the range space.
  • label (string) – Label for the operator.
  • symmetry (string) – Symmetry mode. Possible values are: ‘no_symmetry’, ‘symmetric’, ‘hermitian’.
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.

Notes

The spaces for this operator must be spaces of continuous functions.

bempp.api.operators.boundary.sparse.multitrace_identity(grid, parameters=None, spaces='linear')

Return the multitrace identity operator.

Parameters:
  • grid (bempp.api.grid.Grid) – The underlying grid for the multitrace operator
  • parameters (bempp.api.common.ParameterList) – Parameters for the operator. If none given the default global parameter object bempp.api.global_parameters is used.
  • spaces (string) – Choose ‘linear’ to assemble the operator with continuous linear function spaces for the Dirichlet and Neumann component (default). For a dual pairing of a linear space for the Dirichlet data and piecewise constant space for the Neumann data choose ‘dual’. For the proper dual spaces in Maxwell problems choose ‘maxwell’.