tmatrix Package

The tmatrix package provides methods for calculating how particles scatter light via the T-matrix method. The main class in this package is Tmatrix which provides methods for manipulating the T-matrix data. The data is stored internally as a matrix, allowing the type to be changed to any valid Matlab matrix type (such as gpuArray or sparse).

Other classes in this package provide methods for calculating T-matrices using different approximations, these methods are summarised in Fig. 8. Most of these classes have a FromShape method which can be used to calculate the T-matrix from a geometric shape description. Once the T-matrix is calculated, it can be multiplied by a ott.bsc.Bsc object to calculate the scattered or internal fields.

The package also contains two sub-packages: smatries contains components of SMARTIES which are used in Smarties; and the dda sub-package contains components of a discrete dipole approximation implementation that will eventually move into OTTv2 (currently used only by Dda).

Graphical (overview) table of contents for drag package

Fig. 8 Graphical display of different particle shapes and the corresponding T-matrix methods that might be suitable for modelling them. The accuracy of the resulting T-matrix often depends on the size and aspect ratio of the particle as well as the choosen method.

Tmatrix

class ott.tmatrix.Tmatrix(varargin)

Class representing T-matrix of a scattering particle or lens. This class can either be instantiated directly or used as a base class for defining custom T-matrix types.

This class is the base class for all other T-matrix objects, you should inherit from this class when defining your own T-matrix creation methods. This class doesn’t inherit from double or single, instead the internal array type can be set at creation allowing the use of different data types such as sparse or gpuArray.

Properties
  • data – The T-matrix this class encapsulates

  • type – Type of T-matrix (total, scattered or internal)

  • Nmax – Size of the T-matrix data (number of multipoles)

  • total – Total-field instance of the T-matrix

  • scattered – Scattered-field instance of the T-matrix

Methods
  • Tmatrix – Class constructor

  • issparse – Returns true if the internal data is sparse

  • full – Convert internal data to full

  • sparse – Convert internal data to sparse

  • makeSparse – Make the data sparse (with additional options)

  • setNmax – Resize data to desired Nmax

  • shrinkNmax – Reduce Nmax while preserving column/row power

  • gpuArray – Make T-matrix a gpuArray

  • gather – Apply gather to T-matrix data

  • setType – Set the T-matrix type property (doesn’t change data)

  • columnCheck – Calculate and check T-matrix column power

  • mergeCols – Merge columns of tmatrices

Mathematical and matrix operations
  • times – Scalar multiplication of data

  • mtimes – Scalar and matrix multiplication of data

  • rdivide – Scalar division of data

  • mrdivide – Scalar division of data

  • uminus – Negation of data

  • minus – Subtraction of data

  • plus – Addition of data

  • real – Extract real part of T-matrix

  • imag – Extract imaginary part of T-matrix

  • diag – Extract the diagonal of the T-matrix

Static methods
  • FromShape – Take a guess at a suitable T-matrix method

  • SmartCylinder – Smart method selection for cylindrical particles

Casts
  • ott.bsc.Bsc – Convert each T-matrix column to beam vector

  • ott.tmatrix.Tmatrix – Downcast T-matrix superclass to base class

See also Tmatrix() and ott.tmatrix.Mie.

static FromShape(shape, varargin)

Take a guess at a suitable T-matrix method.

The T-matrix can only be calculated easily and accurately for a few very specific cases; as such, this function defaults to DDA for most particles which may result in very slow calculations (or failures due to memory limitations). The resulting T-matrices may not be accurate and it is recommended to inspect the fields and compare results with another method.

This method does the following
  • spheres – Uses Mie.

  • spheroids – Uses Smarties

  • cylinders – Uses SmartCylinder(). This method may also work well for other semi-elongated rotationally symmetry shapes.

  • rotationally symmetric – Uses Ebcm.

  • star shaped – Uses Pointmatch.

  • otherwise – Uses Dda.

For many types of particles it would be better to use another method (such as geometric optics for very large particles). This method may change in future releases when other methods are added or when limits of existing methods are further explored.

Usage

tmatrix = ott.tmatrix.Tmatrix.FromShape(shape, index_relative)

Parameters
  • shape (ott.shape.Shape) – Shape to generate T-matrix for. Shape dimensions should be in units of wavelength.

  • index_relative (numeric) – Relative refractive index.

  • internal (logical) – If the T-matrix should be an internal T-matrix. Default: false.

Tmatrix(varargin)

Construct a new T-matrix object.

Usage

tmatrix = Tmatrix(…) New empty T-matrix. Leaves the data uninitialised.

tmatrix = Tmatrix(data, …) Initializes the data with the matrix data.

Parameters
  • data (NxM numeric | cell) – The T-matrix data. Typically a sparse or full matrix. Data must be empty or valid T-matrix size. If cell, describes T-matrix array and elements must be a cell array of NxM matrices.

Optional named arguments
  • type (enum) – Type of T-matrix. Must be ‘internal’, ‘scattered’ or ‘total’. Default: 'scattered'.

Example

The following example creates an identity T-matrix which represents a particle which doesn’t scatter light:

