API Reference

SphericalScattering.DielectricSphereThinImpedanceLayerMethod
DielectricSphereThinImpedanceLayer(
    radius    = error("missing argument `radius`"),
    thickness = error("missing argument `thickness` of the coating"),
    thinlayer = error("missing argument `thinlayer`"),
    filling   = error("missing argument `filling`")
)

Constructor for the dielectric sphere with a thin impedance layer. For this model, it is assumed that the displacement field is only radial direction in the layer, which requires a small thickness and low conductivity. For details, see for example T. B. Jones, Ed., “Models for layered spherical particles,” in Electromechanics of Particles, Cambridge: Cambridge University Press, 1995, pp. 227–235. doi: 10.1017/CBO9780511574498.012.

source
SphericalScattering.FitzgeraldDipoleMethod
ex = FitzgeraldDipole(;
        embedding   = Medium(ε0, μ0),
        frequency   = error("missing argument `frequency`"),
        amplitude   = 1.0,
        position    = error("missing argument `position`"),
        orientation = SVector{3,typeof(frequency)}(0.0, 0.0, 1.0),
)
source
SphericalScattering.HertzianDipoleMethod
ex = HertzianDipole(
        embedding   = Medium(ε0, μ0),
        frequency   = error("missing argument `frequency`"),
        amplitude   = 1.0,
        position    = error("missing argument `position`"),
        orientation = SVector{3,typeof(frequency)}(0.0, 0.0, 1.0),
)
source
SphericalScattering.LayeredSphereMethod
LayeredSphere( 
    radii   = error("Missing argument `radii`"), 
    filling = error("`missing argument `filling`")
)

Constructor for the layered dielectric sphere.

source
SphericalScattering.SphericalModeTEMethod
ex = SphericalModeTE(;
        embedding   = Medium(ε0, μ0),
        frequency   = error("missing argument `frequency`"),
        amplitude   = 1.0,
        m           = 0,
        n           = 1,
        c           = 1,
        center      = SVector{3,typeof(frequency)}(0.0, 0.0, 0.0),
        orientation = SVector{3,typeof(frequency)}(0.0, 0.0, 1.0),
)
source
SphericalScattering.SphericalModeTMMethod
ex = SphericalModeTM(;
        embedding   = Medium(ε0, μ0),
        frequency   = error("missing argument `frequency`"),
        amplitude   = 1.0,
        m           = 0,
        n           = 1,
        c           = 1,
        center      = SVector{3,typeof(frequency)}(0.0, 0.0, 0.0),
        orientation = SVector{3,typeof(frequency)}(0.0, 0.0, 1.0),
)
source
SphericalScattering.amplitudeMethod
amplitude(sphere, excitation::PlaneWave, quantity::MagneticField, r)

Returns $E₀$, where $E₀$ is the electric field of the incident plane wave. For PEC layers, it returns $0$.

source
SphericalScattering.amplitudeMethod
amplitude(sphere, excitation::PlaneWave, quantity::MagneticField, r)

Returns $E₀$, where $E₀$ is the electric field of the incident plane wave.

source
SphericalScattering.amplitudeMethod
amplitude(sphere, excitation::PlaneWave, quantity::MagneticField, r)

Returns $H₀/ηᵢ$, where $H₀$ is the magnetic field of the incident plane wave. For PEC layers, it returns $0$.

source
SphericalScattering.associatedLegendreMethod
associatedLegendre(n::T, m::T, ϑ::F)

Compute the normalized associated Legendre polynomials according to the definition by Hansen (A1.25) times m divided by sin(ϑ), but including the Condon-Shortley phase. By using the recursion relations (A1.34a) all corner cases are accounted for correctly.

m * bar(P)n|m|(cos(ϑ)) / sin(ϑ)

ϑ ∈ [0, π] assumed.

source
SphericalScattering.convertSpherical2CartesianMethod
convertSpherical2Cartesian(F_sph, point_sph)

Takes a 3D (3 entry) vector F_sph and converts it from its spherical basis to its Cartesian basis representation.

The location of the vector has to be provided in spherical coordinates by point_sph ordered as $(r, ϑ, φ)$.

source
SphericalScattering.derivatieAssociatedLegendreMethod
derivatieAssociatedLegendre(n::T, m::T, ϑ::F) where {T<:Integer,F<:Real}

