BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
Classes | Functions
Helmholtz equation in 3D

Classes

class  Bempp::Helmholtz3dDoubleLayerPotentialOperator< BasisFunctionType_ >
 Double-layer potential operator for the Helmholtz equation in 3D. More...
 
class  Bempp::Helmholtz3dFarFieldDoubleLayerPotentialOperator< BasisFunctionType_ >
 Part of the far-field pattern of a radiating solution of the Helmholtz equation in 3D. More...
 
class  Bempp::Helmholtz3dFarFieldSingleLayerPotentialOperator< BasisFunctionType_ >
 Part of the far-field pattern of a radiating solution of the Helmholtz equation in 3D. More...
 
class  Bempp::Helmholtz3dPotentialOperatorBase< Impl, BasisFunctionType_ >
 Base class for potential operators for the Helmholtz equation in 3D. More...
 
class  Bempp::Helmholtz3dSingleLayerPotentialOperator< BasisFunctionType_ >
 Single-layer potential operator for the Helmholtz equation in 3D. More...
 
class  Bempp::ModifiedHelmholtz3dDoubleLayerPotentialOperator< BasisFunctionType_ >
 Double-layer potential operator for the modified Helmholtz equation in 3D. More...
 
class  Bempp::ModifiedHelmholtz3dPotentialOperatorBase< Impl, BasisFunctionType_ >
 Base class for potential operators for the modified Helmholtz equation in 3D. More...
 
class  Bempp::ModifiedHelmholtz3dSingleLayerPotentialOperator< BasisFunctionType_ >
 Single-layer potential operator for the modified Helmholtz equation in 3D. More...
 

Functions

template<typename BasisFunctionType >
BoundaryOperator
< BasisFunctionType, typename
ScalarTraits
< BasisFunctionType >
::ComplexType > 
Bempp::helmholtz3dAdjointDoubleLayerBoundaryOperator (const shared_ptr< const Context< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > > &context, const shared_ptr< const Space< BasisFunctionType > > &domain, const shared_ptr< const Space< BasisFunctionType > > &range, const shared_ptr< const Space< BasisFunctionType > > &dualToRange, typename ScalarTraits< BasisFunctionType >::ComplexType waveNumber, const std::string &label="", int symmetry=NO_SYMMETRY, bool useInterpolation=false, int interpPtsPerWavelength=DEFAULT_HELMHOLTZ_INTERPOLATION_DENSITY)
 Construct a BoundaryOperator object representing the adjoint double-layer boundary operator associated with the Helmholtz equation in 3D. More...
 
template<typename BasisFunctionType >
BoundaryOperator
< BasisFunctionType, typename
ScalarTraits
< BasisFunctionType >
::ComplexType > 
Bempp::helmholtz3dDoubleLayerBoundaryOperator (const shared_ptr< const Context< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > > &context, const shared_ptr< const Space< BasisFunctionType > > &domain, const shared_ptr< const Space< BasisFunctionType > > &range, const shared_ptr< const Space< BasisFunctionType > > &dualToRange, typename ScalarTraits< BasisFunctionType >::ComplexType waveNumber, const std::string &label="", int symmetry=NO_SYMMETRY, bool useInterpolation=false, int interpPtsPerWavelength=DEFAULT_HELMHOLTZ_INTERPOLATION_DENSITY)
 Construct a BoundaryOperator object representing the double-layer boundary operator associated with the Helmholtz equation in 3D. More...
 
template<typename BasisFunctionType >
BoundaryOperator
< BasisFunctionType, typename
ScalarTraits
< BasisFunctionType >
::ComplexType > 
Bempp::helmholtz3dHypersingularBoundaryOperator (const shared_ptr< const Context< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > > &context, const shared_ptr< const Space< BasisFunctionType > > &domain, const shared_ptr< const Space< BasisFunctionType > > &range, const shared_ptr< const Space< BasisFunctionType > > &dualToRange, typename ScalarTraits< BasisFunctionType >::ComplexType waveNumber, const std::string &label="", int symmetry=NO_SYMMETRY, bool useInterpolation=false, int interpPtsPerWavelength=DEFAULT_HELMHOLTZ_INTERPOLATION_DENSITY)
 Construct a BoundaryOperator object representing the hypersingular boundary operator associated with the Helmholtz equation in 3D. More...
 
