Base

Physical mediums

The types and functions shared between all different types of physical mediums.

MultipleScattering.PhysicalMediumType
PhysicalMedium{Dim,FieldDim}

An abstract type used to represent the physical medium, the dimension of the field, and the number of spatial dimensions. We expect every sub-type of PhysicalMedium{Dim,1} to have a field c which is the complex wave speed.

source
MultipleScattering.ScalarMediumType
ScalarMedium{Dim}

An type used to represent a scalar wave which satisfies a Helmoholtz equations. Dim represents the number of spatial dimensions.

source
MultipleScattering.boundary_dataFunction

A tuples of vectors of the field close to the boundary of the shape. The field is calculated from sim::FrequencySimulation, but the PhysicalMedium inside and outside of the shape are assumed to be given by insidemedium and outsidemedium.

source
MultipleScattering.complex_legendre_arrayMethod

complex_legendre_array(l_max::Int, z::Complex{T})

returns a vector of all associated legendre functions with degree l <= l_max. The degree and order (indices) of the elements of the vector are given by spherical_harmonics_indices(l_max::Int).

The associated legendre polynomials are taken from the package AssociatedLegendrePolynomials.jl.

source
MultipleScattering.diffsbesselMethod
diffsbessel(f::Function,m,x)

Differentiates the spherical bessel function 'f' (for any spherical bessel). The order is 'm' and the argument is 'x'. Note 'x' can be a complex number.

source
MultipleScattering.gaunt_coefficientMethod
gaunt_coefficient(l1,m1,l2,m2,l3,m3)

A version of the Gaunt coefficients which are used to write the product of two spherical harmonics. If Y_{l,m} is a complex spherical harmonic, with the typical phase conventions from quantum mechanics, then:

gaunt_coefficient(l1,m1,l2,m2,l3,m3) = 4*π*im^{l2+l3-l1} Integral[Y_{l1,m1}*conj(Y_{l2,m2})*conj(Y_{l3,m3})]

where the integral is over the solid angle.

The most standard gaunt coefficients G(l1,m1;l2,m2;l3) are related through the identity:

4pi * G(l1,m1;l2,m2;l3) = im^(l1-l2-l3) * (-1)^m2 * gaunt_coefficient(l1,m1,l2,-m2,l3,m1+m2)
source
MultipleScattering.shankelh1Method
shankelh1(m,x)

Returns the spherical hankel function of the first kind. The order is 'm' and the argument is 'x'. Note 'x' can be a complex number.

source
MultipleScattering.spherical_harmonicsMethod

spherical_harmonics(l_max::Int, θ::Complex{T}, φ::Union{T,Complex{T}})

returns a vector of all spherical harmonics with degree l <= l_max with complex arguments. The degree and order (indices) of the elements of the vector are given by spherical_harmonics_indices(l_max::Int).

source
MultipleScattering.spherical_harmonicsMethod
spherical_harmonics(l_max::Int, θ::T, φ::T)

Returns a vector of all the spherial harmonics $Y_{(l,m)}(\theta,\phi)$ for all the degrees l and orders m. The degree and order (indices) of the elements of the vector are given by spherical_harmonics_indices(l_max::Int).

source

Shapes and domains

Shapes are used to define the shape of particles, and to define containers (or configurations) where all particles are placed. It is also used in plotting.

MultipleScattering.ShapeMethod
Shape(shape::Shape; addtodimensions = 0.0, vector_translation = zeros(...), kws...)

Alter the shape depending on the keywords.

source
Base.inMethod
in(vector, shape)::Bool

Returns true if vector is in interior of shape, false otherwise.

source
Base.issubsetMethod
issubset(shape1, shape2)::Bool

Returns true if shape1 is entirely contained within shape2, false otherwise (also works with particles).

source
MultipleScattering.:≅Method
iscongruent(p1::Shape, p2::Shape)::Bool
≅(p1::Shape, p2::Shape)::Bool

True if shapes are the same but in different positions (origins), standard mathematical definition.

source
MultipleScattering.iscongruentMethod
iscongruent(p1::Shape, p2::Shape)::Bool
≅(p1::Shape, p2::Shape)::Bool

True if shapes are the same but in different positions (origins), standard mathematical definition.

source
MultipleScattering.points_in_shapeMethod
points_in_shape(Shape; resolution = 20, xres = resolution, yres = resolution,
         exclude_region = EmptyShape(region), kws...)

returns (x_vec, region_inds) where x_vec is a vector of points that cover a box which bounds Shape, and region_inds is an array of linear indices such that x_vec[region_inds] are points contained Shape. For 3D we use zres instead of yres.

source
MultipleScattering.HalfspaceType
Halfspace(normal::AbstractVector{T} [, origin::AbstractVector{T}=zeros()])

