Base

    Defining the material

    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.

    MultipleScattering.AcousticType
    Acoustic{T<:AbstractFloat,Dim}(ρ::T, c::Complex{T})
    Acoustic(ρ::T, c::Union{T,Complex{AbstractFloat}}, Dim::Integer)

    Physical properties for a homogenous isotropic acoustic medium with wavespeed (c) and density (ρ)

    Simulations in this medium produce scalar (1D) fields in Dim dimensions.

    ParticleCorrelations.SpecieType

    Specie

    Represents a set of particles which are all the same. The type of particle is given by Specie.particle and the volume fraction this specie occupies is given by Specie.volume_fraction.

    We can use Specie.numberofparticles to specify the number of particles, otherwise for an infinite Specie.numberofparticles = Inf.

    The minimum distance between any two particles will equal outer_radius(Specie) * Specie.separation_ratio.

    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.

    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)

    Currently the physical medium can be Acoustics. Both Acoustic, PlaneSource, and RegularSource are all imported from the package MultipleScattering. In the future these will be moved to a new package WaveScatteringBase. Much of the code is dispatched based on the underlying symmetries of the problem

    The symmetry of the material and source

    The symmetry shared between the material shape and source are used to specialise the form of the wavemode, see Theoretical background.

    MultipleScattering.WithoutSymmetryType
    WithoutSymmetry

    A type used to describe materials and incident waves which share no common symmetry. This will lead to the most general regular wave expansion for the eignvectors.

    MultipleScattering.RadialSymmetryType
    RadialSymmetry

    A type used to describe materials and incident waves in that are both radially symmetric. That is, the material is a sphere, and the incident wave is radially symmetric around the center of this spherical material.

    Types of waves

    There are two main types used.

    EffectiveWaves.EffectivePlaneWaveModeType
    EffectivePlaneWaveMode{T<:AbstractFloat,Dim} <: AbstractWaveMode

    Is a struct that represents a mode of the average scattering coefficients for plane wave symmetry. That is, when using an incident plane wave and a plate or halfspace for the material geometry.

    This wave mode has frequency $\omega$ and has has the value $A e^{i k \mathbf v \cdot \mathbf x }$ at the point $\mathbf x$, where $A$ are the eigenvectors, $\mathbf v$ is the direction and $k$ is the effective wavenumber. We represent these quantities as

    • $\omega =$EffectivePlaneWaveMode.ω
    • $k =$EffectivePlaneWaveMode.wavenumber
    • $v =$EffectivePlaneWaveMode.direction
    • $A =$EffectivePlaneWaveMode.eigenvectors
    • $M =$EffectivePlaneWaveMode.basis_order

    To recover the average scattering coefficients $F$ from a particle use:

    $F_{ns} = i^m A_{n s} e^{-i n θ}$

    where $n$ is multi-index for 3 dimensions and $s$ ranges from 1 to the number of species.

    In 2D, inconvieniently, an extra factor of $e^{i k \mathbf v \cdot \mathbf x}$ needs to be multiplied on the right of the above equation where θ is calculated from transmission_angle.

    EffectiveWaves.EffectiveRegularWaveModeType
    EffectiveRegularWaveMode{T,P<:PhysicalMedium,S<:AbstractSymmetry} <: AbstractRegularWaveMode

    Is a struct that represents a mode of the average scattering coefficients $F$ in the form

    $F_{n s} (\mathbf r) = \sum_{n_1 p} A_{p n n_1 s} \mathrm v_{n_1}(k \mathbf r)$

    where $\mathbf r$ is a position vector in space, $k$ the wavenumber, $p$ is used to count the number of eigenvectors, $s$ the number of species, $n$ and $n_1$ are multi-indices, and $\mathrm v_{n_1}$ is a regular spherical wave of order $n_1$ which in 3D is given by

    $\mathrm v_{n_1}(\mathbf r) = \mathrm j_{\ell_1} (k r) \mathrm Y_{\ell_1 m_1}(\hat{\mathbf r})$

    and $\mathrm Y_{\ell_1 m_1}$ is a spherical harmonic, $n_1 = (\ell_1,m_1)$ where we always have that $\ell_1 \geq 0$ and $\ell_1 \geq |m_1|$ according to the conventions of spherical harmonics.

    To translate the mathematics to Julia code we use

    • $\omega =$EffectiveRegularWaveMode.ω
    • $k =$EffectiveRegularWaveMode.wavenumber
    • $A =$EffectiveRegularWaveMode.eigenvectors
    • $n =$EffectiveRegularWaveMode.basis_order
    • $n_1 =$EffectiveRegularWaveMode.basis_field_order

    Effective wavenumbers and wavemodes

    EffectiveWaves.wavenumbersFunction
    wavenumbers(ω, medium, specie; kws...)

    Returns all the possible effective wavenumbers with positive imaginary part when using a single type of particle called specie.

    wavenumbers(ω, medium::PhysicalMedium, micro::Microstructure; kws...)

    Returns all the possible effective wavenumbers with positive imaginary part. This function requires significantly numerical optimisation and so can be slow.

    EffectiveWaves.wavenumber_low_volumefractionFunction
    wavenumber_low_volumefraction(ω::T, medium::Acoustic{T,Dim}, micro::Microstructure{Dim}; basis_order::Int = 2)

    Explicit formula for one effective wavenumber based on a low particle volume fraction expansion.

    Effective wavemodes and eignvectors

    EffectiveWaves.WaveModeFunction
    WaveMode(ω::T, wavenumber::Complex{T}, eigenvectors::Array{Complex{T}}, ::SetupSymmetry; kws...)

    Returns a concrete subtype of AbstractWaveMode depending on the SetupSymmetry. The returned type should have all the necessary fields to calculate scattered waves (currently not true for EffectivePlanarWaves).

    Reflection coefficients

    EffectiveWaves.reflection_coefficientFunction
    reflection_coefficient(ω::T, wave_eff::EffectivePlaneWaveMode, psource::PlaneSource{T,2,1,Acoustic}, material::Material{Halfspace}; [x::T = zero(T)])

    The reflection coefficient in 2D for acoustics for just one EffectivePlaneWaveMode"

    Calculates the reflection coefficient from each wave in waves and then sums the results.

    reflection_coefficient(PlaneSource, Acoustic[, Halfspace = Halfspace(-psource.direction)])

    calculates the reflection coefficient from a homogenious halfspace (assumed to direct incidence if not given), which is also the low frequency reflection from a particulate material when using the effective_medium.

    reflection_coefficient(ω::T, m_wave::MatchPlaneWaveMode,
        source::PlaneSource, material::Material{Halfspace{T,2}};

    Calculate the reflection coefficient from a matched wave. This requires using both the discrete part and effective wavemode of the matched wave.

    Transmission

    EffectiveWaves.transmission_directionFunction
    transmission_direction(k_eff::Complex, incident_wavevector::AbstractArray, surface_normal::AbstractArray)

    Gives the effective direction direction where sum(direction .^ 2) = 1.0 and the components of k_eff .* direction and incident_wavevector which are orthogonal to the surface are the same. Note surface_normal is the outward pointing normal and the incident medium's wavenumber k = sqrt(sum(incident_wavevector .^2)). Below we deduce the result.

    Assume that dot(v,w) = conj(v[i])w[i] and surface_normal[i]*surface_normal[i] = 1. Let wnp be orthogonal to surface_normal

    wnp =  incident_wavevector -  dot(surface_normal,incident_wavevector) .* surface_normal

    Then we know that the effective transmission wavevector has to equal

    k_eff_vec = wnp + α .* surface_normal

    where α is determined so that

    sum(x^2 for x in wnp) + α^2  = sum(x^2 for x in k_eff_vec) = k_eff^2
    EffectiveWaves.transmission_angleMethod
    transmission_angle(wavevector, surface_normal)

    Calculates the transmission angle for a wave propagting inside a material which has the outward normal surface_normal. The wave propagates in the direction of wavevector with a wavenumber given by k = sqrt(sum(wavevector .^2)), where k should have a positive imaginary part.

    EffectiveWaves.transmission_angleMethod
    transmission_angle(wavevector::StaticVector{3}, surface_normal::StaticVector{3})

    Some care is needed when wavevector is a complex vector in 3 dimensions.

    The maths

    Let's define the normal projection vn:

    n = - surface_normal
    vn = dot(n,wavevector) .* n

    remembering that in Julia, and mathematics, dot(v,w) = sum(conj(v[i])*w[i] for i in eachindex(v)).

    We now define the orthogonal component vo and its direction no:

    vo = wavevector - vn
    no = vo / sqrt(sum(vo.^2))

    Note this assumes that vo is a real vector, which it should be due to symmetry.

    To determine the transmission angle θT we calculate the wavenumber

    k = ± sqrt(sum(wavevector .^2))

    where we choose the sign in ± such that imag(k) >= 0.

    From the above we can write the wavevector in the form

    wavevector = k .* (n .* cos(θ) + no .* sin(θ))

    which implies that

    k * cos(θ) = dot(conj(n), wavevector)
    k * sin(θ) = dot(conj(no),wavevector) = dot(conj(no),vo) = sqrt(sum(vo.^2))

    From these two we can robustly calculate θT as

    θT = atan(sqrt(sum(vo.^2)),dot(n,wavevector))

    Statistics

    All types

    EffectiveWaves.MicrostructureMethod
    Microstructure(medium::PhysicalMedium, sps::Vector{Specie})

    When no pair-correlation is specified for the species, the microstructure will use the default that assumes that particles can not overlap, but, otherwise, their positions are uncorrelated. This is often called "Hole Correction"

    EffectiveWaves.ParticulateMicrostructureType
    ParticulateMicrostructure

    Represents a microstructure filled with multiply species of particles. ParticulateMicrostructure.paircorrelations specifies the pair correlation between each of the species. That is, how the particles are distributed on average.

    EffectiveWaves.EffectivePlaneWaveModeMethod
    EffectivePlaneWaveMode(ω, wavenumber::Complex, direction::Array{Complex}, eigenvectors::Array{Complex})

    A convienient way to form the type EffectivePlaneWaveMode.

    EffectiveWaves.ScatteringCoefficientsFieldMethod
    ScatteringCoefficientsField(ω::T, medium::P, material::Material, coefficient_field::Function; symmetry::AbstractSymmetry{Dim} = WithoutSymmetry{Dim})

    A type to hold the results of methods which produce a function of the scattering field