Fast Operators
For some theory look here.
Helper Functions
BoundaryIntegralEquations.partial_assemble_parallel!
— Functionpartial_assemble_parallel!(mesh::Mesh3d,k,in_sources,shape_function::Triangular;
fOn=true,gOn=true,n=3,progress=false,depth=1)
Assembles only local contributions in the BEM computations. This is used for singularity extraction when using the Fast Mulitpole Method (FMM) for BEM.
BoundaryIntegralEquations.unroll_interpolations
— Functionunroll_interpolations(interpolations)
Merges the element interpolations from interpolate_elements
into vectors.
BoundaryIntegralEquations.create_coefficient_map
— Functioncreate_coefficient_map(weights,physics_topology,physics_function,n_gauss)
Returns the mapping from nodal values to coefficients for the FMM/H-matrix.
This mapping is represented by a (sparse) matrix $C$ with rows given as
\[ \underbrace{\text{jacobian}(\mathbf{u}_j)w_j\mathbf{T}^{e(j)}(\mathbf{u}_j)\mathbf{L}^{e(j)}}_{j\text{th row of } \mathbf{C}},\]
where each of the $j$ rows corresponds to the Gaussian points from all elements and $e(j)$ is a function that returns the element number that Gaussian point $j$ is located on.
BoundaryIntegralEquations.scale_columns!
— Functionscale_columns!(dipvecs,normals,scalings)
Saves the $i$th column of normals
scaled by the $i$th value of scalings
in dipvecs
.
Used when creating the "dipole-vectors" required for using the FMM.
BoundaryIntegralEquations.setup_fast_operator
— Functionsetup_fast_operator(mesh,zk,n_gauss,nearfield,offset,depth;single_layer=true)
Computes sources
, normals
, C_map
and the nearfield_correction
required when setting up a fast operator.
BoundaryIntegralEquations.create_single_layer_matrix
— Functioncreate_single_layer_matrix(X, Y, k)
Creates a KernelMatrix
representing the single layer potential.
BoundaryIntegralEquations.evaluate_targets
— Methodevaluate_targets(A::FMMFOperator,x,targets)
Evaluates the G-operator with coefficients x
at the targets
positions. The wavenumber used is defined by the FMMGOperaetor
.
BoundaryIntegralEquations.evaluate_targets
— Methodevaluate_targets(A::FMMFOperator,k,x,targets)
Evaluates the F-operator with coefficients x
at the targets
positions. The wavenumber used, k
, is supplied by the user.
BoundaryIntegralEquations.evaluate_targets
— Methodevaluate_targets(A::FMMGOperator,x,targets)
Evaluates the G-operator with coefficients x
at the targets
positions. The wavenumber used is defined by the FMMGOperaetor
.
BoundaryIntegralEquations.evaluate_targets
— Methodevaluate_targets(A::FMMGOperator,k,x,targets)
Evaluates the G-operator with coefficients x
at the targets
positions. The wavenumber used, k
, is supplied by the user.
Fast Multipole Operators
BoundaryIntegralEquations.FMMGOperator
— TypeFMMGOperator(k,tol,targets,sources,C,coefficients,nearfield_correction)
FMMGOperator(mesh,k;tol=1e-6,n_gauss=3,nearfield=true,offset=0.2,depth=1)
A LinearMap
that represents the BEM $\mathbf{G}$ matrix through the FMM. This matrix has $k$th row given by $\mathbf{z}=\mathbf{z}_k$ in the following
\[\begin{aligned} \left(\int_{\Gamma} G(\mathbf{x},\mathbf{z})\mathbf{T}(\mathbf{x}) \mathrm{d}S_\mathbf{x}\right)\mathbf{y} &\approx \left(\sum_{e=1}^{N}\left(\sum_{i=1}^{Q}G(\mathbf{x}^e(\mathbf{u}_i),\mathbf{z})\text{jacobian}(\mathbf{u}_i)w_i\mathbf{T}^e(\mathbf{u}_i)\right)\mathbf{L}^e\right)\mathbf{y} \newline &= \left(\sum_{j=1}^{NQ}G(\mathbf{x}_j,\mathbf{z})\underbrace{\text{jacobian}(\mathbf{u}_j)w_j\mathbf{T}^{e(j)}(\mathbf{u}_j)\mathbf{L}^{e(j)}}_{j\text{th row of } \mathbf{C}}\right)\mathbf{y} \newline &= \begin{bmatrix} G(\mathbf{x}_1,\mathbf{z}) & G(\mathbf{x}_2,\mathbf{z}) & \dots & G(\mathbf{x}_{NQ},\mathbf{z}) \end{bmatrix} \mathbf{C}\mathbf{y}, \end{aligned}\]
where the subscript $j$ refers to an ordering of the collection of Gaussian points from all elements and $e(j)$ is a function that returns the element number that Gaussian point $j$ is located on. Note that $\mathbf{C}$ is the same same for all $k$'s.
The remaining multiplication with the Green's functions is done utilizing the Flatiron Institute Fast Multipole libraries.
BoundaryIntegralEquations.FMMFOperator
— TypeFMMFOperator(k,tol,targets,sources,normals,C,coefficients,dipvecs,nearfield_correction)
FMMFOperator(mesh,k;n_gauss=3,tol=1e-6,nearfield=true,offset=0.2,depth=1,)
A LinearMap
that represents the BEM $\mathbf{H}$ matrix through the FMM. This matrix has $k$th row given by $\mathbf{z}=\mathbf{z}_k$ in the following
\[\begin{aligned} \left(\int_{\Gamma}\frac{\partial G(\mathbf{x}, \mathbf{z})}{\partial\mathbf{n}(\mathbf{x})}\mathbf{T}(\mathbf{x}) \mathrm{d}S_{\mathbf{x}}\right)\mathbf{y} &\approx \left(\sum_{e=1}^{N}\left(\sum_{i=1}^{Q}\frac{\partial G(\mathbf{x}^e(\mathbf{u}_i), \mathbf{z})}{\partial\mathbf{n}(\mathbf{x})}\text{jacobian}(\mathbf{u}_i)w_i\mathbf{T}^e(\mathbf{u}_i)\right)\mathbf{L}^e\right)\mathbf{y}\newline &= \left(\sum_{j=1}^{NQ}\nabla G(\mathbf{x}_j, \mathbf{z})\cdot \mathbf{n}(\mathbf{x}_j)\underbrace{\text{jacobian}(\mathbf{u}_j)w_j\mathbf{T}^{e(j)}(\mathbf{u}_j)\mathbf{L}^{e(j)}}_{j\text{th row of }\mathbf{C}}\right)\mathbf{y}\newline &= \begin{bmatrix} \nabla G(\mathbf{x}_1, \mathbf{z}) \cdot \mathbf{n}(\mathbf{x}_1) & \dots & \nabla G(\mathbf{x}_{NQ}, \mathbf{z})\cdot \mathbf{n}(\mathbf{x}_{NQ}) \end{bmatrix} \mathbf{C}\mathbf{y}, \end{aligned}\]
where $\mathbf{C}$ is coefficient mapping (the same for all $k$).
The remaining multiplication with the Green's functions is performed utilizing the Flatiron Institute Fast Multipole libraries.
H-Matrix Operators
BoundaryIntegralEquations.HGOperator
— TypeHGOperator(k,G,C,nearfield_correction,coefficients)
HGOperator(mesh,k;tol=1e-4,n_gauss=3,nearfield=true,offset=0.2,depth=1)
A LinearMap
that represents the BEM $\mathbf{G}$ matrix through the H-matrix approach. This matrix has $k$th row given by $\mathbf{z}=\mathbf{z}_k$ in the following
\[\begin{aligned} \left(\int_{\Gamma} G(\mathbf{x},\mathbf{z})\mathbf{T}(\mathbf{x}) \mathrm{d}S_\mathbf{x}\right)\mathbf{y} &\approx \left(\sum_{e=1}^{N}\left(\sum_{i=1}^{Q}G(\mathbf{x}^e(\mathbf{u}_i),\mathbf{z})\text{jacobian}(\mathbf{u}_i)w_i\mathbf{T}^e(\mathbf{u}_i)\right)\mathbf{L}^e\right)\mathbf{y} \newline &= \left(\sum_{j=1}^{NQ}G(\mathbf{x}_j,\mathbf{z})\underbrace{\text{jacobian}(\mathbf{u}_j)w_j\mathbf{T}^{e(j)}(\mathbf{u}_j)\mathbf{L}^{e(j)}}_{j\text{th row of } \mathbf{C}}\right)\mathbf{y} \newline &= \begin{bmatrix} G(\mathbf{x}_1,\mathbf{z}) & G(\mathbf{x}_2,\mathbf{z}) & \dots & G(\mathbf{x}_{NQ},\mathbf{z}) \end{bmatrix} \mathbf{C}\mathbf{y}, \end{aligned}\]
where the subscript $j$ refers to an ordering of the collection of Gaussian points from all elements and $e(j)$ is a function that returns the element number that Gaussian point $j$ is located on. Note that $\mathbf{C}$ is the same same for all $k$'s.
The remaining multiplication with the Green's functions is done utilizing the HMatrices.jl
package.
BoundaryIntegralEquations.HFOperator
— TypeHFOperator(k,H,C,nearfield_correction,coefficients)
HFOperator(mesh,k;tol=1e-4,n_gauss=3,nearfield=true,offset=0.2,depth=1)
A LinearMap
that represents the BEM $\mathbf{H}$ matrix through the H-matrix approach. This matrix has $k$th row given by $\mathbf{z}=\mathbf{z}_k$ in the following
\[\begin{aligned} \left(\int_{\Gamma}\frac{\partial G(\mathbf{x}, \mathbf{z})}{\partial\mathbf{n}(\mathbf{x})}\mathbf{T}(\mathbf{x}) \mathrm{d}S_{\mathbf{x}}\right)\mathbf{y} &\approx \left(\sum_{e=1}^{N}\left(\sum_{i=1}^{Q}\frac{\partial G(\mathbf{x}^e(\mathbf{u}_i), \mathbf{z})}{\partial\mathbf{n}(\mathbf{x})}\text{jacobian}(\mathbf{u}_i)w_i\mathbf{T}^e(\mathbf{u}_i)\right)\mathbf{L}^e\right)\mathbf{y}\newline &= \left(\sum_{j=1}^{NQ}\nabla G(\mathbf{x}_j, \mathbf{z})\cdot \mathbf{n}(\mathbf{x}_j)\underbrace{\text{jacobian}(\mathbf{u}_j)w_j\mathbf{T}^{e(j)}(\mathbf{u}_j)\mathbf{L}^{e(j)}}_{j\text{th row of }\mathbf{C}}\right)\mathbf{y}\newline &= \begin{bmatrix} \nabla G(\mathbf{x}_1, \mathbf{z}) \cdot \mathbf{n}(\mathbf{x}_1) & \dots & \nabla G(\mathbf{x}_{NQ}, \mathbf{z})\cdot \mathbf{n}(\mathbf{x}_{NQ}) \end{bmatrix} \mathbf{C}\mathbf{y}, \end{aligned}\]
where $\mathbf{C}$ is coefficient mapping (the same for all $k$). The remaining multiplication with the Green's functions is performed utilizing the HMatrices.jl
package.