Compute the derivatives of the normalized associated Legendre polynomials according to the definition by Hansen (A1.25), but including the Condon-Shortley phase. By using the recursion relations (A1.34b) all corner cases are accounted for correctly.

d bar(P)n|m|(cos(ϑ)) / dϑ

ϑ ∈ [0, π] assumed.

source
SphericalScattering.electricRingCurrentMethod
ex = electricRingCurrent(;
        embedding   = Medium(ε0, μ0),
        frequency   = error("missing argument `frequency`"),
        amplitude   = 1.0,
        radius      = error("missing argument `radius`"),
        center      = error("missing argument `center`"),
        orientation = SVector{3,typeof(frequency)}(0.0, 0.0, 1.0),
)
source
SphericalScattering.expansionMethod
expansion(sphere::Sphere, excitation::PlaneWave, quantity::Field, r, plm, cosϑ, sinϑ, n::Int)

Compute functional dependencies of the Mie series for a plane wave travelling in +z-direction with polarization in x-direction.

source
SphericalScattering.fieldMethod
field(excitation::FitzgeraldDipole, point, quantity::FarField; parameter::Parameter=Parameter())

Compute the electric far field radiated by a Fitzgerald dipole at given position and orientation at point.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.fieldMethod
field(excitation::HertzianDipole, point, quantity::FarField; parameter::Parameter=Parameter())