data = eye(16);
tmatrix = ott.scat.vswf.Tmatrix(data, 'type', 'total');
columnCheck(tmatrix, varargin)

Check the power in each column

For a non-absorbing total-field T-matrix, the power in each column should add up to unity (power should be conserved).

Usage

tmatrix.columnCheck(…) Raises a warning if the power drops bellow a threshold.

column_power = tmatrix.columnCheck(…) Returns the power in each column.

Optional named arguments
  • threshold (numeric) – Threshold for power loss warnings/errors. Default: 1.0e-3.

  • action (enum) – Action to take if power lost/gained. Can be ‘warn’, ‘error’ or ‘none’. If no outputs present, the default behaviour is ‘warn’. Otherwise ‘none’.

diag(tmatrix, k)

Extract the T-matrix diagonal

Usage

d = diag(tmatrix) Returns the diagonal in vector format

d = diag(tmatrix, K) Extracts the K-th diagonal of the T-matrix.

[dA, dB] = diag(…) Returns the diagonal of the upper and lower parts separately.

full(tmatrix)

Convert the data to a full matrix

Usage

tmatrix = full(tmatrix)

gather(tmatrix)

Apply gather to data.

If the data is a gpuArray, returns a copy of the data in the local workspace with data transferred from the GPU.

Usage

tmatrix = gather(tmatrix)

gpuArray(tmatrix)

Copies the tmatrix data to the GPU

Usage

tmatrix = gpuArray(tmatrix)

imag(tmatrix)

Extract imaginary part of T-matrix

Usage

tmatrix = imag(tmatrix);

issparse(tmatrix)

Returns true if the data is sparse

Usage

b = issparse(tmatrix)

makeSparse(tmatrix, varargin)

Make the T-matrix data sparse by removing near-zero power elements

Treats each column as a beam shape vector and applies ott.bsc.Bsc.makeSparse() to each column.

Usage

tmatrix = tmatrix.makeSparse(…)

Optional named arguments
  • AbsTol (numeric) – Absolute tolerance for removing elements. Default: [].

  • RelTol (numeric) – Relative tolerance for removing elements. Power is relative to power in each column. Default: 1.0e-15.

If both AbsTol and RelTol are specified, only elements satisfying both conditions are kept.

mergeCols(tmatrix, tmatrix2, ci)

Merge columns of two T-matrices

Usage

tmatrix = tmatrix.mergeCols(tmatrix2, ci) Keeps all cols of the first T-matrix except thouse replaced by tmatrix2 (specified by ci).

Parameters
  • tmatrix1, tmatrix2 – T-matrices to merge

  • ci – (N numeric) Combined indices of T-matrix columns to replace. For each ci entry, keeps 2 columns from tmatrix2 (i.e., both TE and TM modes).

minus(a, b)

Minus operation on tmatrix

Usage

tmatrix = tmatrix1 - tmatrix2

tmatrix = other - tmatrix tmatrix = tmatrix - other

Note: Uses the uminus operation of the second argument.

mtimes(a, b)

Matrix and scalar multiplication of T-matrix data

Usage

beam = tmatrix * beam Calculate how a beam is scattered by the T-matrix, increase beam or T-matrix Nmax if required.

tmatrix = tmatrix1 * tmatrix2 Multiply T-matrix data, increasing Nmax if required.

vector = vector * tmatrix vector = tmatrix * vector

tmatrix = scalar * tmatrix tmatrix = tmatrix * scalar

plus(a, b)

Apply plus operation on T-matrix data

Usage

tmatrix = tmatrix1 + tmatrix2;

tmatrix = other + tmatrix; tmatrix = tmatrix + other;

rdivide(tmatrix, o)

Scalar division of T-matrix data

Usage

tmatrix = tmatrix ./ scalar

real(tmatrix)

Extract real part of T-matrix

Usage

tmatrix = real(tmatrix);

setNmax(tmatrix, nmax, varargin)

Resize the T-matrix, with additional options

Usage

tmarix = tmatrix.setNmax(nmax, …) or tmatrix.Nmax = nmax Set the Nmax, a optional warning is issued if truncation occurs.

Parameters
  • Nmax (1 | 2 numeric) – Nmax for both dimensions or vector with Nmax for [rows, cols].

Optional named arguments
  • AbsTol (numeric) – Absolute tolerance for removing elements. Default: [].

  • RelTol (numeric) – Relative tolerance for removing rows. Power is relative to power in each column. Default: 1.0e-15.

  • ColTol (numeric) – Absolute tolerance for removing columns. Default: 1.0e-15.

  • powerloss (enum) – Action to take when column power is lost. Can be one of ‘ignore’, ‘warn’ or ‘error’. Default: 'warn'.

setType(tmatrix, val)

Set the T-matrix type paramter (without raising a warning)

Usage

tmatrix = tmatrix.setType(val);

shrinkNmax(tmatrix, varargin)

Shrink the size of the T-matrix while preserving power

Converts to a scattered or internal T-matrix and then removes columns with no significant power and rows by passing each column to ott.bsc.Bsc.shrinkNmax().

Usage

tmatrix = tmatrix.shrinkNmax(…)

