The HoloPy Scattering Theories
The HoloPy ScatteringTheory classes know how to calculate scattered
fields from detector and scatterer information. Each scattering theory is only
able to work with certain specific scatterers.
There are two broad classes of scattering theories in HoloPy:
lens-free theories which treat the recorded fields as
the magnified image of the fields at the focal plane, and
lens-based theories which use a more detailed
description of the effects of the objective lens. The lens-free theories
usually do not need any additional parameters specified, whereas the
lens theories need the lens’s acceptance angle, which can be specified
as either a fixed number or a Prior object, representing an
unknown value to be determined in an inference calculation.
All scattering theories in HoloPy inherit from the ScatteringTheory class.
Not sure how to choose a scattering theory? See the Which Scattering Theory should I use? section.
ScatteringTheory Methods
HoloPy Scattering theories calculate the scattered fields through one of the following methods.
_raw_fields()Calculates the scattering fields.
_raw_scat_matrs()Calculates the scattering matrices.
_raw_cross_sections()Calculates the cross sections.
_can_handle()Checks if the theory is compatible with a given scatterer.
If a theory is asked for the raw fields, but does not have a _raw_fields
method, the scattering theory attempts to calculate them via the scattering
matrices, as called by _raw_scat_matrs. More than one of these methods may
be implemented for performance reasons.
Be advised that the ScatteringTheory class is under active
development, and these method names may change.
Lens-Free Scattering Theories
DDACan handle every
Scattererobject in HoloPyComputes scattered fields using the discrete dipole approximation, as implemented by ADDA.
Requires the ADDA package to be installed separately, as detailed in the DDA section
Functions in two different ways, as controlled by the
use_indicatorsflag. If theuse_indicatorsflag isTrue, the scatterer is voxelated within HoloPy before passing to ADDA. If the flag isFalse, ADDA’s built-in scatterer geometries are used for things like spheres, cylinders, ellipsoids, etc.
MieCan handle
Sphereobjects,LayeredSphereobjects, orSpheresthrough superposition.Computes scattered fields using Mie theory.
By default, calculates the radial (near-field) component of scattered electric fields, which is nonradiative.
Lens-Based Scattering Theories
AberratedMieLensCan handle
Sphereobjects, orSpheresthrough superposition.Computes scattered fields using Mie theory, but incorporates both diffractive effects of an objective lens and arbitrary-order spherical aberration.
AberratedMieLensandMieLenshave the same computational cost, althoughAberratedMieLensrequires more parameters for fitting.
Which Scattering Theory should I use?
HoloPy chooses a default scattering theory based off the scatterer type,
currently determined by the function
determine_default_theory_for(). If you’re not satisfied with
HoloPy’s default scattering theory selection, you should choose the
scattering theory based on (1) the scatterer that you are modeling,
and (2) whether you want to describe the effect of the lens on the
recorded hologram in detail.
An Individual Sphere
For single spheres, the default is to calculate scattering using Mie
theory, implemented in the class Mie. Mie theory is the exact
solution to Maxwell’s equations for the scattered field from a spherical
particle, originally derived by Gustav Mie and (independently) by Ludvig
Lorenz in the early 1900s.
By default, HoloPy computes the radial components of the scattered fields, which can affect the hologram if the detector is very close to the detector. But if you have a spherical particle that is far (more than a few wavelengths) from the detector, you can save time and avoid numerical issues by telling HoloPy to use the far-field Mie solutions, as in this example:
import holopy as hp
from holopy.scattering import Sphere, Mie, calc_holo
sphere = Sphere(center=(100, 100, 1500), n = 1.59, r = 7.5)
medium_index = 1.33
illum_wavelen = 0.465
illum_polarization = (1, 0)
detector = hp.detector_grid(shape=200, spacing=1.0)
far_field_mie = Mie(compute_escat_radial=False, full_radial_dependence=False)
holo = calc_holo(detector, sphere, medium_index, illum_wavelen,
illum_polarization, theory=far_field_mie)
Multiple Spheres
A scatterer composed of multiple spheres can exhibit multiple scattering
and coupling of the near-fields of neighboring particles. Mie theory
doesn’t include these effects, so Spheres objects are by
default calculated using the Multisphere theory, which
accounts for near- and far-field coupling of the scattered fields
through the SCSMFO package from Daniel Mackowski. This calculation relies on the
exact solution to Maxwell’s equation for the scattering from an
arbitrary arrangement of non-overlapping spheres.
If you want to reduce computation time, or your spheres are widely
separated (such that optical coupling between the spheres is
negligible), you can approximate the scattering from a cluster of
spheres by adding the far fields. To perform such a calculation, which
does not account for multiple scattering, you can specify Mie theory
manually when calling the calc_holo() function, as the following
code snippet shows:
from holopy.core.io import get_example_data_path
from holopy.scattering import Spheres
s1 = Sphere(center=(5, 5, 5), n = 1.59, r = .5)
s2 = Sphere(center=(4, 4, 5), n = 1.59, r = .5)
collection = Spheres([s1, s2])
imagepath = get_example_data_path('image0002.h5')
exp_img = hp.load(imagepath)
holo = calc_holo(exp_img, collection, theory=Mie)
Note that the multisphere theory does not work with collections of multi-layered spheres; in this case HoloPy defaults to using Mie theory with superposition (that is, neglecting multiple scattering).
Non-spherical particles
HoloPy also includes scattering theories that can calculate scattering
from non-spherical particles. For cylindrical or spheroidal particles,
by default HoloPy calculates scattering from cylindrical or spheroidal
particles by using the Tmatrix theory, which uses the T-matrix
code from Michael Mishchenko.
from holopy.scattering.theory import Tmatrix
from holopy.scattering.scatterer import Spheroid
spheroid = Spheroid(n=1.59, r=(1., 2.), center=(4, 4, 5))
theory = Tmatrix()
holo = calc_holo(exp_img, spheroid, theory=theory)
Holopy can also access a discrete dipole approximation (DDA) theory to model arbitrary non-spherical objects. See the Scattering from Arbitrary Structures with DDA tutorial for more details.
Including the effect of the lens
Most of the scattering theories in HoloPy treat the fields on the detector as a (magnified) image of the fields at the focal plane. While these theories usually provide a good description of holograms of particles far above the focus, when the particle is near the focus subtle optical effects can cause deviations between the recorded hologram and theories which do not specifically describe the effects of the lens. To deal with this, HoloPy currently offers two scattering theories which describe the effects of a perfect lens on the recorded hologram. Both of these scattering theories need information about the lens to make predictions, specifically the acceptance angle of the lens. The acceptance angle \(\beta\) is related to the numerical aperture or NA of the lens by \(\beta = \arcsin(\text{NA} / n_f)\), where \(n_f\) is the refractive index of the immersion fluid. For more details on the effect of the lens on the recorded hologram, see [Leahy2020] and [Martin2021].
The Lens theory allows HoloPy to include the effects of a
perfect objective lens with any scattering theory. The Lens theory works
by wrapping a normal scattering theory. For instance, to calculate the
image of a sphere in an objective lens with an acceptance angle of 1.0,
do
from holopy.scattering.theory import Lens, Mie
lens_angle = 1.0
theory = Lens(lens_angle, Mie())
This theory can then be passed to calc_holo() just like any other
scattering theory. However, calculations with the Lens theory
are very slow, orders of magnitude slower than calculations without the
lens.
To get around the slow speed of the Lens theory, HoloPy offers
an additional theory, MieLens, specifically for spherical
particles imaged with a perfect lens. For spherical particles, some
analytical simplifications are possible which greatly speed up the
description of the objective lens – in fact, the MieLens
theory’s implementation is slightly faster than Mie theory’s.
The following code creates a MieLens theory, which can be
passed to calc_holo() just like any other scattering theory:
from holopy.scattering.theory import MieLens
lens_angle = 1.0
theory = MieLens(lens_angle)
In addition, HoloPy supports the calculation of holograms of spherical
particles when the imaging objective lens has spherical aberrations of
arbitrary order. Currently only spherical aberrations are supported, and
only for the image of spherical scatterers. The following code creates
a AberratedMieLens theory with aberrations up to 8th order in
the phase. This theory can be passed to calc_holo() just like any
other scattering theory:
from holopy.scattering.theory import AberratedMieLens
aberration_coefficients = [1.0, 2.0, 3.0]
lens_angle = 1.0
theory = AberratedMieLens(
spherical_aberration=aberration_coefficients,
lens_angle=lens_angle)
Why isn’t my scattering theory here?!
If you don’t see your favorite scattering theory, you can add it to HoloPy! See Adding a new scattering theory for details. If you think your new scattering theory may be useful for other users, please consider submitting a pull request.