Title: | Coupled-Dipole Approximation for Electromagnetic Scattering by Three- Dimensional Clusters of Sub-Wavelength Particles |
---|---|
Description: | Coupled-dipole simulations for electromagnetic scattering of light by sub-wavelength particles in arbitrary 3-dimensional configurations. Scattering and absorption spectra are simulated by inversion of the interaction matrix, or by an order-of-scattering approximation scheme. High-level functions are provided to simulate spectra with varying angles of incidence, as well as with full angular averaging. |
Authors: | Baptiste Auguie [aut, cre] |
Maintainer: | Baptiste Auguie <[email protected]> |
License: | GPL-3 |
Version: | 2.0.0 |
Built: | 2025-01-07 05:46:33 UTC |
Source: | https://github.com/nano-optics/cda |
Coupled-dipole approximation for electromagnetic scattering by three-dimensional clusters of subwavelength particles
Coupled-dipole simulations for electromagnetic scattering of light by sub-wavelength particles in arbitrary 3-dimensional configurations. Scattering and absorption spectra are simulated by inversion of the interaction matrix, or by an order-of-scattering approximation scheme. High-level functions are provided to simulate spectra with varying angles of incidence, as well as with full angular averaging.
baptiste Auguie
Draine BT. The discrete-dipole approximation and its application to interstellar graphite grains. Astrophysical Journal. 1988. ## Any one of the following references should be used to cite and acknowledge the use of this package.
Circular dichroism:
B. Auguie, J.L. Alonso-Gomez, A. Guerrero-Martinez, L.M. Liz-Marzan. Fingers crossed: circular dichroism with a dimer of plasmonic nanorods. J. Phys. Chem. Lett. 2, (2011)
Linear extinction:
B. Auguie, W.L. Barnes. Diffractive coupling in gold nanoparticle arrays and the effect of disorder. Optics Letters (2009)
Array factor (infinite case):
B. Auguie, W.L. Barnes. Collective resonances in gold nanoparticle arrays. Physical Review Letters (2008)
Bare (intrinsic) polarizability of a dye in vacuum
alpha_bare( wavelength = seq(300, 800), alpha_inf = 9.6e-39, alpha_k = c(5.76e-38), lambda_k = c(526), mu_k = c(10000) )
alpha_bare( wavelength = seq(300, 800), alpha_inf = 9.6e-39, alpha_k = c(5.76e-38), lambda_k = c(526), mu_k = c(10000) )
wavelength |
wavelength in nm |
alpha_inf |
scalar real offset |
alpha_k |
vector of oscillator strengths |
lambda_k |
vector of oscillator wavelengths |
mu_k |
vector of oscillator damping terms |
Sum of lorentz oscillators
data.frame
baptiste Auguie
Other user_level polarisability:
alpha_dye()
,
alpha_ellipsoid()
Principal polarisability components for a dye molecule
alpha_dye(sizes, wavelength, medium, ...)
alpha_dye(sizes, wavelength, medium, ...)
sizes |
matrix of particle sizes (scaling factors for polarisability tensor) |
wavelength |
wavelength in nm |
medium |
refractive index of incident medium |
... |
further parameters passed to the Lorentzian function |
The dye is modelled as a sum of Lorentz oscillators
matrix of polarisability
baptiste Auguie
Other user_level polarisability:
alpha_bare()
,
alpha_ellipsoid()
Principal polarisability components for an ellipsoidal particle
alpha_ellipsoid(sizes, material, medium)
alpha_ellipsoid(sizes, material, medium)
sizes |
matrix of cluster sizes in nm |
material |
data.frame with wavelength and epsilon |
medium |
refractive index of surrounding medium |
This long-wavelength polarisability approximation uses the Kuwata prescription
The Kuwata prescription includes semi-empirical terms of radiative correction and dynamic depolarisation to better match the fully retarded dipolar response in a reasonable range of (subwavelength) sizes and aspect ratios.
matrix of polarisability
baptiste Auguie
Kuwata et al. Resonant light scattering from metal nanoparticles: Practical analysis beyond Rayleigh approximation Appl. Phys. Lett. 83, 22 (2003)
Other user_level polarisability:
alpha_bare()
,
alpha_dye()
polarizability
alpha_kuwata(wavelength, epsilon, V, axis, L, medium = 1.33)
alpha_kuwata(wavelength, epsilon, V, axis, L, medium = 1.33)
wavelength |
wavelength |
epsilon |
permittivity |
V |
volume |
axis |
semi-axis along incident field |
L |
shape factor |
medium |
refractive index |
prescription from Kuwata
polarizability
baptiste Auguie
Kuwata et al. Resonant light scattering from metal nanoparticles: Practical analysis beyond Rayleigh approximation Appl. Phys. Lett. 83, 22 (2003)
Other user_level polarizability:
depolarisation()
C++ calculation of the array factor
array_factor(wavelength, N, pitch)
array_factor(wavelength, N, pitch)
wavelength |
wavelength in nm |
N |
half the number of dipoles along one side |
pitch |
pitch in nm |
Brute-force numerical evaluation of the truncated 2D sum of dipole fields in a finite square array
complex array factor
baptiste Auguie
S <- array_factor(seq(400, 600), 10, 500) str(S)
S <- array_factor(seq(400, 600), 10, 500) str(S)
Square array of particles
cluster_array(N, pitch = 500, a = 50, b = 50, c = b)
cluster_array(N, pitch = 500, a = 50, b = 50, c = b)
N |
number of particles |
pitch |
center-to-center distance |
a |
semi-axis along x |
b |
semi-axis along y |
c |
semi-axis along z |
A cluster describing a 2D square array of identical particles
list of class cluster with fields: positions, sizes, angles
baptiste Auguie
Other user_level cluster:
cluster_ball()
,
cluster_chain()
,
cluster_dimer()
,
cluster_helix()
,
cluster_pumpkin()
,
cluster_shell()
,
cluster_single()
A ball of particles on a cubic lattice
cluster_ball(N, R0 = 15, a = 1, b = 1, c = b)
cluster_ball(N, R0 = 15, a = 1, b = 1, c = b)
N |
number of particles |
R0 |
ball radius |
a |
semi-axis along x |
b |
semi-axis along y |
c |
semi-axis along z |
Identical particles fill a sphere with a cubic lattice
list of class cluster with fields: positions, sizes, angles
baptiste Auguie
Other user_level cluster:
cluster_array()
,
cluster_chain()
,
cluster_dimer()
,
cluster_helix()
,
cluster_pumpkin()
,
cluster_shell()
,
cluster_single()
b = cluster_ball(100)
b = cluster_ball(100)
Linear chain of particles
cluster_chain(N, pitch = 500, a = 50, b = 30, c = b)
cluster_chain(N, pitch = 500, a = 50, b = 30, c = b)
N |
number of particles |
pitch |
center-to-center distance |
a |
semi-axis along x |
b |
semi-axis along y |
c |
semi-axis along z |
A cluster describing a linear chain of identical particles
list of class cluster with fields: positions, sizes, angles
baptiste Auguie
Other user_level cluster:
cluster_array()
,
cluster_ball()
,
cluster_dimer()
,
cluster_helix()
,
cluster_pumpkin()
,
cluster_shell()
,
cluster_single()
A dimer of two particles
cluster_dimer( d = 100, a = 35, b = 12, c = b, dihedral = pi/4, alpha1 = 0, alpha2 = 0 )
cluster_dimer( d = 100, a = 35, b = 12, c = b, dihedral = pi/4, alpha1 = 0, alpha2 = 0 )
d |
center-to-center distance |
a |
semi-axis along x |
b |
semi-axis along y |
c |
semi-axis along z |
dihedral |
dihedral angle |
alpha1 |
angle first rod |
alpha2 |
angle second rod |
A cluster describing two particles, with dimer axis along z
list of class cluster with fields: positions, sizes, angles
baptiste Auguie
Other user_level cluster:
cluster_array()
,
cluster_ball()
,
cluster_chain()
,
cluster_helix()
,
cluster_pumpkin()
,
cluster_shell()
,
cluster_single()
Particles arranged along a helix
cluster_helix( N = 5, a = 10, b = 10, c = 20, R0 = 100, pitch = 200, delta = pi/5, delta0 = 0, right = TRUE, angles = c("helix", "random", "fixed"), seed = 123, ... )
cluster_helix( N = 5, a = 10, b = 10, c = 20, R0 = 100, pitch = 200, delta = pi/5, delta0 = 0, right = TRUE, angles = c("helix", "random", "fixed"), seed = 123, ... )
N |
number of particles |
a |
semi-axis along x |
b |
semi-axis along y |
c |
semi-axis along z |
R0 |
radius of helix |
pitch |
pitch of helix |
delta |
angle between particles |
delta0 |
initial angle |
right |
logical, helicity |
angles |
type of angular orientation |
seed |
random seed for reproducibility |
... |
extra arguments (ignored) |
Cluster describing a helical assembly of particles
list of class cluster with fields: positions, sizes, angles
baptiste Auguie
Other user_level cluster:
cluster_array()
,
cluster_ball()
,
cluster_chain()
,
cluster_dimer()
,
cluster_pumpkin()
,
cluster_shell()
,
cluster_single()
Sparse shell of nanoparticles around a spherical core
cluster_pumpkin( N = 50, R0 = 30, cone = 2 * pi, d = 1, a = 1, b = 1, c = 1, tilt = 0, position = c("fibonacci", "hc", "random", "landings"), exclusion = 0.7, seed = 123, ... )
cluster_pumpkin( N = 50, R0 = 30, cone = 2 * pi, d = 1, a = 1, b = 1, c = 1, tilt = 0, position = c("fibonacci", "hc", "random", "landings"), exclusion = 0.7, seed = 123, ... )
N |
number of particles |
R0 |
radius of core |
cone |
type of angular orientation |
d |
distance from core |
a |
semi-axis along x |
b |
semi-axis along y |
c |
semi-axis along z |
tilt |
type of angular orientation |
position |
type of random coverage |
exclusion |
minimum exclusion distance for 'hc' positions |
seed |
random seed for reproducibility |
... |
extra arguments (ignored) |
A cluster describing a discrete shell of nanoparticles in a spherical geometry
list of class cluster with fields: positions, sizes, angles
baptiste Auguie
Other user_level cluster:
cluster_array()
,
cluster_ball()
,
cluster_chain()
,
cluster_dimer()
,
cluster_helix()
,
cluster_shell()
,
cluster_single()
Sparse shell of nanoparticles around a spherical core
cluster_shell( N = 50, R0 = 30, d = 1, a = 1, b = 1, c = 1, orientation = c("radial", "flat", "random"), position = c("fibonacci", "hc", "random", "landings"), exclusion = 5 * N^(-1/4), seed = 123, ... )
cluster_shell( N = 50, R0 = 30, d = 1, a = 1, b = 1, c = 1, orientation = c("radial", "flat", "random"), position = c("fibonacci", "hc", "random", "landings"), exclusion = 5 * N^(-1/4), seed = 123, ... )
N |
number of particles |
R0 |
radius of core |
d |
distance from core |
a |
semi-axis along x |
b |
semi-axis along y |
c |
semi-axis along z |
orientation |
type of angular orientation |
position |
type of random coverage |
exclusion |
minimum exclusion distance for 'hc' positions |
seed |
random seed for reproducibility |
... |
extra arguments (ignored) |
A cluster describing a discrete shell of nanoparticles in a spherical geometry
list of class cluster with fields: positions, sizes, angles
baptiste Auguie
Other user_level cluster:
cluster_array()
,
cluster_ball()
,
cluster_chain()
,
cluster_dimer()
,
cluster_helix()
,
cluster_pumpkin()
,
cluster_single()
Trivial cluster
cluster_single(a, b = a, c = b, phi = 0, theta = 0, psi = 0)
cluster_single(a, b = a, c = b, phi = 0, theta = 0, psi = 0)
a |
semi-axis along x |
b |
semi-axis along y |
c |
semi-axis along z |
phi |
first Euler angle |
theta |
second Euler angle |
psi |
third Euler angle |
A single particle cluster
list of class cluster with fields: positions, sizes, angles
baptiste Auguie
Other user_level cluster:
cluster_array()
,
cluster_ball()
,
cluster_chain()
,
cluster_dimer()
,
cluster_helix()
,
cluster_pumpkin()
,
cluster_shell()
cl = cluster_single(10)
cl = cluster_single(10)
Depolarisation factor for an ellipsoid
depolarisation(x1, x2 = x1, x3 = x2)
depolarisation(x1, x2 = x1, x3 = x2)
x1 |
semi-axis in nm |
x2 |
semi-axis in nm |
x3 |
semi-axis in nm |
calculates the 3 depolarisation factors for a general ellipsoid
shape factor along x1
baptiste Auguie
Other user_level polarizability:
alpha_kuwata()
dye_coverage
dye_coverage(rho, R)
dye_coverage(rho, R)
rho |
surface demsity |
R |
radius |
baptiste Auguie
Other user_level cda utility:
equal_angles()
,
equal_sizes()
,
spheroid_ar()
Utility function to create clusters
equal_angles(phi, theta, gamma, N)
equal_angles(phi, theta, gamma, N)
phi |
Euler angle |
theta |
Euler angle |
gamma |
Euler angle |
N |
number of particles |
Identical particles
3xN matrix
baptiste Auguie
Other user_level cda utility:
dye_coverage()
,
equal_sizes()
,
spheroid_ar()
Utility function to create clusters
equal_sizes(a, b, c, N)
equal_sizes(a, b, c, N)
a |
semi-axis along x |
b |
semi-axis along y |
c |
semi-axis along z |
N |
number of particles |
Identical particles
3xN matrix
baptiste Auguie
Other user_level cda utility:
dye_coverage()
,
equal_angles()
,
spheroid_ar()
Exact calculation of the array factor using code from Javier Garcia de Abajo (part of the pxtal program for multiple scattering calculations in infinite layered 2D arrays)
G0
G0
A data frame with 1000 rows and 3 variables:
normalised wavelength lambda/pitch
in-plane component of the wavevector (0, since normal incidence)
complex value of the array factor
Javier Garcia de Abajo
Exact calculation of the array factor using code from Javier Garcia de Abajo (part of the pxtal program for multiple scattering calculations in infinite layered 2D arrays)
gfun
gfun
A list of two interpolation functions:
real part of G0
imaginary part of G0
Javier Garcia de Abajo
Positions along a helix
helix( R0 = 500, pitch = 600, N = 5, delta = pi/8, delta0 = pi/2, n.smooth = 100 * N, right = TRUE )
helix( R0 = 500, pitch = 600, N = 5, delta = pi/8, delta0 = pi/2, n.smooth = 100 * N, right = TRUE )
R0 |
radius of helix |
pitch |
pitch of helix |
N |
number of particles |
delta |
angle between particles |
delta0 |
initial angle |
n.smooth |
number of points for a finer helix (useful for display) |
right |
logical, helicity |
3D points following a helix
list of positions and angles
baptiste Auguie
Quadrature points on a sphere
quadrature_sphere( Nq = 30, quadrature = c("qmc", "gl", "cheap", "random", "grid", "grid2", "fibonacci", "fibonacci2"), init = TRUE )
quadrature_sphere( Nq = 30, quadrature = c("qmc", "gl", "cheap", "random", "grid", "grid2", "fibonacci", "fibonacci2"), init = TRUE )
Nq |
number of integration points |
quadrature |
quadrature method, using either Gauss Legendre quadrature (default), Quasi Monte Carlo, regular grid, or "cheap" (3 axes) |
init |
(qmc method only) logical: restart, or continue from previous call |
Numerical integration points for angular averaging
baptiste Auguie
creates an rgl ellipsoid
rgl.ellipsoid( x = 0, y = 0, z = 0, a = 1, b = 1, c = 1, phi = 0, theta = 0, psi = 0, subdivide = 3, smooth = TRUE, ... )
rgl.ellipsoid( x = 0, y = 0, z = 0, a = 1, b = 1, c = 1, phi = 0, theta = 0, psi = 0, subdivide = 3, smooth = TRUE, ... )
x |
x |
y |
y |
z |
z |
a |
axis |
b |
axis |
c |
axis |
phi |
phi |
theta |
theta |
psi |
psi |
subdivide |
subdivision |
smooth |
smoothing |
... |
additional params |
deforms, rotate, and translate a sphere
an rgl mesh
baptiste Auguie
Other user_level rgl:
rgl.ellipsoids()
## Not run: require(rgl) ; ee <- rgl.ellipsoid() shapelist3d(ee) ## End(Not run)
## Not run: require(rgl) ; ee <- rgl.ellipsoid() shapelist3d(ee) ## End(Not run)
Create a list of rgl ellipsoids oriented in space
rgl.ellipsoids(positions, sizes, angles, colour = "red", ...)
rgl.ellipsoids(positions, sizes, angles, colour = "red", ...)
positions |
matrix of positions |
sizes |
matrix of axis lengths |
angles |
matrix of Euler angles |
colour |
colour of each ellipsoid |
... |
additional params |
each ellipsoid is specified by its position, dimensions, and Euler angles
rgl mesh
baptiste Auguie
Other user_level rgl:
rgl.ellipsoid()
cl <- helix(0.5, 1, 36, delta=pi/6, n.smooth=1e3) sizes <- equal_sizes(0.04,0.02,0.02,NROW(cl$positions)) ## Not run: require(rgl) ; rgl.ellipsoids(cl$positions, sizes, cl$angles, col="gold")
cl <- helix(0.5, 1, 36, delta=pi/6, n.smooth=1e3) sizes <- equal_sizes(0.04,0.02,0.02,NROW(cl$positions)) ## Not run: require(rgl) ; rgl.ellipsoids(cl$positions, sizes, cl$angles, col="gold")
Random sample
Random sample with minimum exlusion zone ("hard-core process")
Random sample with minimum exlusion zone enforced
Fibonacci coverage of a sphere
sample_random(N) sample_hc(N, exclusion = 0.1, maxiter = 200L, k = 30L) sample_landings(N, exclusion = 0.1) sample_fibonacci(N = 301)
sample_random(N) sample_hc(N, exclusion = 0.1, maxiter = 200L, k = 30L) sample_landings(N, exclusion = 0.1) sample_fibonacci(N = 301)
N |
number of points |
exclusion |
minimum distance allowed between points |
maxiter |
maximum number of iterations |
k |
number of extra new points to try at each iteration |
Produces a set of points that covers rather uniformly the unit sphere with N points with a spiral-like pattern based on a Fibonacci sequence
3xN matrix
3xN matrix
3xN matrix
sample_random
: random sample
sample_hc
: random sample with exclusion zone
sample_landings
: random sample with exclusion zone
sample_fibonacci
:
baptiste Auguie
sample_random(10) sample_hc(10) sample_landings(10)
sample_random(10) sample_hc(10) sample_landings(10)
dispersion spectrum
spectrum_dispersion( cluster, material, medium = 1.33, Incidence = 0, Axes = "z", polarisation = c("linear", "circular"), method = c("solve", "cg", "oos"), Nsca = 50, maxiter = 30, tol = 1e-04, progress = FALSE )
spectrum_dispersion( cluster, material, medium = 1.33, Incidence = 0, Axes = "z", polarisation = c("linear", "circular"), method = c("solve", "cg", "oos"), Nsca = 50, maxiter = 30, tol = 1e-04, progress = FALSE )
cluster |
list describing a cluster |
material |
list |
medium |
medium refractive index |
Incidence |
angular directions of incident field |
Axes |
incident field rotation axis |
polarisation |
linear or circular polarisation |
method |
linear system (solve), conjugate-gradient (cg), or order-of-scattering (oos) |
Nsca |
number of quadrature points in calculation of csca |
maxiter |
integer termination of iterative solver |
tol |
double, tolerance of iterative solver |
progress |
logical, display progress bar |
dispersion spectrum
data.frame
The incident wavevector is along the z direction.
baptiste Auguie
Orientation-averaged spectrum
spectrum_oa( cluster, material, medium = 1.33, quadrature = c("gl", "qmc", "random", "cheap"), Nq = 100, iterative = FALSE, precision = 0.001, Qmax = 10000, dN = Nq, method = c("solve", "cg", "oos"), Nsca = 50, maxiter = 30, tol = 1e-04, progress = FALSE, verbose = TRUE )
spectrum_oa( cluster, material, medium = 1.33, quadrature = c("gl", "qmc", "random", "cheap"), Nq = 100, iterative = FALSE, precision = 0.001, Qmax = 10000, dN = Nq, method = c("solve", "cg", "oos"), Nsca = 50, maxiter = 30, tol = 1e-04, progress = FALSE, verbose = TRUE )
cluster |
cluster (list) |
material |
material |
medium |
refractive index medium |
quadrature |
quadrature method, using either Gauss Legendre quadrature (default), Quasi Monte Carlo, random, or "cheap" (3 Axes) |
Nq |
number of integration points |
iterative |
logical, increase N until convergence (QMC only) |
precision |
relative diff between two runs (QMC only) |
Qmax |
maximum N if convergence not attained (QMC only) |
dN |
iterative increase in N (QMC only) |
method |
linear system (solve), conjugate-gradient (cg), or order-of-scattering (oos) |
Nsca |
quadrature points for scattering cross-section |
maxiter |
integer termination of iterative solver |
tol |
double, tolerance of iterative solver |
progress |
print progress lines |
verbose |
display messages |
OA spectrum
baptiste Auguie
Y. Okada, Efficient numerical orientation quadrature of light scattering properties with a quasi-Monte-Carlo method, Journal of Quantitative Spectroscopy and Radiative Transfer, Volume 109, Issue 9, June 2008, Pages 1719-1742.
Spheroid described by effective radius and aspect ratio
spheroid_ar(rv, h, type = c("prolate", "oblate"))
spheroid_ar(rv, h, type = c("prolate", "oblate"))
rv |
equivolume sphere radius |
h |
aspect ratio |
type |
class of spheroid |
Describe a spheroid by the aspect ratio and effective radius of an equi-volume sphere V = 4/3 pi rv^3 = 4/3 pi a^2 c c = h * a
baptiste Auguie
Other user_level cda utility:
dye_coverage()
,
equal_angles()
,
equal_sizes()
Visualise a cluster of particles
visualise(x, type, outfile = NULL, ...)
visualise(x, type, outfile = NULL, ...)
x |
cluster |
type |
type of visualisation (rgl or povray output) |
outfile |
optional output file for the results |
... |
additional arguments passed to the visualise method |
Helper function for rapid visualisation of cluster geometries.
baptiste Auguie