Tmatrix classes¶
This section describes the Tmatrix classes currently
implemented in the toolbox.
Tmatrix¶
- class ott.Tmatrix(data, type)¶
Class representing the 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 object, you should inherit from this class when defining your own T-matrix creation methods. This class doesn’t inherit from
doubleorsingle, instead the internal array type can be set at creation allowing the use of different data types such assparseorgpuArray.This class is not a handle class, therefore, when using the class methods you need to store the resulting T-matrix output, for example:
tmatrix = ott.Tmatrix(); new_tmatrix = tmatrix.scattered();
- Properties
data – The T-matrix this class encapsulates
type (enum) – Type of T-matrix (total or scattered)
- Methods
total() – Convert to a total-field T-matrix
scattered() – Convert to a scattered-field T-matrix
real – Extract real part of T-matrix
imag – Extract imaginary part of T-matrix
- Static methods
simple() – Construct a simple particle T-matrix
See also Tmatrix, simple,
ott.TmatrixMie.- Tmatrix(data, type)¶
Construct a new T-matrix object.
- Usage
TMATRIX() leaves the data uninitialised.
TMATRIX(data, type) initializes the data with the matrix data.
- Parameters
data (numeric) – The T-matrix data. Typically a sparse or full matrix.
type (enum) – Type of T-matrix. Must be ‘internal’, ‘scattered’ or ‘total’.
- Example
The following example creates an identity T-matrix which represents a particle which doesn’t scatter light:
data = eye(16); tmatrix = ott.Tmatrix(data, 'total');
- static simple(shape, varargin)¶
Constructs a T-matrix for different simple particle shapes. This method creates an instance of one of the other T-matrix classes and is here only as a helper method.
- Usage
SIMPLE(shape) constructs a new simple T-matrix for the given
ott.shapes.Shapeobject.SIMPLE(name, parameters) constructs a new T-matrix for the shape described by the name and parameters.
- Supported shape names [parameters]
‘sphere’ – Spherical (or layered sphere) [ radius ]
‘cylinder’ – z-axis aligned cylinder [ radius height ]
‘ellipsoid’ – Ellipsoid [ a b c]
‘superellipsoid’ – Superellipsoid [ a b c e n ]
‘cone-tipped-cylinder’ – [ radius height cone_height ]
‘cube’ – Cube [ width ]
‘axisym’ – Axis-symetric particle [ rho(:) z(:) ]
- Optional named arguments
method (enum) – Allows you to choose the preferred method to use for T-matrix calculation. Supported methods are ‘mie’, smarties’, ‘dda’, ‘ebcm’, and ‘pm’. Default: ‘’.
method_tol (numeric) – Specifies the error tolerances, a number between (0, 1] to use for method selection. Smaller values correspond to more accurate methods. Default: [].
- Example
The following example creates a T-matrix for a cube with side length of 1 micron using a guess at the best available method. Illumination wavelength is 1064 nm, relative index 1.5/1.33:
tmatrix = ott.Tmatrix.simple('cube', 1.0e-6, ... 'wavelength0', 1064e-9, ... 'index_medium', 1.33, 'index_particle', 1.5);
TmatrixEbcm¶
Constructs a T-matrix using extended boundary conditions method.
- class ott.TmatrixEbcm(rtp, normals, ds, varargin)¶
Constructs a T-matrix using extended boundary conditions method. Inherits from
ott.Tmatrix.- TmatrixEbcm properties:
k_medium Wavenumber in the surrounding medium k_particle Wavenumber of the particle
This class is based on tmatrix_ebcm_axisym.m from ottv1.
See also TmatrixEbcm
TmatrixMie¶
Construct T-matrix from Mie scattering coefficients.
- class ott.TmatrixMie(radius, varargin)¶
TmatrixMie construct T-matrix from Mie scattering coefficients
- TmatrixMie properties:
radius The radius of the sphere the T-matrix represents k_medium Wavenumber in the trapping medium k_particle Wavenumber of the particle
This class is based on tmatrix_mie.m and tmatrix_mie_layered.m from ottv1.
TmatrixPm¶
Constructs a T-matrix using the point matching method.
- class ott.TmatrixPm(rtp, normals, varargin)¶
TmatrixPm constructs a T-matrix using the point matching method
- TmatrixPm properties:
k_medium Wavenumber in the surrounding medium k_particle Wavenumber of the particle
- TmatrixPm methods:
getInternal Get the internal T-matrix
This class is based on tmatrix_pm.m from ottv1.
TmatrixSmarties¶
- class ott.TmatrixSmarties(a, c, varargin)¶
Constructs a T-matrix using SMARTIES. Inherits from
ott.TmatrixSMARTIES is a method for calculating T-matrices for spheroids.
This class requires the SMARTIES toolbox has been loaded onto the path. Details about the toolbox can be found in:
Somerville, Auguié, Le Ru. JQSRT, Volume 174, May 2016, Pages 39-55. https://doi.org/10.1016/j.jqsrt.2016.01.005
See also TmatrixSmarties
TmatrixDda¶
- class ott.TmatrixDda(xyz, varargin)¶
Constructs a T-matrix using discrete dipole approximation. Inherits from
ott.Tmatrix.To construct a T-matrix with DDA, either the simple interface or the class constructor can be used. Using the simple interface, the following should produce something similar to
TmatrixMie:Tmatrix = ott.TmatrixDda.simple('sphere', 0.1, 'index_relative', 1.2);
The DDA method requires a lot of memory to calculate the T-matrix. 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). For these particles, consider using Geometric Optics or Finite Difference Time Domain method.
See also TmatrixDda, simple.
- TmatrixDda(xyz, varargin)¶
Calculates T-matrix using discrete dipole approximation.
- Usage
TmatrixDda(xyz, …) calculates the T-matrix for the particle described by voxels xyz. xyz is a 3xN matrix of coordinates for each voxel.
The method supports homogenous and inhomogenous particles. For homogeneous particles, specify the material as a scalar, 3x1 vector or 3x3 polarizability matrix. For inhomogeneous particles use a N, 3xN or 3x3N vector/matrix.
- Optional named parameters
Nmax [r,c] Size of the T-matrix to generate. Default:
ott.utils.ka2nmax(max_radius*k_medium)k_medium (numeric) – Wavenumber in medium
wavelength_medium (numeric) – Wavelength in medium
index_medium (numeric) – Refractive index in medium. Default:
k_medium = 2*pik_particle (numeric) – Wavenumber in particle
wavelength_particle (numeric) – Wavelength in particle
index_particle (numeric) – Refractive index in particle. Default:
k_particle = 2*pi*index_relativepolarizability (enum|numeric) – Polarizability or method name to use to calculate from relative refractive index. Default: ‘LDR’. Supported methods: ‘LDR’, ‘FCD’, ‘CM’.
index_relative (numeric) – Relative refractive index. Default: 1.0
wavelength0 (numeric) – Wavelength in vacuum. Default: 1.0
spacing (numeric) – spacing for estimating Nmax and calculating the polarizability. Only required when polarizability is non-numeric. Default:
[]z_mirror_symmetry (logical) – If z-mirror symmetry should be used. All voxels less than 0 are ignored. Default:
false.z_rotational_symmetry (numeric) – z-rotational symmetry. Degree of rotational symmetry. Objects with no rotational symetry should set this to 1. Default:
1.low_memory (logical) – If true, the DDA implementation favours additional calculations over additional memory use allowing simualtion of larger particles. Default: false. Only applicable with z_mirror_symmetry and z_rotational_symmetry.
modes (numeric) – Specifies the modes to include in the calculated T-matrix. Can either be a Nx2 matrix or a N element vector with the (n, m) modes or combined index modes respectively. Default:
[].use_nearfield (logical) – If true, uses near-field point matching. Default: false.
use_iterative (logical) – If true, uses an iterative solver. Default:
false.verbose (logical) – Display additional information. Doesn’t affect the display of the progress callback. Default: false.
- static simple(shape, varargin)¶
Construct a T-matrix using DDA for simple shapes.
- Usage
simple(shape, …) constructs a new simple T-matrix for the given
ott.shapes.Shapeobject.simple(name, parameters, …) constructs a new T-matrix for the shape described by the name and parameters. For supported shape names, see
ott.shapes.Shape.simple.- Optional named arguments:
spacing (numeric) – Spacing between dipoles. Default:
wavelength_particle/20
For other named parameters, see
TmatrixDda().