A halfspace defined by all the points $\mathbf x$ that satify $(\mathbf x - \mathbf o) \cdot \mathbf n < 0$ where $\mathbf n$ is the unit normal and $\mathbf o$ is the origin.

source
MultipleScattering.PlateType
Plate(normal::AbstractVector{T}, width::T, [, origin::AbstractVector{T}=zeros()])

A plate defined by all the points $\mathbf x$ that satify $|(\mathbf x - \mathbf o) \cdot \mathbf n| < w /2$ where $\mathbf n$ is the unit normal, $\mathbf o$ is the origin, and $w$ is the width.

source
MultipleScattering.BoxType
Box([origin::AbstractVector{T}=zeros(),] dimensions::AbstractVector{T})

A Box for 2D and 3D with axis aligned sides, defined by dimensions and origin (at the center).

source
MultipleScattering.BoxMethod
Box(points::Vector{v} where v <: AbstractVector)

A Box for any dimension with axis aligned sides, that is a minimal cover for the points.

source
MultipleScattering.TimeOfFlightPlaneWaveToPointType

A shape which contains all particles,$x > x_0$, necessary to simulate a plane-wave scattering from an infinite medium, for a reciever at the focal point, for time $t < D / c$, where $c$ is the sound speed of the background medium, and $D$ is some chosen focal distance. The plane-wave travels towards the positive direction of the $x$ axis.

More precisely, if the focal point is at $(x_f,y_f)$ then the interior of the shape is defined as $(y - y_f)^2 < (D + x_0)^2 - x_f^2 - 2(D + x_0 - x_f)x$ and $x > min(x_0, x_f)$ where $D$ is the focal distance.

source

Particle

Particle types and functions.

MultipleScattering.iscongruentMethod
iscongruent(p1::AbstractParticle, p2::AbstractParticle)::Bool
≅(p1::AbstractParticle, p2::AbstractParticle)::Bool

Returns true if medium and shape of particles are the same, ignoring origin, false otherwise.

source

RegularSource

RegularSource types and functions for all physical mediums.

MultipleScattering.PlaneSourceType
PlaneSource(medium::P, amplitude::SVector, direction::SVector)

Is a struct type which describes a plane-wave source that drives/forces the whole system. It has three fields: a physical medium, an amplitude of the source, and the direction the propagate in direction.

For any given angular frequency ω, the PlaneSource has the value $e^{i ω/c \mathbf v \cdot \mathbf x }$ at the point $\mathbf x$, where $c$ is the medium wavespeed and $\mathbf v$ is the direction.

source
MultipleScattering.RegularSourceType
RegularSource(medium::P, field::Function, coef::Function)

Is a struct type which describes the source field that drives/forces the whole system. It is also known as an incident wave. It has three fields RegularSource.medium, RegularSource.field, and RegularSource.coef.

The source field at the position 'x' and angular frequency 'ω' is given by

x = [1.0,0.0]
ω = 1.0
RegularSource.field(x,ω)

The field RegularSource.coef regularbasisfunction(medium::Acoustic{T,2}, ω::T)

source
MultipleScattering.regular_spherical_sourceMethod
regular_spherical_source(PhysicalMedium,regular_coefficients; amplitude = one(T), position = zeros(T,Dim))

$c_n =$regular_coefficients[n], $x_o=$position, and let $v_n(kx)$ represent the regular spherical basis with wavenumber $k$ at position $x$. The above returns a source where source.field(x,ω) =$\sum_n c_n v_n(kx -k x_0)$

source
MultipleScattering.source_expandMethod
source_expand(RegularSource, centre; basis_order = 4)

Returns a function of (x,ω) which approximates the value of the source at (x,ω). That is, the source is written in terms of a regular basis expansion centred at centre.

source

Impulse

For convert to the time domain.

MultipleScattering.ContinuousImpulseType

See also: DiscreteImpulse, frequency_to_time

ContinuousImpulse{T<:AbstractFloat}

A struct used to represent an analytic impulse function. Has two fields: in_time a function of time t, and in_freq a function of the angular frequency ω. in_freq should be the Fourier transform of in_time, though this is not enforced.

We use the Fourier transform convention: F(ω) = ∫ f(t)exp(imωt) dt, f(t) = (2π)^(-1) * ∫ F(ω)exp(-imωt) dω.

An impluse f(t) is convoluted in time with the field u(t), however we avoid the convolution by working with the fourier transform F(ω) of the impulse f(t), which results in

frequency to time: (2π)^(-1) * ∫ F(ω)U(ω)exp(-imωt) dω

source
MultipleScattering.FreqDiracImpulseMethod

Dirac Delta function of unit area in the frequency domain centred at ω=ω0.

Warning: in frequency space this is a singuarity and so may lead to unexpected behaviour.

source
MultipleScattering.TimeDiracImpulseMethod
TimeDiracImpulse(t0::T)