Optional named arguments
  • AbsTol (numeric) – Absolute tolerance for removing elements. Default: [].

  • RelTol (numeric) – Relative tolerance for removing rows. Power is relative to power in each column. Default: 1.0e-15.

  • ColTol (numeric) – Absolute tolerance for removing columns. Default: 1.0e-15.

sparse(tmatrix)

Convert the data to a sparse matrix

This function doesn’t change the data. For a method that removes near-zeros elements, see makeSparse().

Usage

tmatrix = sparse(tmatrix)

times(a, b)

Element-wise multiplication of T-matrix data

Usage

tmatrix = tmatrix1 .* tmatrix2;

tmatrix = other .* tmatrix; tmatrix = tmatrix .* other;

Dda

class ott.tmatrix.Dda(dda, varargin)

Construct a T-matrix using the discrete dipole approximation. Inherits from ott.tmatrix.Tmatrix.

This method calculates how each VSWF mode is scattered and then uses point-matching to calculate each column of the T-matrix. The field calculations are done using the discrete dipole approximation, but the same approach could be used for T-matrix calculation via any other field calculation method.

The current DDA implementation requires a lot of memory. Most small desktop computers will be unable to calculate T-matrices for large particles (i.e., particles larger than a couple of wavelengths in diameter using 20 dipoles per wavelength). The aim of the next OTT release is to include geometric optics and finite difference time domain as alternative methods for force calculation with these larger particles.

Properties
  • dda – DDA instance used for field calculation

  • pmrtp – (3xN numeric) Locations for point matching

Static methods
  • FromShape – Construct from a geometric shape

  • DefaultPmrtp – Build default pmrtp locatiosn for point matching

  • DefaultProgressCallback – Default progress callback for method

See also Dda(), ott.tmatrix.dda.

Dda(dda, varargin)

Construct a T-matrix using a DDA simulation for field calculation

Usage

tmatrix = Dda(dda, …)

[tmatrix, incData, pmData] = Dda(dda, …) Also returns the VSWF data structures used for calculating the incident field and the point matching field.

Parameters
  • dda – (ott.tmatrix.dda.Dda instance) The DDA instance used for field calculations.

Optional named arguments
  • Nmax – (numeric) Size of the VSWF expansion used for the T-matrix point matching (determines T-matrix rows). Default: ott.utils.ka2nmax(2*pi*shape.maxRadius) (may need different values to give convergence for some shapes).

  • ci – (N numeric) Number of modes to calculate scattering for. This determines number of T-matrix columns. Default: 1:ott.utils.combined_index(Nmax, Nmax) (all modes).

  • pmrtp – (3xN numeric) Coordinatse for point matching. Radial coordinate must all be finite or all be Inf. Default uses ott.tmatrix.Dda.DefaultPmrtp.

  • incData (ott.utils.VswfData) – Data structure for repeated incident field calculations. Default: ott.utils.VswfData().

  • pmData (ott.utils.VswfData) – Data structure for repeated point matching field calculations. Default: ott.utils.VswfData().

Unmatched parameters passed to calculate_columns().

static DefaultPmrtp(Nmax, varargin)

Build default grid of points for point matching

Usage

pmrtp = DefaultPmrtp(Nmax, …)

Parameters
  • Nmax – (numeric) Row Nmax for generated T-matrix.

Optional named parameters
  • radius – (numeric) Radius for point matching locations. Must be positive scalar or Inf for far-field point matching. Default: Inf.

  • angulargrid – ({theta, phi}) Angular grid of points for calculation of radii. Default is equally spaced angles with the number of points determined by Nmax.

  • xySymmetry – (logical) If the generated grid should be for mirror symmetry DDA. Default: false.

  • zRotSymmetry – (numeric) If the generated grid should be for rotationally symmetric DDA. Default: 1.

static DefaultProgressCallback(data)

Default progress callback for Dda

Prints the progress to the terminal.

Usage

DefaultProgressCallback(data)

Parameters
  • data (struct) – Structure with two fields: index and total.

static FromShape(shape, varargin)

Construct a T-matrix from a geometric shape

Usage

tmatrix = ott.tmatrix.Dda.FromShape(shape, …)

[tmatrix, incData, pmData] = ott.tmatrix.Dda.FromShape(…) Returns the field calculation data for repeated calculations.

Optional named arguments
  • spacing – (numeric) – Dipole spacing in wavelength units. Default: 1/20.

  • polarizability – (function_handle | 3x3 numeric) Method to calculate polarizability or 3x3 tensor for homogeneous material. Default: @(xyz, spacing, ri) polarizability.LDR(spacing, ri)

  • index_relative – (function_handle | numeric) Method to calculate relative refractive index or homogeneous value. Ignored if polarizability is a 3x3 tensor.

  • low_memory – (logical) If we should use the low memory DDA implementation. Default: false.

Additional parameters passed to class constructor.

Ebcm

class ott.tmatrix.Ebcm(varargin)

Constructs a T-matrix using extended boundary conditions method. Inherits from ott.tmatrix.Tmatrix.

Implements the extended boundary conditions methods for rotationally symmetric homogeneous particles.