template<typename BasisFunctionType >
BoundaryOperator
< BasisFunctionType, typename
ScalarTraits
< BasisFunctionType >
::ComplexType > 
Bempp::helmholtz3dSingleLayerBoundaryOperator (const shared_ptr< const Context< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > > &context, const shared_ptr< const Space< BasisFunctionType > > &domain, const shared_ptr< const Space< BasisFunctionType > > &range, const shared_ptr< const Space< BasisFunctionType > > &dualToRange, typename ScalarTraits< BasisFunctionType >::ComplexType waveNumber, const std::string &label="", int symmetry=NO_SYMMETRY, bool useInterpolation=false, int interpPtsPerWavelength=DEFAULT_HELMHOLTZ_INTERPOLATION_DENSITY)
 Construct a BoundaryOperator object representing the single-layer boundary operator associated with the Helmholtz equation in 3D. More...
 

Detailed Description

This submodule contains classes implementing kernels, boundary operators and potential operators related to the Helmholtz equation in 3D,

\[ \biggl(\frac{\partial^2}{\partial x^2} + \frac{\partial^2}{\partial y^2} + \frac{\partial^2}{\partial z^2} + k^2\biggr) u(x, y, z) = 0. \]

The number $k$ is referred to as the wave number.

Note
The term wave number refers to different physical quantities in the Helmholtz equation and in the modified Helmholtz equation.

Function Documentation

template<typename BasisFunctionType >
BoundaryOperator< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > Bempp::helmholtz3dAdjointDoubleLayerBoundaryOperator ( const shared_ptr< const Context< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > > &  context,
const shared_ptr< const Space< BasisFunctionType > > &  domain,
const shared_ptr< const Space< BasisFunctionType > > &  range,
const shared_ptr< const Space< BasisFunctionType > > &  dualToRange,
typename ScalarTraits< BasisFunctionType >::ComplexType  waveNumber,
const std::string &  label = "",
int  symmetry = NO_SYMMETRY,
bool  useInterpolation = false,
int  interpPtsPerWavelength = DEFAULT_HELMHOLTZ_INTERPOLATION_DENSITY 
)

Construct a BoundaryOperator object representing the adjoint double-layer boundary operator associated with the Helmholtz equation in 3D.

Parameters
[in]contextA Context object that will be used to build the weak form of the boundary operator when necessary.
[in]domainFunction space being the domain of the boundary operator.
[in]rangeFunction space being the range of the boundary operator.
[in]dualToRangeFunction space dual to the the range of the boundary operator.
[in]waveNumberWave number. See Helmholtz equation in 3D for its definition.
[in]labelTextual label of the operator. If empty, a unique label is generated automatically.
[in]symmetrySymmetry of the weak form of the operator. Can be any combination of the flags defined in the enumeration type Symmetry.
[in]useInterpolationIf set to false (default), the standard exp() function will be used to evaluate the exponential factor occurring in the kernel. If set to true, the exponential factor will be evaluated by piecewise-cubic interpolation of values calculated in advance on a regular grid. This normally speeds up calculations, but might result in a loss of accuracy. This is an experimental feature: use it at your own risk.
[in]interpPtsPerWavelengthIf useInterpolation is set to true, this parameter determines the number of points per "effective wavelength" (defined as $2\pi/|k|$, where $k$ = waveNumber) used to construct the interpolation grid. The default value (5000) is normally enough to reduce the relative or absolute error, whichever is smaller, below 100 * machine precision. If useInterpolation is set to false, this parameter is ignored.

None of the shared pointers may be null and the spaces range and dualToRange must be defined on the same grid, otherwise an exception is thrown.

If local-mode ACA assembly is requested (see AcaOptions::mode), after discretization, the weak form of this operator is stored as the product

\[ P A_{\textrm{d}} Q, \]

where $A_{\textrm{d}}$ is the weak form of this operator discretized with test and trial functions being the restrictions of the basis functions of domain and range to individual elements; $Q$ is the sparse matrix representing the expansion of the basis functions of domain in the just mentioned single-element trial functions; and $P$ is the sparse matrix whose transpose represents the expansion of the basis functions of dualToRange in the single-element test functions.

Template Parameters
BasisFunctionTypeType of the values of the basis functions into which functions acted upon by the operator are expanded. It can take the following values: float, double, std::complex<float> and std::complex<double>.

References Bempp::modifiedHelmholtz3dAdjointDoubleLayerBoundaryOperator().

template<typename BasisFunctionType >
BoundaryOperator< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > Bempp::helmholtz3dDoubleLayerBoundaryOperator ( const shared_ptr< const Context< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > > &  context,
const shared_ptr< const Space< BasisFunctionType > > &  domain,
const shared_ptr< const Space< BasisFunctionType > > &  range,
const shared_ptr< const Space< BasisFunctionType > > &  dualToRange,
typename ScalarTraits< BasisFunctionType >::ComplexType  waveNumber,
const std::string &  label = "",
int  symmetry = NO_SYMMETRY,
bool  useInterpolation = false,
int  interpPtsPerWavelength = DEFAULT_HELMHOLTZ_INTERPOLATION_DENSITY 
)