Dirac Delta function of unit area in the time domain centred at t=t0.

Warning: in the time domain this is a singuarity and so may lead to unexpected behaviour.

source
MultipleScattering.frequency_to_timeMethod

Convert a FrequencySimulationResult into a TimeSimulationResult by using the inverse fourier transform. Assumes only positive frequencies and a real time signal

source
MultipleScattering.frequency_to_timeMethod
frequency_to_time(field_mat::AbstractArray, ω_vec::AbstractVector,
    t_vec::AbstractArray = ω_to_t(ω_vec);
    method = :DFT,
    impulse::ContinuousImpulse = TimeDiracImpulse(zero(T)),
    discrete_impulse::DiscreteImpulse = continuous_to_discrete_impulse(impulse, t_vec, ω_vec))

See also: DiscreteImpulse, ContinuousImpulse

Calculates the time response from the frequency response by approximating an inverse Fourier transform. The time signal is assumed to be real and the frequenices ω_vec are assumed to be positive (can include zero) and sorted. The result is convoluted in time ωith the user specified impulse.

We use the Fourier transform convention: F(ω) = ∫ f(t)exp(imωt) dt f(t) = (2π)^(-1) * ∫ F(ω)exp(-imωt) dt

To easily sample any time, the default is not FFT, but a discrete version of the transform above.

source
MultipleScattering.time_to_frequencyMethod
function time_to_frequency(field_mat::AbstractArray, t_vec::AbstractVector,
    ω_vec::AbstractArray = t_to_ω(t_vec);
    method = :DFT,
    impulse::ContinuousImpulse = TimeDiracImpulse(zero(T)),
    discrete_impulse::DiscreteImpulse = continuous_to_discrete_impulse(impulse, t_vec, ω_vec)
)

The inverse of the function frequencytotime (only an exact inverse when using :DFT integration). We use the Fourier transform convention: F(ω) = ∫ f(t)exp(imω*t) dt

source
MultipleScattering.ω_to_tMethod

returns an array of time from the frequency array ωvec. Uses the same convention for sampling the time as the discrete Fourier transfrom. Assumes ωvec is ordered and non-negative.

source

Simulation

Simulation types and functions.

Base.runMethod
run(sim::FrequencySimulation, x, ω; basis_order=5)

Run the simulation sim for the position x and angular frequency ω.

Position can be an SVector or Vector{SVector} and frequency can be a float or vector of floats.

source
Base.runMethod
run(sim::FrequencySimulation, region::Shape;
    res=20, xres=res, yres=res, basis_order=5)

Run the simulation sim for a grid of positions in region and for angular frequency ω.

Frequency can be a float or vector of floats. The resolution of the grid points is defined by xres and yres.

source
MultipleScattering.basis_coefficientsFunction
basis_coefficients(sim::FrequencySimulation, ω::AbstractFloat; basis_order::Int=5)::Matrix{Complex}

Return coefficients for bases around each particle for a given simulation and angular frequency (ω).

source
MultipleScattering.fieldFunction
field(result::SimulationResult, [i::Integer, j::Integer])

Get field from result, optionally specifying indices.

Returns single value of/matrix of complex SVectors() if vector field, and complex float if scalar field.

source
MultipleScattering.t_matrixFunction
t_matrix(particle, medium, ω, order)

Returns a finite T-matrix, with size depending on order, for a specific particle within a medium with specific physical properties.

source
t_matrix(Particle{2,Acoustic{T,2},Sphere{T,2}}, Acoustic{T,2}, ω, order)

The T-matrix for a 2D circlular acoustic particle in a 2D acoustic medium.

source
t_matrix(Particle{2,Acoustic{T,2},Sphere{T,2}}, Acoustic{T,2}, ω, order)

The T-matrix for a 2D isotropic Helmholtz resonator in a 2D acoustic medium.

source
t_matrix(Particle{2,Acoustic{T,2},Sphere{T,2}}, Acoustic{T,2}, ω, order)

The T-matrix for a 2D non-isotropic Helmholtz resonator in a 2D acoustic medium.

source
t_matrix(CapsuleParticle{2,Acoustic{T,2},Sphere{T,2}}, Acoustic{T,2}, ω, order)

The T-matrix for a 2D circlular capsule particle in an acoustic medium.

source
t_matrix(Particle{3,Acoustic{T,3},Sphere{T,3}}, Acoustic{T,3}, ω, order)

The T-matrix for a 3D spherical acoustic particle in a 3D acoustic medium.

source
MultipleScattering.get_t_matricesFunction
get_t_matrices(PhysicalMedium, Vector{Particles}, ω, basis_order::Integer)

Returns vector of T-matrices from a vector of particles in a specific domain. Can save computation if multiple of the same kind of particle are present in the vector.

source