Properties
  • points – (2xN numeric) Surface points [r; theta]

  • normals – (2xN numeric) Normals at points [nr; ntheta]

  • areas – (1xN numeric) Conic section surface areas

  • index_relative – (numeric) Relative refractive index of particle.

  • xySymmetry – (logical) True if using xy-mirror optimisations

  • invMethod – Inversion method used for T-matrix calculation

Static methods
  • FromShape – Construct from geometric shape object.

Additional methods/properties inherited from Tmatrix.

This class is based on tmatrix_ebcm_axisym.m from OTTv1.

Ebcm(varargin)

Calculates T-matrix using extended boundary condition method

Usage

tmatrix = Ebcm(points, normals, area, index_relative, …) Calculate external T-matrix (unless internal is true)

[external, internal, data] = Ebcm(…) Calculate both the internal and external T-matrices, and return the VswfData structure used during calculations.

Parameters
  • points (2xN numeric) – Coordinates describing surface. Spherical coordinates (omitting azimuthal angle: [r; theta]).

  • normals (2xN numeric) – Normals at surface points.

  • areas (N numeric) – Area elements at surface points.

  • index_relative (numeric) – Particle relative refractive index.

Optional named parameters
  • xySymmetry (logical) – If calculation should use xy-mirror symmetry optimisations. Default: false. Doesn’t check points describe valid mirror symmetric shape.

  • invMethod (enum|function handle) – Inversion method for T-matrix calculation. Currently supported methods are ‘forwardslash’, ‘backslash’, ‘inv’, ‘pinv’. A custom inversion method can be specified using a function handle with the format inv_method(RgQ, Q) which returns the T-matrix data. Default: forwardslash. Only used for external T-matrix calculation, may change in future.

  • internal (logical) – If true, the returned T-matrix is an internal T-matrix. Ignored for two outputs. Default: false.

  • Nmax (numeric) – Size of the VSWF expansion used for the T-matrix calculation. In some cases it can be reduced after construction. Default: ott.utis.ka2nmax(2*pi*shape.maxRadius) (may need different values to give convergence for some shapes).

  • verbose (logical) – If true, outputs the condition number of the Q and RgQ matrices. Default: false.

  • data (ott.utils.VswfData) – Field data for repeated field calculation. Default is an empty VswfData structure.

static FromShape(varargin)

Construct a T-matrix using EBCM from a shape object.

If the shape object is not of type ott.shape.AxisymLerp, first casts the shape to this type using the default cast (if supported, otherwise raises an error).

Usage

tmatrix = Ebcm.FromAxisymInterpShape(shape, index_relative, …) Calculate external T-matrix (unless internal is true)

[external, internal] = Ebcm.FromAxisymInterpShape(…) Calculate both the internal and external T-matrices.

Parameters
  • shape (ott.shape.Shape) – Description of shape geometry. Object must be a AxisymInterp or be castable to AxisymInterp.

  • index_relative (numeric) – Particle relative refractive index.

See Ebcm() for additional named parameters.

Mie

class ott.tmatrix.Mie(varargin)

Construct T-matrix with Mie scattering coefficients. Inherits from ott.tmatrix.Homogeneous.

The Mie coefficients describe the scattering of a sphere. They can also be used to give a reasonable estimate of the force for non-spherical particles when no other suitable method is available.

This class supports both homogeneous dielectric (conductive and non-conductive) and magnetic isotropic materials. For other materials, consider Dda or MieLayered (for layered spheres).

Properties
  • radius – Radius of sphere

  • index_relative – Relative refractive index

  • relative_permeability – Relative permeability

Static methods
  • FromShape – Uses ShapeMaxRadius but raises a warning

  • ShapeVolume – Construct with radius set from particle volume

  • ShapeMaxRadius – Construct with radius set from particle max radius

See base class for additional methods/properties.

This class is based on tmatrix_mie.m and tmatrix_mie_layered.m from ottv1.

static FromShape(shape, varargin)

Construct a T-matrix from a shape object.

Uses ShapeMaxRadius() but raises a warning if the shape isn’t a sphere.

Usage

tmatrix = Mie.FromShape(shape, …)

Parameters
  • shape (ott.shape.Shape) – The shape input.

All other parameters passed to constructor.

Mie(varargin)

Construct a new Mie T-matrix for a sphere.

Usage

tmatirx = Mie(radius, index_relative, …) Calculate the external T-matrix (unless internal = true).

[external, internal] = Mie(radius, index_relative, …) Calculate both the internal and external T-matrices.

Parameters
  • radius (numeric) – Radius of the sphere (in wavelength units).

  • index_relative (numeric) – The relative refractive index of the particle compared to the surrounding medium.

Optional named parameters
  • relative_permeability (numeric) – Relative permeability of the particle compared to surrounding medium. Default: 1.0.

  • Nmax (numeric) – Size of the VSWF expansion used for the T-matrix calculation. Default: ott.utis.ka2nmax(2*pi*radius) (external) or ott.utis.ka2nmax(2*pi*radius*index_relative) (internal).

  • internal (logical) – If true, the returned T-matrix is an internal T-matrix. Ignored for two outputs. Default: false.

