holopy.scattering package
Scattering calculations
The scattering package provides objects and methods to define scatterer geometries, and theories to compute scattering from specified geometries. Scattering depends on holopy.core, and certain scattering theories may require external scattering codes.
The HoloPy scattering module is used to:
Describe geometry as a
scattererobjectDefine the result you want as a xarray.DataArray xarray.DataArray
Calculate scattering quantities with an
theoryappropriate for yourscatterer
Subpackages
- holopy.scattering.scatterer package
- Submodules
- holopy.scattering.scatterer.bisphere module
- holopy.scattering.scatterer.capsule module
- holopy.scattering.scatterer.composite module
- holopy.scattering.scatterer.csg module
- holopy.scattering.scatterer.cylinder module
- holopy.scattering.scatterer.ellipsoid module
- holopy.scattering.scatterer.janus module
- holopy.scattering.scatterer.scatterer module
- holopy.scattering.scatterer.sphere module
- holopy.scattering.scatterer.spherecluster module
- holopy.scattering.scatterer.spheroid module
- holopy.scattering.theory package
- Subpackages
- Submodules
- holopy.scattering.theory.dda module
- holopy.scattering.theory.lens module
- holopy.scattering.theory.mie module
- holopy.scattering.theory.mielens module
- holopy.scattering.theory.mielensfunctions module
- holopy.scattering.theory.multisphere module
- holopy.scattering.theory.scatteringtheory module
- holopy.scattering.theory.tmatrix module
Submodules
holopy.scattering.errors module
Exceptions used in scatterpy module. These are separated out from the other exceptions in other parts of HoloPy to keep things modular.
- exception OverlapWarning(scatterer, overlaps)
Bases:
UserWarning
holopy.scattering.imageformation module
- class ImageFormation(scattering_theory)
Bases:
HoloPyObjectCalculates fields, holograms, intensities, etc.
- calculate_cross_sections(scatterer, medium_wavevec, medium_index, illum_polarization)
- calculate_scattered_field(scatterer, schema)
- Parameters:
scatterer (
scattererobject) – (possibly composite) scatterer for which to compute scattering- Returns:
e_field – scattered electric field
- Return type:
VectorGrid
- calculate_scattering_matrix(scatterer, schema)
Compute scattering matrices for scatterer
- Parameters:
scatterer (
holopy.scattering.scattererobject) – (possibly composite) scatterer for which to compute scattering- Returns:
scat_matr – Scattering matrices at specified positions
- Return type:
Marray
- get_wavevec_from(schema)
- select_scatterer_by_illumination(scatterer, illum)
holopy.scattering.interface module
Base class for scattering theories. Implements python-based calc_intensity and calc_holo, based on subclass’s calc_field
- calc_cross_sections(scatterer, medium_index=None, illum_wavelen=None, illum_polarization=None, theory='auto')
Calculate scattering, absorption, and extinction cross sections, and asymmetry parameter <cos heta>.
- Parameters:
scatterer (
scattererobject) – (possibly composite) scatterer for which to compute scatteringmedium_index (float or complex) – Refractive index of the medium in which the scatter is imbedded
illum_wavelen (float or ndarray(float)) – Wavelength of illumination light. If illum_wavelen is an array result will add a dimension and have all wavelengths
theory (
theoryobject (optional)) – Scattering theory object to use for the calculation. This is optional if there is a clear choice of theory for your scatterer. If there is not a clear choice, calc_cross_sections will error out and ask you to specify a theory
- Returns:
cross_sections – Dimensional scattering, absorption, and extinction cross sections, and <cos theta>
- Return type:
array (4)
- calc_field(detector, scatterer, medium_index=None, illum_wavelen=None, illum_polarization=None, theory='auto')
Calculate the scattered fields from a scatterer illuminated by a reference wave.
- Parameters:
detector (xarray object) – The detector points and calculation metadata used to calculate the scattered fields.
scatterer (
scattererobject) – (possibly composite) scatterer for which to compute scatteringmedium_index (float or complex) – Refractive index of the medium in which the scatter is imbedded
illum_wavelen (float or ndarray(float)) – Wavelength of illumination light. If illum_wavelen is an array result will add a dimension and have all wavelengths
theory (
theoryobject (optional)) – Scattering theory object to use for the calculation. This is optional if there is a clear choice of theory for your scatterer. If there is not a clear choice, calc_field will error out and ask you to specify a theory
- Returns:
e_field – Calculated hologram from the given distribution of spheres
- Return type:
Vectorobject
- calc_holo(detector, scatterer, medium_index=None, illum_wavelen=None, illum_polarization=None, theory='auto', scaling=1.0)
Calculate hologram formed by interference between scattered fields and a reference wave
- Parameters:
detector (xarray object) – The detector points and calculation metadata used to calculate the hologram.
scatterer (
scattererobject) – (possibly composite) scatterer for which to compute scatteringmedium_index (float or complex) – Refractive index of the medium in which the scatter is imbedded
illum_wavelen (float or ndarray(float)) – Wavelength of illumination light. If illum_wavelen is an array result will add a dimension and have all wavelengths
theory (
theoryobject (optional)) – Scattering theory object to use for the calculation. This is optional if there is a clear choice of theory for your scatterer. If there is not a clear choice, calc_holo will error out and ask you to specify a theoryscaling (scaling value (alpha) for amplitude of reference wave)
- Returns:
holo – Calculated hologram from the given distribution of spheres
- Return type:
xarray.DataArray
- calc_intensity(detector, scatterer, medium_index=None, illum_wavelen=None, illum_polarization=None, theory='auto')
Calculate intensity from the scattered field at a set of locations
- Parameters:
detector (xarray object) – The detector points and calculation metadata used to calculate the intensity.
scatterer (
scattererobject) – (possibly composite) scatterer for which to compute scatteringmedium_index (float or complex) – Refractive index of the medium in which the scatter is imbedded
illum_wavelen (float or ndarray(float)) – Wavelength of illumination light. If illum_wavelen is an array result will add a dimension and have all wavelengths
theory (
theoryobject (optional)) – Scattering theory object to use for the calculation. This is optional if there is a clear choice of theory for your scatterer. If there is not a clear choice, calc_intensity will error out and ask you to specify a theory
- Returns:
inten – scattered intensity
- Return type:
xarray.DataArray
- calc_scat_matrix(detector, scatterer, medium_index=None, illum_wavelen=None, theory='auto')
Compute farfield scattering matrices for scatterer
- Parameters:
detector (xarray object) – The detector points and calculation metadata used to calculate the scattering matrices.
scatterer (
holopy.scattering.scattererobject) – (possibly composite) scatterer for which to compute scatteringmedium_index (float or complex) – Refractive index of the medium in which the scatter is imbedded
illum_wavelen (float or ndarray(float)) – Wavelength of illumination light. If illum_wavelen is an array result will add a dimension and have all wavelengths
theory (
theoryobject (optional)) – Scattering theory object to use for the calculation. This is optional if there is a clear choice of theory for your scatterer. If there is not a clear choice, calc_scat_matrix will error out and ask you to specify a theory
- Returns:
scat_matr – Scattering matrices at specified positions
- Return type:
Marray
- determine_default_theory_for(scatterer)
- finalize(detector, result)
- interpret_theory(scatterer, theory='auto')
- prep_schema(detector, medium_index, illum_wavelen, illum_polarization)
- scattered_field_to_hologram(scat, ref)
Calculate a hologram from an E-field
- Parameters:
scat (
VectorGrid) – The scattered (object) fieldref (xarray[vector]]) – The reference field
- validate_scatterer(scatterer)