|
BEM++
2.0
|
Classes | |
| class | Bempp::Maxwell3dDoubleLayerPotentialOperator< BasisFunctionType_ > |
| Double-layer potential operator for the Maxwell equations in 3D. More... | |
| class | Bempp::Maxwell3dFarFieldDoubleLayerPotentialOperator< BasisFunctionType_ > |
| Far-field limit of the double-layer potential operator for the Maxwell equations in 3D. More... | |
| class | Bempp::Maxwell3dFarFieldSingleLayerPotentialOperator< BasisFunctionType_ > |
| Far-field limit of the single-layer potential operator for the Maxwell equations in 3D. More... | |
| class | Bempp::Maxwell3dIdentityOperator< BasisFunctionType_, ResultType_ > |
| "Identity operator" for Maxwell equations in 3D. More... | |
| class | Bempp::Maxwell3dSingleLayerPotentialOperator< BasisFunctionType_ > |
| Single-layer potential operator for the Maxwell equations in 3D. More... | |
Functions | |
| template<typename BasisFunctionType > | |
| BoundaryOperator < BasisFunctionType, typename ScalarTraits < BasisFunctionType > ::ComplexType > | Bempp::maxwell3dDoubleLayerBoundaryOperator (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 double-layer boundary operator for Maxwell equations in 3D. More... | |
| template<typename BasisFunctionType > | |
| BoundaryOperator < BasisFunctionType, typename ScalarTraits < BasisFunctionType > ::ComplexType > | Bempp::maxwell3dSingleLayerBoundaryOperator (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 single-layer boundary operator associated with the Maxwell equations in 3D. More... | |
This submodule contains classes relevant to the implementation of boundary operators and potential operators related to the time-harmonic Maxwell equations in 3D. Inside a region with a uniform scalar permittivity and permeability, these equations can be reduced to the form
where
stands either for the electric or the magnetic field and the number
is referred to as the wave number.
The treatment of Maxwell equations in BEM++ closely follows that of A. Buffa and R. Hiptmair, "Galerkin Boundary Element Methods for Electromagnetic Scattering", in "Computational Methods in Wave Propagation", M. Ainsworth, ed., Springer, New York 2003, pp. 85-126. In the following we cite the principal definitions relevant to the implementation of BEM in the Maxwell case.
The starting point is the Stratton-Chu theorem (Theorem 6 in Buffa and Hiptmair), which plays a part analogous to the Green's representation theorem for the Laplace equation. It states that any solution
of the Maxwell equations in a bounded domain
with boundary
has the representation
(all the symbols will be defined shortly). Furthermore, any solution
that satisfies the Maxwell equations in a domain
extending to infinity and meets the Silver-Mueller radiation conditions at infinity has the representation
The action of the single-layer potential operator is defined by
where
is the Green's function of the Helmholtz equation. The action of the double-layer potential operator, in turn, is defined by
The interior Dirichlet trace
at a point
is defined as
where
is the outward unit vector normal to
at
and
is the limit of
as
approaches
from within
. The interior Neumann trace
at
is defined as
Note that compared to Buffa and Hiptmair, we include an additional
factor in the denominator of the Neumann trace and in the numerator of the single-layer potential operator, as this makes both quantities real-valued in the special case of
purely imaginary (and
purely real).
The exterior Dirichlet and Neumann traces are defined analogously, with the relevant quantities assumed to approach the point
from within the complement of
.
The single-layer boundary operator
is defined as the average of the interior and exterior Dirichlet traces of the single-layer potential operator:
The double-layer boundary operator
is defined, analogously, as the average of the interior and exterior Dirichlet traces of the double-layer potential operator.
To exploit the symmetry between the electric and magnetic field in Maxwell equations, it is advantageous to use the following antisymmetric pseudo-inner product in the Galerkin discretisation of the above boundary operators:
where
denotes complex conjugation. The weak forms of the single- and double-layer boundary operators with respect to this pseudo-inner product are given by
and
Solution of particular Maxwell problems relies on the following boundary integral equations imposed on
:
Here, the symbol
denotes the identity operator. It is crucial to remember that the weak forms of
and
are defined with respect to the pseudo-inner product rather than the standard sesquilinear inner product. Thus, in practical implementation an operator whose weak form represents the identity operator under the pseudo-inner product is indispensable. Such an operator is provided by the Maxwell3dIdentityOperator class.
| BoundaryOperator< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > Bempp::maxwell3dDoubleLayerBoundaryOperator | ( | 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 double-layer boundary operator for Maxwell equations in 3D.
This function constructs a BoundaryOperator object representing the double-layer boundary operator
for the Maxwell equations in 3D, as defined in Maxwell equations in 3D.
| [in] | context | A Context object that will be used to build the weak form of the boundary operator when necessary. |
| [in] | domain | Function space being the domain of the boundary operator. |
| [in] | range | Function space being the range of the boundary operator. |
| [in] | dualToRange | Function space dual to the the range of the boundary operator. |
| [in] | waveNumber | Wave number. See Maxwell equations in 3D for its definition. |
| [in] | label | Textual label of the operator. If empty, a unique label is generated automatically. |
| [in] | symmetry | Symmetry of the weak form of the operator. Can be any combination of the flags defined in the enumeration type Symmetry. |
| [in] | useInterpolation | If 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] | interpPtsPerWavelength | If useInterpolation is set to true, this parameter determines the number of points per "effective wavelength" (defined as , where = 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.
| BoundaryOperator< BasisFunctionType, typename ScalarTraits< BasisFunctionType >::ComplexType > Bempp::maxwell3dSingleLayerBoundaryOperator | ( | 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 single-layer boundary operator associated with the Maxwell equations in 3D.
This function constructs a BoundaryOperator object representing the single-layer boundary operator
for the Maxwell equations in 3D, as defined in Maxwell equations in 3D.
| [in] | context | A Context object that will be used to build the weak form of the boundary operator when necessary. |
| [in] | domain | Function space being the domain of the boundary operator. |
| [in] | range | Function space being the range of the boundary operator. |
| [in] | dualToRange | Function space dual to the the range of the boundary operator. |
| [in] | waveNumber | Wave number. See Maxwell equations in 3D for its definition. |
| [in] | label | Textual label of the operator. If empty, a unique label is generated automatically. |
| [in] | symmetry | Symmetry of the weak form of the operator. Can be any combination of the flags defined in the enumeration type Symmetry. |
| [in] | useInterpolation | If 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] | interpPtsPerWavelength | If useInterpolation is set to true, this parameter determines the number of points per "effective wavelength" (defined as , where = 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.
This operator supports local-mode ACA assembly. It does not support hybrid-mode assembly, though.
If local-mode ACA assembly is requested (see AcaOptions::mode), after discretization, the weak form of this operator is stored as
where:
is the wave number
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;
is the sparse matrix whose transpose represents the expansion of the
th component of the basis functions of dualToRange in the single-element test functions mentioned above;
is the sparse matrix representing the expansion of the
th component of the basis functions of domain in the single-element trial functions;
is the sparse matrix whose transpose represents the expansion of the surface divergence of the basis functions of dualToRange in the single-element test functions;
is the sparse matrix representing the expansion of the surface divergence of the basis functions of domain in the single-element trial functions.| BasisFunctionType | Type 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::AssemblyOptions::ACA, Bempp::AssemblyOptions::acaOptions(), Bempp::AssemblyOptions::assemblyMode(), and Bempp::AcaOptions::mode.
1.8.5