static ShapeMaxRadius(shape, varargin)

Construct Mie T-matrix with radius from shape max radius

Usage

tmatrix = Mie.ShapeMaxRadius(shape, …)

Parameters
  • shape (ott.shape.Shape) – The shape input.

All other parameters passed to constructor.

static ShapeVolume(shape, varargin)

Construct Mie T-matrix with radius from shape volume

Usage

tmatrix = Mie.ShapeVolume(shape, …)

Parameters
  • shape (ott.shape.Shape) – The shape input.

All other parameters passed to constructor.

MieLayered

class ott.tmatrix.MieLayered(varargin)

Construct T-matrices for a layered sphere. Inherits from ott.tmatrix.Tmatrix.

This class implements the first part of

“Improved recursive algorithm for light scattering by a multilayered sphere”, Wen Yang, Applied Optics 42(9), 2003

and can be used to model layered spherical particles.

Properties
  • radii – Radii of sphere layers (inside to outside)

  • relative_indices – Relative refractive indices (inside to outside)

Static methods
  • FromShape – Construct layered sphere from array of shapes

See base class for additional methods/properties.

This class is based on tmatrix_mie_layered.m from ottv1.

static FromShape(shape, relative_indices, varargin)

Construct a T-matrix from a shape array

Uses the maxRadius property of each shape in the shape array for the sphere radii.

Shapes radii must be in ascending order. If the shapes aren’t all centred, raises a warning.

Usage

tmatrix = MieLayered.FromShape(shape, relative_indices, …)

All other parameters are passed to MieLayered.

MieLayered(varargin)

Construct a new Mie layered T-matrix

Usage

tmatirx = MieLayered(radii, relativeMedium, …) Calculate the external T-matrix (unless internal = true).

[external, internal] = Mie(radii, relativeMedium, …) Calculate both the internal and external T-matrices.

Parameters
  • radii (numeric) – Radius of sphere layers (in relative units unless optional parameter wavelength is specified). Radii should be in ascending order.

  • relativeMedium (ott.beam.medium.RelativeMedium) – The relative medium describing the particle material layers. Particle must be isotropic homogeneous magnetic or dielectric. The number of layers should match number of radii.

Optional named parameters
  • wavelength (numeric) – Used to convert radius input to relative units, i.e. radius_rel = radius ./ wavelength. This parameter not used for setting the T-matrix material. Default: 1.0 (i.e., radius is already in relative units).

  • Nmax (numeric) – Size of the VSWF expansion used for the T-matrix calculation. Can be reduced after construction. Default: 100 (should be numerically stable).

  • internal (logical) – If true, the returned T-matrix is an internal T-matrix. Ignored for two outputs. Default: false.

Pointmatch

class ott.tmatrix.Pointmatch(varargin)

Constructs a T-matrix using the point matching method. Inherits from ott.tmatrix.Homogeneous.

The point matching method is described in

T. A. Nieminen, H. Rubinsztein-Dunlop, N. R. Heckenberg JQSRT 79-80, 1019-1029 (2003), 10.1016/S0022-4073(02)00336-9

The method can be used to construct T-matrices for star shaped particles with aspect ratios close to unity. Supports homogeneous isotropic materials. This implementation includes symmetry optimisations for rotationally symmetric and mirror symmetric particles.

Properties
  • index_relative – Relative refractive index of particle

  • rtp – Locations used for point matching

  • nrtp – Surface normals at rtp-locations (spherical coords.)

  • zRotSymmetry – Z-axis rotational symmetry (1 = no symmetry)

  • xySymmetry – XY mirror symmetry (logical)

Static methods
  • DefaultProgressCallback – Default progress call-back method

  • FromShape – Construct T-matrix from star shape

This class is based on tmatrix_pm.m from ottv1.

static DefaultProgressCallback(data)

Default progress callback for Pointmatch

Prints the progress to the terminal.

Usage

DefaultProgressCallback(data)

Parameters
  • data (struct) – Structure with three fields: stage (either ‘setup’ or ‘inv’), iteration (numeric), and total (numeric).

static FromShape(shape, varargin)

Construct T-matrix using point matching from a star shaped object

Usage

tmatrix = Pointmatch.FromShape(shape, index_relative, …) Calculate external T-matrix.

[external, internal] = Pointmatch.FromShape(…) Calculate external and internal T-matrices.

Parameters
  • shape (ott.shape.Shape) – A star-shaped object describing the geometry (must have a valid starRadii method).

  • index_relative (numeric) – Relative refractive index.

Optional named parameters
  • Nmax (numeric) – Size of the VSWF expansion used for the T-matrix calculation. In some cases it can be reduced after construction. Default: ott.utils.ka2nmax(2*pi*shape.maxRadius) (may need different values to give convergence for some shapes).

  • angulargrid ({theta, phi}) – Angular grid of points for calculation of radii. Default is equally spaced angles with the number of points determined by Nmax.

Additional parameters are passed to the class constructor.

Pointmatch(varargin)

Calculates T-matrix using the point matching method.

Usage

tmatrix = Pointmatch(rtp, nrtp, index_relative, …) Calculate external T-matrix.