Construct a BoundaryOperator object representing the double-layer boundary operator associated with the Helmholtz equation in 3D.

Parameters
[in]contextA Context object that will be used to build the weak form of the boundary operator when necessary.
[in]domainFunction space being the domain of the boundary operator.
[in]rangeFunction space being the range of the boundary operator.
[in]dualToRangeFunction space dual to the the range of the boundary operator.
[in]waveNumberWave number. See Helmholtz equation in 3D for its definition.
[in]labelTextual label of the operator. If empty, a unique label is generated automatically.
[in]symmetrySymmetry of the weak form of the operator. Can be any combination of the flags defined in the enumeration type Symmetry.
[in]useInterpolationIf set to false (default), the standard exp() function will be used to evaluate the exponential factor occurring in the kernel. If set to true, the exponential factor will be evaluated by piecewise-cubic interpolation of values calculated in advance on a regular grid. This normally speeds up calculations, but might result in a loss of accuracy. This is an experimental feature: use it at your own risk.
[in]interpPtsPerWavelengthIf useInterpolation is set to true, this parameter determines the number of points per "effective wavelength" (defined as $2\pi/|k|$, where $k$ = waveNumber) used to construct the interpolation grid. The default value (5000) is normally enough to reduce the relative or absolute error, whichever is smaller, below 100 * machine precision. If useInterpolation is set to false, this parameter is ignored.

None of the shared pointers may be null and the spaces range and dualToRange must be defined on the same grid, otherwise an exception is thrown.

If local-mode ACA assembly is requested (see AcaOptions::mode), after discretization, the weak form of this operator is stored as the product

\[ P A_{\textrm{d}} Q, \]

where $A_{\textrm{d}}$ is the weak form of this operator discretized with test and trial functions being the restrictions of the basis functions of domain and range to individual elements; $Q$ is the sparse matrix representing the expansion of the basis functions of domain in the just mentioned single-element trial functions; and $P$ is the sparse matrix whose transpose represents the expansion of the basis functions of dualToRange in the single-element test functions.

Template Parameters
BasisFunctionTypeType of the values of the basis functions into which functions acted upon by the operator are expanded. It can take the following values: float, double, std::complex<float> and std::complex<double>.

References Bempp::modifiedHelmholtz3dDoubleLayerBoundaryOperator().

template<typename BasisFunctionType >
BoundaryOperator< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > Bempp::helmholtz3dHypersingularBoundaryOperator ( const shared_ptr< const Context< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > > &  context,
const shared_ptr< const Space< BasisFunctionType > > &  domain,
const shared_ptr< const Space< BasisFunctionType > > &  range,
const shared_ptr< const Space< BasisFunctionType > > &  dualToRange,
typename ScalarTraits< BasisFunctionType >::ComplexType  waveNumber,
const std::string &  label = "",
int  symmetry = NO_SYMMETRY,
bool  useInterpolation = false,
int  interpPtsPerWavelength = DEFAULT_HELMHOLTZ_INTERPOLATION_DENSITY 
)

Construct a BoundaryOperator object representing the hypersingular boundary operator associated with the Helmholtz equation in 3D.

Parameters
[in]contextA Context object that will be used to build the weak form of the boundary operator when necessary.
[in]domainFunction space being the domain of the boundary operator.
[in]rangeFunction space being the range of the boundary operator.
[in]dualToRangeFunction space dual to the the range of the boundary operator.
[in]waveNumberWave number. See Helmholtz equation in 3D for its definition.
[in]labelTextual label of the operator. If empty, a unique label is generated automatically.
[in]symmetrySymmetry of the weak form of the operator. Can be any combination of the flags defined in the enumeration type Symmetry.
[in]useInterpolationIf set to false (default), the standard exp() function will be used to evaluate the exponential factor occurring in the kernel. If set to true, the exponential factor will be evaluated by piecewise-cubic interpolation of values calculated in advance on a regular grid. This normally speeds up calculations, but might result in a loss of accuracy. This is an experimental feature: use it at your own risk.
[in]interpPtsPerWavelengthIf useInterpolation is set to true, this parameter determines the number of points per "effective wavelength" (defined as $2\pi/|k|$, where $k$ = waveNumber) used to construct the interpolation grid. The default value (5000) is normally enough to reduce the relative or absolute error, whichever is smaller, below 100 * machine precision. If useInterpolation is set to false, this parameter is ignored.