Compute the electric far field radiated by a Hertzian dipole at given position and orientation at point.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.fieldMethod
field(excitation::PlaneWave, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field of a plane wave.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.fieldMethod
field(excitation::PlaneWave, point, quantity::FarField; parameter::Parameter=Parameter())

Throw error since the far-field of a plane wave is not defined.

source
SphericalScattering.fieldMethod
field(excitation::PlaneWave, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the magnetic field of a plane wave.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.fieldMethod
field(excitation::PlaneWave, quantity::Field; parameter::Parameter=Parameter())

Compute the electric field of a plane wave.

source
SphericalScattering.fieldMethod
field(excitation::SphericalMode, quantity::Field; parameter::Parameter=Parameter())

Compute the electric field of a spherical mode.

source
SphericalScattering.fieldMethod
field(excitation::SphericalModeTE, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field of a TE spherical mode.

source
SphericalScattering.fieldMethod
field(excitation::SphericalModeTE, point, quantity::FarField; parameter::Parameter=Parameter())

Compute the electric far-field of a TE spherical mode.

source
SphericalScattering.fieldMethod
field(excitation::SphericalModeTE, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field of a TM spherical mode.

source
SphericalScattering.fieldMethod
field(excitation::SphericalModeTE, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric far-field of a TM spherical mode.

source
SphericalScattering.fieldMethod
field(excitation::Dipole, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field radiated by a Hertzian dipole at given position and orientation at point.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.fieldMethod
field(excitation::Dipole, point, quantity::MagneticField; parameter::Parameter=Parameter())

Compute the magnetic field radiated by a Hertzian dipole at given position and orientation at point.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.fieldMethod
field(excitation::Dipole, quantity::Field; parameter::Parameter=Parameter())

Compute the electric field radiated by a magnetic/electric ring current at some position and with some orientation.

source
SphericalScattering.fieldMethod
field(excitation::ElectricRingCurrent, r, ϑ, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field radiated by an electric ring current placed in origin at point (r, ϑ)

source
SphericalScattering.fieldMethod
field(excitation::ElectricRingCurrent, r, ϑ, quantity::FarField; parameter::Parameter=Parameter())

Compute the (electric) far-field radiated by an electric ring current placed in origin at point (r, ϑ)

source
SphericalScattering.fieldMethod
field(excitation::ElectricRingCurrent, r, ϑ, quantity::MagneticField; parameter::Parameter=Parameter())

Compute the magnetic field radiated by an electric ring current placed in origin at point (r, ϑ)

source
SphericalScattering.fieldMethod
field(excitation::RingCurrent, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field radiated by a magnetic/electric ring current at some position and with some orientation.

source
SphericalScattering.fieldMethod
field(sphere::Sphere, excitation::Excitation, quantity::Field; parameter::Parameter=Parameter())

Compute the total field in the presence of a sphere for a given excitation.

source
SphericalScattering.fieldMethod
field(sphere::Sphere, excitation::PlaneWave, quantity::Field; parameter::Parameter=Parameter())

Descriptive error for the total far-field in the presence of a sphere for an incident plane wave.

source
SphericalScattering.fieldMethod
field(sphere::Sphere, excitation::SphericalMode, quantity::Field; parameter::Parameter=Parameter())

Descriptive error for the total far-field in the presence of a sphere for an incident spherical mode.

source
SphericalScattering.fieldMethod
field(excitation::UniformField, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field of a uniform field.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.fieldMethod
field(excitation::UniformField, point, quantity::ScalarPotential; parameter::Parameter=Parameter())

Compute the scalar potential of a uniform field.

source
SphericalScattering.fieldMethod
field(excitation::UniformField, quantity::Field; parameter::Parameter=Parameter())

Compute the field or potential of a uniform field.

source
SphericalScattering.getFieldTypeMethod
getFieldType(excitation::FitzgeraldDipole, quantity::Field)

Exchange electric and magnetic field for a magnetic current + apply duality relations.

source
SphericalScattering.impedanceMethod
impedance(sp::Sphere, ex::Excitation, r)

Returns the wavenumber at radius r in the sphere sp. If this part is PEC, a zero medium is returned.

source
SphericalScattering.magneticRingCurrentMethod
ex = magneticRingCurrent(;
        embedding   = Medium(ε0, μ0),
        frequency   = error("missing argument `frequency`"),
        amplitude   = 1.0,
        radius      = error("missing argument `radius`"),
        center      = error("missing argument `center`"),
        orientation = SVector{3,typeof(frequency)}(0.0, 0.0, 1.0),
)
source
SphericalScattering.mediumMethod
medium(sp::Sphere, ex::Excitation, r)

Returns the medium of the sphere sp at radius r. If this part is PEC, a zero medium is returned. The argument r may be a vector of positions.

source
SphericalScattering.permeabilityMethod
permeability(sp::Sphere, ex::Excitation, r)

Returns the permeability of the sphere sp at radius r. If this part is PEC, a zero permeability is returned. The argument r may be a vector of positions.

source
SphericalScattering.permittivityMethod
permittivity(sp::Sphere, ex::Excitation, r)

Returns the permittivity of the sphere sp at radius r. If this part is PEC, a zero permittivity is returned. The argument r may be a vector of positions.

source
SphericalScattering.phiCutPointsMethod
phiCutPoints(ϕ; r=1.0, resolution=1)

Convenience function returning points (in Cartesian coordinates) on a phi cut at a distance 'r' and a resolution in degrees.

source
SphericalScattering.planeWaveMethod
ex = planeWave(;
        embedding    = Medium(ε0, μ0),
        frequency    = error("missing argument `frequency`"),
        amplitude    = 1.0,
        direction    = SVector{3,typeof(frequency)}(0.0, 0.0, 1.0),
        polarization = SVector{3,typeof(frequency)}(1.0, 0.0, 0.0),
)
source
SphericalScattering.prefacMethod
prefac(m::T, n::T)

Prefactor for both types of spherical vector wave functions.

Note the factor (-m / abs(m))^m is now part of the associated Legendre polynomials via the Condon-Shortley phase.

source
SphericalScattering.rcsMethod
rcs(sphere::Sphere, excitation::Excitation, point_cart; parameter::Parameter=Parameter())

RCS only defined for plane waves, so far.

source
SphericalScattering.rcsMethod
rcs(sphere::Sphere, excitation::Excitation; parameter::Parameter=Parameter())

RCS only defined for plane waves, so far.

source
SphericalScattering.rcsMethod
rcs(sphere::Sphere, excitation::PlaneWave, point_cart; parameter::Parameter=Parameter())

Compute the bistatic radar cross-section (RCS).

source
SphericalScattering.rcsMethod
rcs(sphere::Sphere, excitation::PlaneWave; parameter::Parameter=Parameter())

Compute the monostatic radar cross-section (RCS): the bistatic RCS solely for the incident direction of the plane wave.

source
SphericalScattering.rotate!Method
rotate!(excitation::Excitation, vectors_list; inverse=false)

Determine rotation matrix and perform rotation for a general excitation.

The points are assumed to be in Cartesian coordinates.

The vectors_list IS modified (overwritten).

source
SphericalScattering.rotateMethod
rotate(excitation::Excitation, vectors_list; inverse=false)

Determine rotation matrix and perform rotation for general excitations.

The points are assumed to be in Cartesian coordinates.

The vectors_list is NOT modified.

source
SphericalScattering.scatterCoeffMethod
scatterCoeff(sp::DielectricSphereThinImpedanceLayer, ex::UniformField)

Compute the expansion coefficients for the thin impedance layer case and a uniform static field excitation.

source
SphericalScattering.scatterCoeffMethod
scatterCoeff(sphere::PECSphere, excitation::PlaneWave, n::Int)

Compute scattering coefficients for a plane wave travelling in +z-direction with polarization in x-direction.

source
SphericalScattering.scatterCoeffMethod
scatterCoeff(sphere::PECSphere, excitation::SphericalModeTE, n::Int, ka)

Compute scattering coefficients for a spherical TE mode travelling towards the origin.

source
SphericalScattering.scatterCoeffMethod
scatterCoeff(sphere::PECSphere, excitation::SphericalModeTM, n::Int, ka)

Compute scattering coefficients for a spherical TM mode travelling towards the origin.

source
SphericalScattering.scatterCoeffMethod
scatterCoeff(sphere::LayeredSpherePEC{LN,LD,LR,LC}, excitation::UniformField{FC,FT,FR}, r) where {LN,LD,LR,LC,FC,FT,FR}

Scatter coefficients according to Sihvola and Lindell. However, the radii of the shells are ordered from inside to outside as depicted in the documentation.

source
SphericalScattering.scatterCoeffMethod
scatterCoeff(sphere::LayeredSphere{LN,LR,LC}, excitation::UniformField{FC,FT,FR}, r) where {LN,LR,LC,FC,FT,FR}

Scatter coefficients according to Sihvola and Lindell. However, the radii of the shells are ordered from inside to outside as depicted in the documentation.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::DielectricSphere, excitation::UniformField, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field scattered by a Dielectric sphere, for an incident uniform field with polarization in given direction.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::DielectricSphere, excitation::UniformField, point, quantity::ScalarPotential; parameter::Parameter=Parameter())

Compute the scalar potential scattered by a dielectric sphere, for an incident uniform field with polarization in given direction.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::DielectricSphereThinImpedanceLayer, excitation::UniformField, point, quantity::DisplacementField; parameter::Parameter=Parameter())

Compute the displacement field D = ε * E.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::DielectricSphereThinImpedanceLayer, excitation::UniformField, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field scattered by a dielectric sphere with a thin coating, where the displacement field in the coating is only in radial direction. We assume an an incident uniform field with polarization in the given direction.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::DielectricSphereThinImpedanceLayer, excitation::UniformField, point, quantity::ScalarPotentialJump; parameter::Parameter=Parameter())

Compute the jump of the scalar potential for a dielectric sphere with a thin coating, where the displacement field in the coating is only in radial direction. We assume an an incident uniform field with polarization in the given direction.

More precisely, we compute the difference Δ = Φi - Φe, where Φi is the potential on the inner side and ϕe the exterior potential.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::DielectricSphereThinImpedanceLayer, excitation::UniformField, point, quantity::ScalarPotential; parameter::Parameter=Parameter())

Compute the scalar potential scattered by a dielectric sphere with a thin coating, where the displacement field in the coating is only in radial direction. We assume an an incident uniform field with polarization in the given direction.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::PECSphere, excitation::ElectricRingCurrent, quantity::FarField; parameter::Parameter=Parameter())

Compute the electric far-field scattered by the PEC sphere, where the ring current is placed along the z-axis.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::PECSphere, excitation::FitzgeraldDipole, quantity::FarField; parameter::Parameter=Parameter())

Compute the electric far-field scattered by a PEC sphere, where the dipole is placed along the z-axis at z0.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::PECSphere, excitation::HertzianDipole, quantity::FarField; parameter::Parameter=Parameter())

Compute the electric far-field scattered by a PEC sphere, where the dipole is placed along the z-axis at z0.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::PECSphere, excitation::MagneticRingCurrent, quantity::FarField; parameter::Parameter=Parameter())

Compute the electric far-field scattered by the PEC sphere, where the ring current is placed along the z-axis.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::PECSphere, excitation::SphericalModeTE, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field scattered by a PEC sphere, when excited by a spherical mode travelling towards the origin.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::PECSphere, excitation::SphericalMode, quantity::Field; parameter::Parameter=Parameter())

