FMM2D.FMM2D
— Modulemodule FMM2D
Wrappers for the Flatiron Institute's FMM2D library.
List of wrappers
All N-body codes return output in an FMMVals
structure. See documentation of N-body codes for details.
N-body interactions with the Helmholtz kernel
FMM2D.cfmm2d
— Method vals = cfmm2d(eps,zk,sources;
charges=nothing,dipstr=nothing,targets=nothing,pg=0,pgt=0,nd=1)
This subroutine computes the N-body Laplace interactions in two dimensions where the interaction kernel is given by log(r) and its gradients.
\[ u(z) = \sum_{j=1}^{N} c_{j} * log\(\|z - \epsilon_j\|\) - \frac{v_j}{z - \epsilon_j},\]
where $c_{j}$ are the charge strengths, $v_{j}$ are the dipole strengths, $z_{j} = x_{j}(1) + ix_{j}(2)$ are the complex number corresponding to $x_{j}$. when $z=\epsilon_{j}$, the jth term is dropped from the sum
Input
eps::Float64
precision requestedsources::Array{Float64}
size (2,n) source locations ($x_{j}$)
Optional Keyword Input
charges::Array{ComplexF64}
size (nd,n)dipstr::Array{ComplexF64}
size (nd,n)targets::Array{Float64}
size (2,nt) target locations ($x$)pg::Integer
source eval flag.- Do not compute any quantities at sources if
pg == 0
- Potential ($u$) at sources evaluated if
pg == 1
. - Potential and gradient ($\nabla u$) at sources evaluated if
pg == 2
- Potential, gradient, and Hessian ($\nabla \nabla u$) at sources evaluated if
pg == 3
- Do not compute any quantities at sources if
pgt::Integer
target eval flag.- Do not compute any quantities at targets if
pgt == 0
- Potential at targets evaluated if
pgt == 1
. - Potential and gradient at targets evaluated if
pgt == 2
- Potential, gradient, and Hessian at targets evaluated if
pgt == 3
- Do not compute any quantities at targets if
nd::Integer
number of densities
Note: if all default values are used for optional input, nothing is computed.
Output
vals::FMMVals
with the fields Note that the Hessian is returned as a length 4 vector at each point with the second derivatives in the order: $\partial_{xx}$, $\partial_{yy}$, $\partial_{xy}$, $\partial_{yx}$.
vals.pot::Array{ComplexF64}
size (nd,n) or (n) potential at source locations if requestedvals.grad::Array{ComplexF64}
size (nd,2,n) or (2,n) gradient at source locations if requestedvals.hess::Array{ComplexF64}
size (nd,3,n) or (3,n) Hessian at source locations if requestedvals.pottarg::Array{ComplexF64}
size (nd,nt) or (nt) potential at target locations if requestedvals.gradtarg::Array{ComplexF64}
size (nd,2,nt) or (2,nt) gradient at target locations if requestedvals.hesstarg::Array{ComplexF64}
size (nd,3,nt) or (3,nt) Hessian at target locations if requestedvals.ier
error flag as returned by FMM2D library. A value of 0 indicates a successful call.
Non-zero values may indicate insufficient memory available. See the documentation for the FMM2D library. If not set (nothing
), then FMM2D library was never called.
FMM2D.h2ddir
— Method vals = h2ddir(zk,sources,targets, charges=nothing,dipstr=nothing,dipvecs=nothing,pgt=0,nd=1,thresh=1e-16)
This function computes the N-body Helmholtz interactions in two dimensions where the interaction kernel is given by $H_0^{(1)}(kr)$ and its gradients. This is the $O(N^2)$ direct evaluation code. By convention this code only computes the effect of sources on targets. If the value at sources is also needed, the routine can be called again with targets equal to the source locations.
\[ u(x) = \sum_{j=1}^{N} c_jH_0^{(1)}(k\|x - x_{j}\|) - v_j\mathbf{d}_j\cdot\nabla H_0^{(1)}(k\|x - x_{j}\|),\]
where $c_{j}$ are the charge densities, $v_{j}$ are the dipole strengths, $d_{j}$ are the dipole orientation vectors and $x_{j}$ are the source locations. When $x=x_{j}$, the jth term is dropped from the sum
Input
zk::ComplexF64
Helmholtz parametersources::Array{Float64}
size (2,n) source locations ($x_{j}$)targets::Array{Float64}
size (2,nt) target locations ($x$)
Optional Keyword Input
charges::Array{ComplexF64}
size (nd,n) or (n) charge densities (c_{j})dipstr::Array{ComplexF64}
size (nd,n) or (n) dipole strengths ($v_{j}$)dipvecs::Array{Float64}
size (nd,2,n) or (2,n) dipole orientation vectors ($d_{j}$)pgt::Integer
target eval flag.- Potential at targets evaluated if
pgt == 1
. - Potential and gradient at targets evaluated if
pgt == 2
- Potential, gradient, and Hessian at targets evaluated if
pgt == 3
- Potential at targets evaluated if
nd::Integer
number of densitiesthresh::Float64
threshold for ignoring interactions when $\|x-x_{j}\| \leq thresh$
Note: if all default values are used for optional input, nothing is computed.
Output
vals::FMMVals
with the fields Note that the Hessian is returned as a length 4 vector at each point with the second derivatives in the order: $\partial_{xx}$, $\partial_{yy}$, $\partial_{xy}$, $\partial_{yx}$.
vals.pottarg::Array{ComplexF64}
size (nd,nt) or (nt) potential at target locations if requestedvals.gradtarg::Array{ComplexF64}
size (nd,2,nt) or (2,nt) gradient at target locations if requestedvals.hesstarg::Array{ComplexF64}
size (nd,4,nt) or (4,nt) Hessian at target locations if requested
FMM2D.hfmm2d
— Method vals = hfmm2d(eps,zk,sources;
charges=nothing,dipstr=nothing,dipvecs=nothing,targets=nothing,pg=0,pgt=0,nd=1)
This function computes the N-body Helmholtz interactions in two dimensions where the interaction kernel is given by $H_0^{(1)}(kr)$ and its gradients. This is the $O(N)$ fast multipole code which computes the interactions to the requested precision.
\[ u(x) = \sum_{j=1}^{N} c_jH_0^{(1)}(k\|x - x_{j}\|) - v_j\mathbf{d}_j\cdot\nabla H_0^{(1)}(k\|x - x_{j}\|),\]
where $c_{j}$ are the charge densities, $v_{j}$ are the dipole strengths, $d_{j}$ are the dipole orientation vectors and $x_{j}$ are the source locations. When $x=x_{j}$, the jth term is dropped from the sum
Input
eps::Float64
precision requestedzk::ComplexF64
Helmholtz parametersources::Array{Float64}
size (2,n) source locations ($x_{j}$)
Optional Keyword Input
charges::Array{ComplexF64}
size (nd,n) or (n) charge densities (c_{j})dipstr::Array{ComplexF64}
size (nd,n) or (n) dipole strengths ($v_{j}$)dipvecs::Array{Float64}
size (nd,2,n) or (2,n) dipole orientation vectors ($d_{j}$)targets::Array{Float64}
size (2,nt) target locations ($x$)pg::Integer
source eval flag.- Do not compute any quantities at sources if
pg == 0
- Potential ($u$) at sources evaluated if
pg == 1
. - Potential and gradient ($\nabla u$) at sources evaluated if
pg == 2
- Potential, gradient, and Hessian ($\nabla \nabla u$) at sources evaluated if
pg == 3
- Do not compute any quantities at sources if
pgt::Integer
target eval flag.- Do not compute any quantities at targets if
pgt == 0
- Potential at targets evaluated if
pgt == 1
. - Potential and gradient at targets evaluated if
pgt == 2
- Potential, gradient, and Hessian at targets evaluated if
pgt == 3
- Do not compute any quantities at targets if
nd::Integer
number of densities
Note: if all default values are used for optional input, nothing is computed.
Output
vals::FMMVals
with the fields Note that the Hessian is returned as a length 4 vector at each point with the second derivatives in the order: $\partial_{xx}$, $\partial_{yy}$, $\partial_{xy}$, $\partial_{yx}$.
vals.pot::Array{ComplexF64}
size (nd,n) or (n) potential at source locations if requestedvals.grad::Array{ComplexF64}
size (nd,2,n) or (2,n) gradient at source locations if requestedvals.hess::Array{ComplexF64}
size (nd,4,n) or (4,n) Hessian at source locations if requestedvals.pottarg::Array{ComplexF64}
size (nd,nt) or (nt) potential at target locations if requestedvals.gradtarg::Array{ComplexF64}
size (nd,2,nt) or (2,nt) gradient at target locations if requestedvals.hesstarg::Array{ComplexF64}
size (nd,4,nt) or (4,nt) Hessian at target locations if requestedvals.ier
error flag as returned by FMM2D library. A value of 0 indicates a successful call.
Non-zero values may indicate insufficient memory available. See the documentation for the FMM2D library. If not set (nothing
), then FMM2D library was never called.
FMM2D.lfmm2d
— Method vals = lfmm2d(eps,zk,sources;
charges=nothing,dipstr=nothing,dipvecs=nothing,targets=nothing,pg=0,pgt=0,nd=1)
This subroutine computes the N-body Laplace interactions in two dimensions where the interaction kernel is given by log(r) and its gradients.
\[ u(x) = \sum_{j=1}^{N} c_{j} * log\(\|x-x_{j}\|\) - d_{j}v_{j} \cdot \nabla( log(\|x-x_{j}\|) ),\]
where $c_{j}$ are the charge densities, $d_{j}$ are the dipole strength vectors $v_{j}$ are the dipole orientation vectors, $x_{j}$ are the source locations. when $x=x_{j}$, the jth term is dropped from the sum
Input
eps::Float64
precision requestedsources::Array{Float64}
size (2,n) source locations ($x_{j}$)
Optional Keyword Input
charges::Array{ComplexF64}
size (nd,n)dipstr::Array{ComplexF64}
size (nd,n)dipvecs::Array{Float64}
size (nd,2,n) or (2,n) dipole orientation vectors ($d_{j}$)targets::Array{Float64}
size (2,nt) target locations ($x$)pg::Integer
source eval flag.- Do not compute any quantities at sources if
pg == 0
- Potential ($u$) at sources evaluated if
pg == 1
. - Potential and gradient ($\nabla u$) at sources evaluated if
pg == 2
- Potential, gradient, and Hessian ($\nabla \nabla u$) at sources evaluated if
pg == 3
- Do not compute any quantities at sources if
pgt::Integer
target eval flag.- Do not compute any quantities at targets if
pgt == 0
- Potential at targets evaluated if
pgt == 1
. - Potential and gradient at targets evaluated if
pgt == 2
- Potential, gradient, and Hessian at targets evaluated if
pgt == 3
- Do not compute any quantities at targets if
nd::Integer
number of densities
Note: if all default values are used for optional input, nothing is computed.
Output
vals::FMMVals
with the fields Note that the Hessian is returned as a length 4 vector at each point with the second derivatives in the order: $\partial_{xx}$, $\partial_{yy}$, $\partial_{xy}$, $\partial_{yx}$.
vals.pot::Array{ComplexF64}
size (nd,n) or (n) potential at source locations if requestedvals.grad::Array{ComplexF64}
size (nd,2,n) or (2,n) gradient at source locations if requestedvals.hess::Array{ComplexF64}
size (nd,3,n) or (3,n) Hessian at source locations if requestedvals.pottarg::Array{ComplexF64}
size (nd,nt) or (nt) potential at target locations if requestedvals.gradtarg::Array{ComplexF64}
size (nd,2,nt) or (2,nt) gradient at target locations if requestedvals.hesstarg::Array{ComplexF64}
size (nd,3,nt) or (3,nt) Hessian at target locations if requestedvals.ier
error flag as returned by FMM2D library. A value of 0 indicates a successful call.
Non-zero values may indicate insufficient memory available. See the documentation for the FMM2D library. If not set (nothing
), then FMM2D library was never called.
FMM2D.rfmm2d
— Method vals = rfmm2d(eps,zk,sources;
charges=nothing,dipstr=nothing,dipvecs=nothing,targets=nothing,pg=0,pgt=0,nd=1)
This subroutine computes the N-body Laplace interactions in two dimensions where the interaction kernel is given by log(r) and its gradients.
\[ u(x) = \sum_{j=1}^{N} c_{j} * log\(\|x-x_{j}\|\) - d_{j}v_{j} \cdot \nabla( log(\|x-x_{j}\|) ),\]
where $c_{j}$ are the charge densities, $d_{j}$ are the dipole strength vectors $v_{j}$ are the dipole orientation vectors, $x_{j}$ are the source locations. when $x=x_{j}$, the jth term is dropped from the sum
Input
eps::Float64
precision requestedsources::Array{Float64}
size (2,n) source locations ($x_{j}$)
Optional Keyword Input
charges::Array{Float64}
size (nd,n)dipstr::Array{Float64}
size (nd,n)dipvecs::Array{Float64}
size (nd,2,n) or (2,n) dipole orientation vectors ($d_{j}$)targets::Array{Float64}
size (2,nt) target locations ($x$)pg::Integer
source eval flag.- Do not compute any quantities at sources if
pg == 0
- Potential ($u$) at sources evaluated if
pg == 1
. - Potential and gradient ($\nabla u$) at sources evaluated if
pg == 2
- Potential, gradient, and Hessian ($\nabla \nabla u$) at sources evaluated if
pg == 3
- Do not compute any quantities at sources if
pgt::Integer
target eval flag.- Do not compute any quantities at targets if
pgt == 0
- Potential at targets evaluated if
pgt == 1
. - Potential and gradient at targets evaluated if
pgt == 2
- Potential, gradient, and Hessian at targets evaluated if
pgt == 3
- Do not compute any quantities at targets if
nd::Integer
number of densities
Note: if all default values are used for optional input, nothing is computed.
Output
vals::FMMVals
with the fields Note that the Hessian is returned as a length 4 vector at each point with the second derivatives in the order: $\partial_{xx}$, $\partial_{yy}$, $\partial_{xy}$, $\partial_{yx}$.
vals.pot::Array{ComplexF64}
size (nd,n) or (n) potential at source locations if requestedvals.grad::Array{ComplexF64}
size (nd,2,n) or (2,n) gradient at source locations if requestedvals.hess::Array{ComplexF64}
size (nd,3,n) or (3,n) Hessian at source locations if requestedvals.pottarg::Array{ComplexF64}
size (nd,nt) or (nt) potential at target locations if requestedvals.gradtarg::Array{ComplexF64}
size (nd,2,nt) or (2,nt) gradient at target locations if requestedvals.hesstarg::Array{ComplexF64}
size (nd,3,nt) or (3,nt) Hessian at target locations if requestedvals.ier
error flag as returned by FMM2D library. A value of 0 indicates a successful call.
Non-zero values may indicate insufficient memory available. See the documentation for the FMM2D library. If not set (nothing
), then FMM2D library was never called.
FMM2D.scalarfmm2dinputcheck
— Method anyfail, n, nt, ifcharge, ifdipole = (
scalarfmm2dinputcheck(sources,charges,dipvecs,targets,pg,pgt,nd) )
Input checking for Helmholtz and Laplace routines. Checks the sizes of these arrays and the compatibility of the various flags, nd, and provided arrays. If something is off, a warning is issued and anyfail
is set to true. Non-fatal mistakes result in a warning but anyfail
remains false.
Output:
anyfail
- boolean is true if any checks fail, false otherwise- n - integer, number of sources
- nt - integer, number of targets
- ifcharge - integer, is 1 if there are charges, 0 otherwise
- ifdipole - integer, is 1 if there are dipoles, 0 otherwise