None of the shared pointers may be null and the spaces range and dualToRange must be defined on the same grid, otherwise an exception is thrown.

If local-mode ACA assembly is requested (see AcaOptions::mode), after discretization, the weak form of this operator is stored as

\[ A = \sum_{i=1}^3 P_i A_{\textrm{d}} Q_i - k^2 \sum_{i=1}^3 R_i A_{\textrm{d}} S_i, \]

where:

  • $A_{\textrm{d}}$ is the weak form of the single-layer modified Helmholtz boundary operator discretised with test and trial functions being the restrictions of the basis functions of domain and range to individual elements;
  • $P_i$ is the sparse matrix whose transpose represents the expansion of the $i$th component of the surface curl of the basis functions of dualToRange in the single-element test functions mentioned above;
  • $Q_i$ is the sparse matrix representing the expansion of the $i$th component of the surface curl of the basis functions of domain in the single-element trial functions;
  • $k$ is the wave number
  • $R_i$ is the sparse matrix whose transpose represents the expansion of the basis functions of dualToRange, multiplied by the $i$th component of the local vector normal to the surface on which functions from dualToRange live, in the single-element test functions;
  • $S_i$ is the sparse matrix representing the expansion of basis functions of domain, multiplied by the $i$th component of the local vector normal to the surface on which functions from domain live, in the single-element trial functions.
Template Parameters
BasisFunctionTypeType of the values of the basis functions into which functions acted upon by the operator are expanded. It can take the following values: float, double, std::complex<float> and std::complex<double>.

References Bempp::modifiedHelmholtz3dHypersingularBoundaryOperator().

template<typename BasisFunctionType >
BoundaryOperator< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > Bempp::helmholtz3dSingleLayerBoundaryOperator ( const shared_ptr< const Context< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > > &  context,
const shared_ptr< const Space< BasisFunctionType > > &  domain,
const shared_ptr< const Space< BasisFunctionType > > &  range,
const shared_ptr< const Space< BasisFunctionType > > &  dualToRange,
typename ScalarTraits< BasisFunctionType >::ComplexType  waveNumber,
const std::string &  label = "",
int  symmetry = NO_SYMMETRY,
bool  useInterpolation = false,
int  interpPtsPerWavelength = DEFAULT_HELMHOLTZ_INTERPOLATION_DENSITY 
)

Construct a BoundaryOperator object representing the single-layer boundary operator associated with the Helmholtz equation in 3D.

Parameters
[in]contextA Context object that will be used to build the weak form of the boundary operator when necessary.
[in]domainFunction space being the domain of the boundary operator.
[in]rangeFunction space being the range of the boundary operator.
[in]dualToRangeFunction space dual to the the range of the boundary operator.
[in]waveNumberWave number. See Helmholtz equation in 3D for its definition.
[in]labelTextual label of the operator. If empty, a unique label is generated automatically.
[in]symmetrySymmetry of the weak form of the operator. Can be any combination of the flags defined in the enumeration type Symmetry.
[in]useInterpolationIf set to false (default), the standard exp() function will be used to evaluate the exponential factor occurring in the kernel. If set to true, the exponential factor will be evaluated by piecewise-cubic interpolation of values calculated in advance on a regular grid. This normally speeds up calculations, but might result in a loss of accuracy. This is an experimental feature: use it at your own risk.
[in]interpPtsPerWavelengthIf useInterpolation is set to true, this parameter determines the number of points per "effective wavelength" (defined as $2\pi/|k|$, where $k$ = waveNumber) used to construct the interpolation grid. The default value (5000) is normally enough to reduce the relative or absolute error, whichever is smaller, below 100 * machine precision. If useInterpolation is set to false, this parameter is ignored.

None of the shared pointers may be null and the spaces range and dualToRange must be defined on the same grid, otherwise an exception is thrown.

If local-mode ACA assembly is requested (see AcaOptions::mode), after discretization, the weak form of this operator is stored as the product

\[ P A_{\textrm{d}} Q, \]

where $A_{\textrm{d}}$ is the weak form of this operator discretized with test and trial functions being the restrictions of the basis functions of domain and range to individual elements; $Q$ is the sparse matrix representing the expansion of the basis functions of domain in the just mentioned single-element trial functions; and $P$ is the sparse matrix whose transpose represents the expansion of the basis functions of dualToRange in the single-element test functions.

Template Parameters
BasisFunctionTypeType of the values of the basis functions into which functions acted upon by the operator are expanded. It can take the following values: float, double, std::complex<float> and std::complex<double>.

References Bempp::modifiedHelmholtz3dSingleLayerBoundaryOperator().