[external, internal] = Pointmatch(rtp, nrtp, index_relative, …) Calculate external and internal T-matrices.

Parameters
  • rtp (3xN numeric) – Locations of surface points to point-match.

  • nrtp (3xN numeric) – Normals at surface locations. Spherical coordinates.

  • index_relative (numeric) – Relative refractive index.

Optional named parameters
  • Nmax (numeric) – Size of the VSWF expansion used for the T-matrix calculation. In some cases it can be reduced after construction. Default: ott.utis.ka2nmax(2*pi*max(rtp(1, :))) (may need different values to give convergence for some shapes).

  • zRotSymmetry (numeric) – Degree of rotational symmetry about the z-axis. Default: 1 (no symmetry).

  • xySymmetry (logical) – If the particle is mirror symmetry about the xy-plane. Default: false.

  • internal (logical) – If true, the returned T-matrix is an internal T-matrix. Ignored for two outputs. Default: false.

  • progress (function_handle) – Function to call for progress updates during method evaluation. Takes one argument, see DefaultProgressCallback() for more information. Default: [] (for Nmax < 20) and @DefaultProgressCallback (otherwise).

Smarties

class ott.tmatrix.Smarties(varargin)

Constructs a T-matrix using SMARTIES. Inherits from ott.tmatrix.Homogeneous.

SMARTIES is a method for calculating T-matrices for spheroids, full details can be found in

W.R.C. Somerville, B. Auguié, E.C. Le Ru, JQSRT, Volume 174, May 2016, Pages 39-55. https://doi.org/10.1016/j.jqsrt.2016.01.005

SMARTIES is distributed with a Creative Commons Attribution-NonCommercial 4.0 International License. Copyright 2015 Walter Somerville, Baptiste Auguié, and Eric Le Ru. This version of OTT includes a minimal version of SMARTIES for T-matrix calculation, if you use the SMARTIES code, please cite the above paper.

Properties
  • ordinary – Ordinary radius

  • extraordinary – Extra-ordinary radius

  • index_relative – Relative refractive index of particle

Static methods
  • FromShape – Construct a T-matrix from a shape description

See also Smarties, ott.tmatrix.Mie and ott.tmatrix.smarties.

static FromShape(shape, varargin)

Construct a T-matrix using SMARTIES/EBCM for spheroids.

Usage

tmatrix = Smarties.FromShape(shape, index_relative, …)

[texternal, tinternal] = Smarties.FromShape(…)

Parameters
  • shape (ott.shape.Shape) – A spheroid with the extraordinary axis aligned to the z-axis.

  • index_relative (numeric) – The relative refractive index of the particle compared to the surrounding medium.

All other parameters passed to class constructor.

Smarties(varargin)

Construct a T-matrix using SMARTIES/EBCM for spheroids.

Usage

tmatrix = Smarties(ordinary, extraordinary, index_relative…)

[texternal, tinternal] = Smarties(…)

Parameters
  • ordinary (numeric) – Ordinary radius.

  • extraordinary (numeric) – Extraordinary radius.

  • index_relative (numeric) – The relative refractive index of the particle compared to the surrounding medium.

Optional named parameters
  • internal (logical) – If true, the returned T-matrix is an internal T-matrix. Ignored for two outputs. Default: false.

  • Nmax (numeric) – Size of the VSWF expansion used for the T-matrix calculation. Default: ott.utis.ka2nmax(2*pi*max(radius)) (external) or ott.utis.ka2nmax(2*pi*max(radius)*index_relative) (internal).

  • npts (numeric) – Number of points for surface integral. Default: Nmax*Nmax.

  • verbose (logical) – Enables additional output from SMARTIES. Default: false.

Sub-packages

These sub-packages are used by T-matrix calculation methods and only minimal documentation is provided. Their location may change in a future releases (for instance, when DDA is extended to include force calculation methods).

SMARTIES

Contains components of the SMARTIES package for easy installation with the optical tweezers toolbox. Users interested in SMARTIES may want to download the full version, for details see

Somerville, Auguié, Le Ru. JQSRT, Volume 174, May 2016, Pages 39-55. https://doi.org/10.1016/j.jqsrt.2016.01.005

DDA

class ott.tmatrix.dda.Dipole(varargin)

Describes an array of radiating dipoles.

This class stores the dipole locations and polarizations. The scatted fields can be calculated using:

Es = F p

where F describes the locations where fields should be calculated and p describes the polarization of each dipole. The class provides methods for calculating F and Es.

Methods
  • setDipoles – Set the dipole data (location/polarization)

  • efield – Calculate E near-fields

  • hfield – Calculate H near-fields

  • efarfield – Calculate E far-fields

  • hfarfield – Calculate H far-fields

  • efarfield_matrix – Calculate far-field matrix for field calculation

  • hfarfield_matrix – Calculate far-field matrix for field calculation

  • enearfield_matrix – Calculate near-field matrix for field calculation

  • hnearfield_matrix – Calculate near-field matrix for field calculation

  • mtimes – Apply field matrix and calculate fields