Compute the electric field scattered by a dipole at some position and orientation.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::PECSphere, excitation::Dipole, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field scattered by a PEC sphere, where the dipole is placed along the z-axis at z0.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::PECSphere, excitation::Dipole, point, quantity::MagneticField; parameter::Parameter=Parameter())

Compute the magnetic field scattered by a PEC sphere, where the dipole is placed along the z-axis at z0.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::PECSphere, excitation::Dipole, quantity::Field; parameter::Parameter=Parameter())

Compute the field scattered by a PEC sphere excited by a dipole at some position and orientation.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::PECSphere, excitation::RingCurrent, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field scattered by the PEC sphere, where the ring current is placed along the z-axis.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::PECSphere, excitation::RingCurrent, quantity::MagneticField; parameter::Parameter=Parameter())

Compute the magnetic field scattered by the PEC sphere, where the ring current is placed along the z-axis.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
field(excitation::ElectricRingCurrent, quantity::Field; parameter::Parameter=Parameter())

Compute the electric field radiated by an electric ring current at some position and orientation

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::PECSphere, excitation::UniformField, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field scattered by a PEC sphere, for an incident uniform field with polarization in the given direction.

The point and returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::PECSphere, excitation::UniformField, point, quantity::ScalarPotential; parameter::Parameter=Parameter())