Properties
  • location – Dipole locations

  • polarization – Dipole polarization

  • xySymmetry – True if using z-mirror symmetry

  • zRotSymmetry – Order of z-rotational symmetry (0 - infinite)

  • parity – Parity of incident beam

  • rorder – Rotational order of incident beam

  • ndipoles – Number of dipoles in the array

  • nbeams – Number of beams in the array

Dipole(varargin)

Construct a new dipole array

Usage

beam = Dipole(locations, polarization, …) Parameters can also be passed as named arguments.

Parameters
  • locations (3xN numeric) – Locations of dipoles

  • polarization (3NxM) – Dipole polarizations sorted packaged in [x1;y1;z1; x2;y2;z2; …] order.

Optional named parameters
  • parity (enum) – Parity of incident beam (even or odd). Only used when using z_mirror. Default: 'even'.

  • rorder (numeric) – Rotational order of incident beam. Only used when using z_rotation. Default: 0.

  • xySymmetry (logical) – If the particle has z-mirror symmetry. Default: false.

  • zRotSymmetry (numeric) – Order of the particle is z-rotational symmetric. Default: 1. If 0, uses fourth order rotational symmetry (might change in a future release).

efarfield(beam, rtp, varargin)

Calculate the E-field

Usage

E = beam.efarfield(rtp)

Parameters
  • rtp – (3xN | 2xN numeric) Spherical coordinates. Either [radius; theta; phi] or [theta; phi]. Radius is ignored.

Unmatched parameters are passed to efarfield_matrix().

efarfield_matrix(beam, rtp, varargin)

Evaluates the electric far-field matrix for a set of points

Can be applied to the beam to evaluate the scattered fields.

Es = F * beam

Usage

F = beam.farfield_matrix(rtp, …)

Parameters
  • rtp (2xN|3xN numeric) – Far-field coordinates. Can either be [theta; phi] or [r; theta; phi] in which case r is ignored.

Optional named arguments
  • low_memory (logical) – If true, evaluates the low-memory version of F. Default: false.

efield(beam, xyz, varargin)

Calculate the E-field

Evaluates:

Es = F * p
Usage

E = beam.efield(xyz)

Parameters
  • xyz – (3xN numeric) Cartesian coordinates.

Unmatched parameters are passed to enearfield_matrix().

mtimes(F, beam)

Matrix multiplication for calculating scattered fields

Evaluates:

Es = F * p
Usage

Es = F * beam Does not reshape the output.

Parameters
  • F (numeric) – Field matrix calculated using enearfield_matrix(), efarfield_matrix(), hnearfield_matrix(), or hfarfield_matrix().

setDipoles(beam, locations, polarization)

Set the dipole position and polarization data

Usage

beam = beam.setDipoles(location, polarization)

Parameters
  • location (3xN numeric) – Locations of dipoles

  • polarization (3NxM) – Dipole polarizations sorted packaged in [x1;y1;z1; x2;y2;z2; …] order.

class ott.tmatrix.dda.Dda(varargin)

Minimal implementation of the discrete dipole approximation.

This class calculates the dipole polarizations from the interaction matrix A and the incident electric field E:

p = A \ E

The class implements methods for calculating the interaction matrix with optimisation for mirror, rotational symmetry and low memory.

Properties
  • locations – Location of each dipole

  • polarizability – Dipole polarizabilities

  • xySymmetry – True if using z-mirror symmetry

  • zRotSymmetry – Order of z-rotational symmetry (0 - infinite)

  • ndipoles – Number of dipoles in the array

Methods
  • solve – Solve for dipole polarizations

  • interaction_matrix – Calculate the interaction matrix

Static methods
  • FromShape – Construct instance from geometric shape

Dda(varargin)

Construct new DDA solver instance.

Usage

dda = Dda(locations, polarizabilities, …)

Parameters
  • voxels – (3xN numeric) Voxel locations in Cartesian coordinates.

  • polarizabilities – (3x3N numeric) Array of 3x3 dipole polarizability tensors.

Optional named arguments
  • xySymmetry (logical) – If the particle has z-mirror symmetry. Default: false. If true, locations with z<0 are removed.

  • zRotSymmetry (numeric) – Order of the particle z-rotational symmetry. Default: 1. If 0, uses fourth order rotational symmetry (might change in a future release). If zRotSymmetry is 2, removes x<0 points. If zRotSymmetry is mod 4, removes x<0 | y<0 points. No other options supported for now.

static FromShape(shape, varargin)

Construct a DDA instance from a geometric shape.

Usage

dda = ott.tmatrix.dda.Dda.FromShape(shape, …)

Optional named arguments
  • spacing – (numeric) – Dipole spacing in wavelength units. Default: 1/20.

  • polarizability – (function_handle | 1x1 | 3x3 numeric) Particle polarizability. Must be a function handle for inhomogeneous materials or either a scalar or 3x3 matrix for homogeneous materials. Default: @(xyz, spacing, ri) polarizability.LDR(spacing, ri)

  • index_relative – (function_handle | numeric) Method to calculate relative refractive index or homogeneous value. Ignored if polarizability is not a function handle.

Unmatched parameters are passed to FromPolarizability().