Compute the scalar potential scattered by a PEC sphere, for an incident uniform field with polarization in the given direction.

The point and returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::Sphere, excitation::PlaneWave, point, quantity::Field; parameter::Parameter=Parameter())

Compute the electric field scattered by a PEC or dielectric sphere, for an incident plane wave travelling in +z-direction with E-field polarization in x-direction.

The point and the returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::Sphere, excitation::PlaneWave, quantity::Field; parameter::Parameter=Parameter())

Compute the electric field scattered by a PEC sphere, for an incident plane wave.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::Sphere, excitation::PlaneWave, quantity::Field; parameter::Parameter=Parameter())

Compute the electric field scattered by a sphere, for an incident uniform field.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::LayeredSpherePEC, excitation::UniformField, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field scattered by a layered dielectric sphere with PEC core, for an incident uniform field with polarization in the given direction using Sihvola and Lindell, 1988, Transmission line analogy for calculating the effective permittivity of mixtures with spherical multilayer scatterers.

In contrast to Sihvola and Lindell the radii of the shells are ordered from inside to outside as depicted in the documentation.

The point and returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::LayeredSpherePEC, excitation::UniformField{FC,FT,FR}, point, quantity::ScalarPotential; parameter::Parameter=Parameter())

Compute the scalar potential scattered by a layered dielectric sphere with PEC core, for an incident uniform field with polarization in the given direction using Sihvola and Lindell, 1988, Transmission line analogy for calculating the effective permittivity of mixtures with spherical multilayer scatterers

In contrast to Sihvola and Lindell the radii of the shells are ordered from inside to outside as depicted in the documentation.

The point and returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::LayeredSphere, excitation::UniformField, point, quantity::ElectricField; parameter::Parameter=Parameter())

Compute the electric field scattered by a layered dielectric sphere, for an incident uniform field with polarization in the given direction using Sihvola and Lindell, 1988, Transmission line analogy for calculating the effective permittivity of mixtures with spherical multilayer scatterers.

In contrast to Sihvola and Lindell the radii of the shells are ordered from inside to outside as depicted in the documentation.

The point and returned field are in Cartesian coordinates.

source
SphericalScattering.scatteredfieldMethod
scatteredfield(sphere::LayeredSphere, excitation::UniformField, point, quantity::ScalarPotential; parameter::Parameter=Parameter())

Compute the scalar potential scattered by a layered dielectric sphere, for an incident uniform field with polarization in the given direction using Sihvola and Lindell, 1988, Transmission line analogy for calculating the effective permittivity of mixtures with spherical multilayer scatterers.

In contrast to Sihvola and Lindell the radii of the shells are ordered from inside to outside as depicted in the documentation.

The point and returned field are in Cartesian coordinates.

source
SphericalScattering.sphericalGridPointsMethod
patternPoints(;r=1.0, resolution=5)

Convenience function returning points (in Cartesian and spherical coordinates) on a spherical grid at a distance 'r' and a resolution in degrees.

source
SphericalScattering.thetaCutPointsMethod
thetaCutPoints(ϑ; r=1.0, resolution=1)

Convenience function returning points (in Cartesian coordinates) on a theta cut at a distance 'r' and a resolution in degrees.

source
SphericalScattering.translateMethod
translate(points, translation::SVector{3,T})

Translate the points in the direction of the translation vector.

All inputs are assumed to be in Cartesian coordinates.

source