interaction_matrix(dda, varargin)

Calculate the interaction matrix assuming memory is limited

This method is rather time consuming and involves a lot of redundant calculations if called repeatedly with different rorder or parity parameters. For a faster but more memory intensive method see DdaHighMemory.

Usage

A = dda.interaction_matrix(rorder, parity)

Optional parameters
  • parity (enum) – Parity of incident beam (even or odd). Only used when using xySymmetry. Default: 'even'.

  • rorder (numeric) – Rotational order of incident beam. Only used when using zRotSymmetry. Default: 0.

solve(dda, Eincident, varargin)

Solve for dipole polarizations

Usage

dipoles = dda.solve(Eincident, …)

Parameters
  • Eincident – (3NxM numeric) Incident electric field at each dipole. Format: [x1;y1;z1;x2;y2;z2;…].

Optional named arguments
  • solver – (function_handle) Solver to use. Good solvers to try include gmres, bicgstab and \. Default: @(A, E) A \ E.

  • multiBeam – (logical) If true, passes the whole beam into the solver, otherwise iterates over each beam. Default: true.

  • parity (enum) – Parity of incident beam (even or odd). Only used when using z_mirror. Default: 'even'.

  • rorder (numeric) – Rotational order of incident beam. Only used when using z_rotation. Default: 0.

class ott.tmatrix.dda.DdaHighMem(varargin)

Slightly faster but more memory intensive version of DDA. Inherits from ott.tmatrix.dda.Dda.

Properties
  • interaction – Stored interaction matrix

Methods
  • update_interaction_matrix – Re-calculate the interaction matrix

  • interaction_matrix – Get a usable interaction matrix

Static methods
  • FromShape – Construct from a geometric shape

For other methods/properties, see Dda.

DdaHighMem(varargin)

Construct DDA instance and pre-compute data for interaction matrix.

Usage

dda = DdaHigMem(dda) Convert an existing DDA instance into a pre-computed instance.

dda = DdaHighMem(locations, interaction, …) Construct a new DDA instance. See base class for parameters.

static FromShape(shape, varargin)

Construct a DDA instance from a geometric shape.

Usage

dda = ott.tmatrix.dda.Dda.FromShape(shape, …)

Optional named arguments
  • spacing – (numeric) – Dipole spacing in wavelength units. Default: 1/20.

  • polarizability – (function_handle | 3x3 numeric) Method to calculate polarizability or 3x3 tensor for homogeneous material. Default: @(xyz, spacing, ri) polarizability.LDR(spacing, ri)

  • index_relative – (function_handle | numeric) Method to calculate relative refractive index or homogeneous value. Ignored if polarizability is a 3x3 tensor.

For further details and options, see Dda.FromShape().

update_interaction_matrix(dda)

Update the interaction matrix data.

This method is called when the class is constructed. If you change the class properties, call this method again.

Usage

dda = dda.update_interaction_matrix()

DDA polarizability methods

ott.tmatrix.dda.polarizability.CM(spacing, index)

Clausius-Mossoti Polarizability

Evaluates:

\alpha = \frac{3}{4\pi} s^3 \frac{n^2 - 1}{n^2 + 2}

Where \(n\) is the refractive index and \(s\) is the dipole spacing. This has units of L^3. To convert to SI units, multiply by \(4\pi\epsilon_0\).

Usage

alpha = CM(spacing, index) Calculates a Nx1 element vector containing the isotropic polarisabilities for N dipoles.

Parameters
  • spacing (numeric scalar) – lattice spacing parameter [L]

  • index (Nx1 numeric) – Relative refractive indices for N dipoles.

ott.tmatrix.dda.polarizability.LDR(spacing, index, varargin)

Lattice dispersion relation polarizablity

Polarizability calculation based on

Draine & Goodman, Beyond Clausius-Mossoti: wave propagation on a polarizable point lattice and the discrete dipole approximation, The Astrophysical Journal, 405:685-697, 1993 March 10

Usage

alpha = LDR(spacing, index, …) Calculates a Nx1 element vector containing the isotropic polarisabilities for N dipoles.

alpha = LDR(spacing, index, kvec, E0, …) As above but specifies the polarisability information for use with plane wave illumination.

Parameters
  • spacing (numeric scalar) – lattice spacing parameter

  • index (Nx1 numeric) – Relative refractive indices for N dipoles.

  • kvec (1x3 numeric) – Wave vector [kx, ky, kz]

  • E0 (1x3 numeric) – E-field polarisation [Ex, Ey, Ez]

Optional named arguments
  • k0 (numeric) – Wavenumber to scale spacing by. Default: 2*pi.

ott.tmatrix.dda.polarizability.FCD(spacing, index, varargin)

Filtered coupled dipole polarizability

Usage

alpha = FCD(spacing, index) Calculates a Nx1 element vector containing the isotropic polarizabilities for N dipoles.

Parameters
  • spacing (numeric scalar) – lattice spacing parameter

  • index (Nx1 numeric) – Relative refractive indices for N dipoles.

Optional named arguments
  • k0 (numeric) – Wavenumber to scale spacing by. Default: 2*pi.