peri.comp

Index of components within peri.comp:

peri.comp.comp

class peri.comp.comp.ParameterGroup(params=None, values=None, ordered=True, category=’param’)

Any object which computes something based on parameters and values can be considered a ParameterGroup. This class provides a common interface since ParameterGroup appears throughout PERI including Components, Priors, States. In the very basic form, a ParameterGroup is a dict or OrderedDict of:

{ parameter_name: parameter_value, ... }

The use of a dictionary is strictly optional – as long as the following methods are provided, the parameters and values may be stored in any format that is convenient:

Parameters:
  • params (string, list of strings) – The names of the parameters, in the proper order
  • values (number, list of numbers) – The values corresponding to the parameter names
  • ordered (boolean (default: True)) – If True, uses an OrderedDict so that parameter order is deterministic independent of number of parameters
  • category (string (default: 'param')) –

    Name of the category associated with this ParameterGroup.

    Warning

    FIXME : should only be a property of Component

get_values(params)

Get the value of a list or single parameter.

Parameters:params (string, list of string) – name of parameters which to retrieve
initargs()

Pickling helper method which returns a dictionary of function parameters which get passed to pickle via __getinitargs__

Returns:arg_dict**kwargs to be passed to the __init__ func after unpickling
Return type:dictionary
nopickle()

Elements of the class that should not be included in pickled objects. If inheriting a new class, should be:

super(Class, self).nopickle() + ['other1', 'other2', ...]
Returns:elements – The name of class member variables which should not be pickled
Return type:list of strings
params

The list of parameters

set_values(params, values)

Directly set the values corresponding to certain parameters. This does not necessarily trigger and update of the calculation,

See also

update() : full update func

update(params, values)

Update the calculation of the class based on a pair or pairs of parameters and associated values.

Parameters:
  • params (string, list of strings) – name of parameters to update
  • values (number, list of numbers) – cooresponding values to update
values

The list of values

class peri.comp.comp.Component(params, values, ordered=True, category=’comp’)

A ParameterGroup which specifically computes over sections of an image for an ImageState. To this end, we require the implementation of several new member functions:

In order to facilitate optimizations such as caching and local updates, we must incorporate tiling in this object.

Parameters:
  • params (string, list of strings) – The names of the parameters, in the proper order
  • values (number, list of numbers) – The values corresponding to the parameter names
  • ordered (boolean (default: True)) – If True, uses an OrderedDict so that parameter order is deterministic independent of number of parameters
  • category (string (default: 'param')) – Name of the category associated with this ParameterGroup.
execute(*args, **kwargs)

Perform its routine, whatever that may be

exports()

Which methods a class wants to expose to parent classes

get()

Return the natural part of the model. In the case of most elements it is the calculated field, for others it is the object itself.

get_padding_size(tile)

Get the amount of padding required for this object when calculating about a tile tile. Padding size is the total size, so half that on each side. For example, if this Component is a Gaussian point spread function, then the padding returned might be:

peri.util.Tile(np.ceil(2*self.sigma))
Parameters:tile (Tile) – A tile defining the region of interest
Returns:pad – A tile corresponding to the required padding size
Return type:Tile
get_update_tile(params, values)

This method returns a Tile object defining the region of a field that has to be modified by the update of (params, values). For example, if this Component is the point-spread-function, it might return a tile of entire image since every parameter affects the entire image:

return self.shape
Parameters:
  • params (single param, list of params) – A single parameter or list of parameters to be updated
  • values (single value, list of values) – The values corresponding to the params
Returns:

tile – A tile corresponding to the image region

Return type:

Tile

initialize()

Begin anew and initialize the component

nopickle()

Class attributes which should not be included in a pickle object

register(obj)

Registery a parent object so that communication maybe happen upwards

set_shape(shape, inner)

Set the overall shape of the calculation area. The total shape of that the calculation can possibly occupy, in pixels. The second, inner, is the region of interest within the image.

set_tile(tile)

Set the currently active tile region for the calculation

trigger_parameter_change()

Notify parents of a parameter change

trigger_update(params, values)

Notify parent of a parameter change

class peri.comp.comp.ComponentCollection(comps, field_reduce_func=None, category=’comp’)

Group a number of components into a single coherent object which a single interface for each function. Obvious reductions are performed in places such as get_update_tile which takes the bounding tile of all constituent get_update_tile. The only reduction which is not straight-forward is get, but by default it adds all fields. This class has the same interface as Component itself.

Parameters:
  • comps (list of Component) – The components to group together
  • field_reduce_func (function(list of ndarrays)) –

    Reduction function for get object of collection, what happens when all comps are reduced to a single field (etc). For example, nested point spread functions may want convolution to be the reduce func:

    def reduce_conv(psfs):
        ffts = np.array([fft.fftn(p) for p in psfs])
        return fft.ifftn(np.prod(ffts, axis=0))
    

    Or, if it is just a simple field that needs to be added together:

    def reduce_add(x):
        return reduce(add, x)
    
  • category (string) – Name of the category associated with this ComponentCollection
get()

Combine the fields from all components

initargs()

Return arguments that are passed to init to setup the class again

set_shape(shape, inner)

Set the shape for all components

set_tile(tile)

Set the current working tile for components

setup_passthroughs()

Inherit some functions from the components that we own. In particular, let’s grab all functions that begin with param_ so the super class knows how to get parameter groups. Also, take anything that is listed under Component.exports and rename with the category type, i.e., SphereCollection.add_particle -> Component.obj_add_particle

split_params(params, values=None)

Split params, values into groups that correspond to the ordering in self.comps. For example, given a sphere collection and slab:

[
    (spheres) [pos rad etc] [pos val, rad val, etc]
    (slab) [slab params] [slab vals]
]
sync_params()

Ensure that shared parameters are the same value everywhere

peri.comp.objs

class peri.comp.objs.PlatonicSpheresCollection(pos, rad, shape=None, zscale=1.0, support_pad=4, method=’exact-gaussian-fast’, alpha=None, user_method=None, exact_volume=True, volume_error=1e-05, max_radius_change=0.01, param_prefix=’sph’, grouping=’particle’, category=’obj’, float_precision=<class ‘numpy.float64’>)

A collection of spheres in real-space with positions and radii, drawn not necessarily on a uniform grid (i.e. scale factor associated with z-direction). There are many ways to draw the sphere, currently supported methods can be one of:

[
    'bool', 'lerp', 'logistic', 'triangle', 'constrained-cubic',
    'exact-gaussian', 'exact-gaussian-trim', 'exact-gaussian-fast',
    'user-defined'
]
Parameters:
  • pos (ndarray [N,3]) – Initial positions of the spheres
  • rad (ndarray [N] or float) – Initial radii of the spheres
  • shape (tuple) – Shape of the field over which to draw the platonic spheres
  • zscale (float) – scaling of z-pixels in the platonic image
  • support_pad (int) – how much to pad the boundary of particles when calculating support so that there is not more contribution
  • method (string) – The sphere drawing function to use, see above.
  • alpha (float) – Parameter supplied to sphere drawing function, set to value to override default value
  • user_method (tuple (function, parameters)) – Provide your own sphere function to the drawing method. First element of tuple is function with call signature func(dr, a, *args) where the second element is the *args that are not the distance to edge (dr) or particles radius (a). method must be set to ‘user-defined’.
  • exact_volume (boolean) – whether to iterate effective particle size until exact volume (within volume_error) is achieved
  • volume_error (float) – relative volume error tolerance in iteration steps
  • max_radius_change (float) – maximum relative radius change allowed during iteration (due to edge particles and other confounding factors)
  • grouping (string) – Either ‘particle’ or ‘parameter’ parameter grouping. If ‘particle’ then grouped by xyza,xyza if ‘parameter’ then xyz,xyz,a,a
  • float_precision (numpy float datatype) – One of numpy.float16, numpy.float32, numpy.float64; precision for precomputed arrays. Default is np.float64; make it 16 or 32 to save memory.
add_particle(pos, rad)

Add a particle or list of particles given by a list of positions and radii, both need to be array-like.

Parameters:
  • pos (array-like [N, 3]) – Positions of all new particles
  • rad (array-like [N]) – Corresponding radii of new particles
Returns:

inds – Indices of the added particles.

Return type:

N-element numpy.ndarray.

param_particle(ind)

Get position and radius of one or more particles

param_particle_pos(ind)

Get position of one or more particles

param_particle_rad(ind)

Get radius of one or more particles

param_radii()

Return params of all radii

remove_particle(inds)

Remove the particle at index inds, may be a list. Returns [3,N], [N] element numpy.ndarray of pos, rad.

update(params, values)

Calls an update, but clips radii to be > 0

class peri.comp.objs.Slab(zpos=0, angles=(0, 0), param_prefix=’slab’, shape=None, float_precision=<class ‘numpy.float64’>, category=’obj’)

A half plane corresponding to a cover-slip.

Parameters:
  • shape (tuple) – field shape over which to calculate
  • zpos (float) – position of the center of the slab in pixels
  • angles (tuple of float (2,), optional) – Euler-like Angles of rotation of the normal with respect to the z-axis, i.e. angles=(0., 0.) gives a slab with a normal along z. The first angle theta is the rotation about the x-axis; the second angle phi is the rotation about the y-axis. Default is (0,0).
  • float_precision (numpy float datatype) – One of numpy.float16, numpy.float32, numpy.float64; precision for precomputed arrays. Default is np.float64; make it 16 or 32 to save memory.
rmatrix()

Generate the composite rotation matrix that rotates the slab normal.

The rotation is a rotation about the x-axis, followed by a rotation about the z-axis.

peri.comp.ilms

class peri.comp.ilms.Polynomial3D(order=(1, 1, 1), tileinfo=None, constval=None, category=’ilm’, shape=None, float_precision=<class ‘numpy.float64’>)

A polynomial 3D class for updating large fields of polys.

Parameters:
  • shape (peri.util.Tile) – shape of the field (z,y,x)
  • order (tuple) – number of terms in each direction
  • tileinfo (tuple of 2 peri.util.Tile) – These objects help in the transfer of fields from different sections of the same image to new fields. tileinfo is a tuple containing the Tile representing the entire image as well as the Tile representing this particular section of field. (typically given by peri.rawimage.tile)
  • constval (float) – The initial value of the entire field, if a constant.
  • float_precision (numpy float datatype) – One of numpy.float16, numpy.float32, numpy.float64; precision for precomputed arrays. Default is np.float64; make it 16 or 32 to save memory.
class peri.comp.ilms.Polynomial2P1D(order=(1, 1, 1), tileinfo=None, constval=None, operation=’*’, category=’ilm’, shape=None, float_precision=<class ‘numpy.float64’>)

A polynomial 2+1D class for updating large fields of polys. The form of these polynomials if P(x,y) () Q(z), separated in the z-direction.

Parameters:
  • shape (tuple) – shape of the field (z,y,x)
  • order (tuple) – number of terms in each direction
  • tileinfo (tuple of 2 peri.util.Tile) – These objects help in the transfer of fields from different sections of the same image to new fields. tileinfo is a tuple containing the Tile representing the entire image as well as the Tile representing this particular section of field. (typically given by peri.rawimage.tile)
  • constval (float) – The initial value of the entire field, if a constant.
  • operation (string) – Type of joining operation between the (x,y) and (z) poly. Can be either ‘*’ or ‘+’
  • float_precision (numpy float datatype) – One of numpy.float16, numpy.float32, numpy.float64; precision for precomputed arrays. Default is np.float64; make it 16 or 32 to save memory.
class peri.comp.ilms.LegendrePoly2P1D(order=(1, 1, 1), **kwargs)
class peri.comp.ilms.BarnesStreakLegPoly2P1D(npts=(40, 20), zorder=7, op=’*’, barnes_dist=1.75, barnes_clip_size=3, local_updates=True, category=’ilm’, shape=None, float_precision=<class ‘numpy.float64’>, donorm=True)

A Barnes interpolant. This one is of the form

\[I = \left[1 + \left(\sum b_k(x) (o) L_k(y)\right)\right] (1 + z q(z)) + c\]

where b_k are independent barnes interpolants and L_k are legendre polynomials. q is a polynomial strictly in z. Additionally, the operation (o) is settable.

Parameters:
  • shape (iterable) – size of the field in pixels, needs to be padded shape
  • npts (tuple of ints, optional) – Number of control points used for the Barnes interpolant b_k in the x-y sum. Default is (40,20)
  • zorder (integer) – Number of orders for the z-polynomial.
  • op (string) – The operation to perform between Barnes and LegPoly, ‘*’ or ‘+’.
  • barnes_dist (float) – Fractional distance to use for the barnes interpolator
  • local_updates (boolean) – Whether to perform local updates on the ILM
  • float_precision (numpy float datatype) – One of numpy.float16, numpy.float32, numpy.float64; precision for precomputed arrays. Default is np.float64; make it 16 or 32 to save memory.
  • donorm (Bool) – Whether or not to normalize the Barnes interpolation (compatibility patch). Use True, i.e. normalize the Barnes interpolant. Old version is False. Default is True.
randomize_parameters(ptp=0.2, fourier=False, vmin=None, vmax=None)

Create random parameters for this ILM that mimic experiments as closely as possible without real assumptions.

peri.comp.psfs

class peri.comp.psfs.PSF(params, values, shape=None)

Point spread function classes must contain the following classes in order to interface with the states class:

get_padding_size(z) : get the psf size at certain z position get_params : return the parameters of the psf set_tile : set the current update tile size update : update the psf based on new parameters execute : apply the psf to an image
class peri.comp.psfs.IdentityPSF

Delta-function PSF; returns the field passed to execute identically. Params is an N-element numpy.ndarray, doesn’t do anything.

class peri.comp.psfs.AnisotropicGaussian(sigmas=(2.0, 1.0), error=0.00392156862745098, shape=None)
class peri.comp.psfs.AnisotropicGaussianXYZ(sigmas=(2.0, 0.5, 1.0), error=0.00392156862745098, shape=None)
class peri.comp.psfs.PSF4D(params, values, shape=None)

4-dimensional Point-Spread-Function (PSF) is implemented by assuming that the there is only z-dependence of parameters (so that it can be separated into x-y and z parts). Therefore, we keep track of 2D FFTs and X-Y psf functions that get convolved with FFTs. Then, we manually convolve in the z direction

The key variables are rpsf (2d) and kpsf (2d) which are used for the x-y convolution. The z-convolution cannot be cached.

rpsf_xy(vecs, z)

Returns the x-y plane real space psf function as a function of z values. This function does not necessarily have to be normalized, it will be normalized in k-space layer by layer later.

rpsf_z(z, zp)

Returns the z dependence of the PSF. This section needs to be noramlized.

class peri.comp.psfs.Gaussian4D(sigmas=(2.0, 0.5, 1.0), order=(1, 1, 1), error=0.00392156862745098, zrange=128, shape=None)
class peri.comp.psfs.FromArray(array, *args, **kwargs)

Only thing to pass is the values of the point spread function (does not need to be normalized) in the form of a numpy ndarray of shape

(z’, z, y, x)

so if there are 50 layers in the image and the psf is at most 16 wide then it must be shaped (50,16,16,16). The values of the psf must be centered in this array. Hint: np.fft.fftfreq provides the correct ordering of values for both even and odd lattices.

peri.comp.exactpsf

class peri.comp.exactpsf.ExactPSF(shape=None, zrange=None, laser_wavelength=0.488, zslab=0.0, zscale=1.0, kfki=0.889, n2n1=0.9486166007905138, alpha=1.173, polar_angle=0.0, pxsize=0.125, support_factor=2, normalize=False, sigkf=0.0, nkpts=None, cutoffval=None, measurement_iterations=None, k_dist=’gaussian’, use_J1=True, sph6_ab=None, global_zscale=False, cutbyval=False, cutfallrate=0.25, cutedgeval=1e-12, pinhole_width=None, do_pinhole=False, *args, **kwargs)

Superclass for all the exact PSFs, i.e. any PSF that is based on physical properties of the imaging system such as the laser wavelength.

This PSF functions by calculating the local PSF for every z layer in the image, and convolving each layer independently (numerically the exact model of image formation).

Parameters:
  • shape (tuple) – Shape of the image in (z,y,x) pixel numbers (to be deprecated)
  • zrange (tuple) – range of z pixels over which we should calculate the psf, good pixels being zrange[0] <= z <= zrange[1]. currently must be set to the interior of the image, so [state.pad, state.image.shape[0] - state.pad]
  • laser_wavelength (float) – wavelength of light in um of the incoming laser light
  • zslab (float) – Pixel position of the optical interface where 0 is the edge of the image in the z direction
  • zscale (float) – Scaling of the z pixel size, so that each pixel is located at zscale*(z - zint), such that the interface does not move with zscale
  • kfki (float) – Ratio of outgoing to incoming light wavevectors, 2pi/lambda
  • n2n1 (float) – Ratio of the index mismatch across the optical interface. For typical glass with glycerol/water 80/20 mixture this is 1.4/1.518
  • alpha (float) – Aperture of the lens in radians, set by arcsin(n2n1)?
  • polar_angle (float) – the angle of the light polarization with respect to the scanning axis
  • pxsize (float) – the size of a xy pixel in um, defaults to cohen group size 0.125 um
  • support_factor (integer) – size of the support
  • normalize (boolean) – if True, use the normalization as calculated by the PSF instead of unit normalization
  • sigkf (float) – Width of wavelengths to use in the polychromatic psf, None is monochromatic. Values for kfki = kfki +- sigkf, unitless
  • nkpts (integer) – number of integration points to use for the polychromatic psf
  • cutoffval (float) – Percentage of peak psf value to cutoff using a 3-axis exp(-(r-r0)**4) function where r0 is determined by cutoffval. A good value for this is the bit depth of the camera, or possibly the SNR, so 1/2**8 or 1/50.
  • measurement_iterations (int) – number of interations used when trying to find the center of mass of the psf in a certain slice
  • k_dist (str) – Eithe [‘gaussian’, ‘gamma’] which control the wavevector distribution for the polychromatic detection psf. Default is Gaussian.
  • use_J1 (boolean) – Which hdet confocal model to use. Set to True to include the J1 term corresponding to a large-NA focusing lens, False to exclude it. Default is True
  • cutbyval (boolean) – If True, cuts the PSF based on the actual value instead of the position associated with the nearest value.
  • cutfallrate (float) – The relative value of the cutoffval over which to damp the remaining values of the psf. 0.3 seems a good fit now.
  • cutedgeval (float) – The value with which to determine the edge of the psf, typically taken around floating point, 1e-12
  • pinhole_width (Float) – The width of the line illumination, in 1/k units. Default is 1.0.
  • do_pinhole (Bool) – Whether or not to include pinhole line width in the sampling. Default is False.
characterize_psf()

Get support size and drift polynomial for current set of params

drift(z)

Give the pixel offset at a given z value for the current parameters

measure_size_drift(z, size=31, zoffset=0.0)

Returns the ‘size’ of the psf in each direction a particular z (px)

pack_args()

Pack the parameters into the form necessary for the integration routines above in psfcalc.

psf_slice(zint, size=11, zoffset=0.0, getextent=False)

Calculates the 3D psf at a particular z pixel height

Parameters:
  • zint (float) – z pixel height in image coordinates , converted to 1/k by the function using the slab position as well
  • size (int, list, tuple) – The size over which to calculate the psf, can be 1 or 3 elements for the different axes in image pixel coordinates
  • zoffset (float) – Offset in pixel units to use in the calculation of the psf
  • cutval (float) – If not None, the psf will be cut along a curve corresponding to p(r) == 0 with exponential damping exp(-d^4)
  • getextent (boolean) – If True, also return the extent of the psf in pixels for example to get the support size. Can only be used with cutval.
psffunc(*args, **kwargs)

The function to evaluate the exact psf. Syntax must be:

psf = psffunc(x, y, z, **kwargs)

and return a [x.size, y.size, z.size] numpy.ndarray, where x,y,z are 1D arrays. Implement in subclass.

class peri.comp.exactpsf.ExactLineScanConfocalPSF(shape=None, zrange=None, laser_wavelength=0.488, zslab=0.0, zscale=1.0, kfki=0.889, n2n1=0.9486166007905138, alpha=1.173, polar_angle=0.0, pxsize=0.125, support_factor=2, normalize=False, sigkf=0.0, nkpts=None, cutoffval=None, measurement_iterations=None, k_dist=’gaussian’, use_J1=True, sph6_ab=None, global_zscale=False, cutbyval=False, cutfallrate=0.25, cutedgeval=1e-12, pinhole_width=None, do_pinhole=False, *args, **kwargs)

PSF for line-scanning confocal microscopes that can be used with the peri framework. Calculates the spatially varying point spread function for confocal microscopes and allows them to be applied to images as a convolution.

This PSF assumes that the z extent is large compared to the image size and so calculates the local PSF for every z layer in the image.

Parameters:
  • shape (tuple) – Shape of the image in (z,y,x) pixel numbers (to be deprecated)
  • zrange (tuple) – range of z pixels over which we should calculate the psf, good pixels being zrange[0] <= z <= zrange[1]. currently must be set to the interior of the image, so [state.pad, state.image.shape[0] - state.pad]
  • laser_wavelength (float) – wavelength of light in um of the incoming laser light
  • zslab (float) – Pixel position of the optical interface where 0 is the edge of the image in the z direction
  • zscale (float) – Scaling of the z pixel size, so that each pixel is located at zscale*(z - zint), such that the interface does not move with zscale
  • kfki (float) – Ratio of outgoing to incoming light wavevectors, 2pi/lambda
  • n2n1 (float) – Ratio of the index mismatch across the optical interface. For typical glass with glycerol/water 80/20 mixture this is 1.4/1.518
  • alpha (float) – Aperture of the lens in radians, set by arcsin(n2n1)?
  • polar_angle (float) – the angle of the light polarization with respect to the scanning axis
  • pxsize (float) – the size of a xy pixel in um, defaults to cohen group size 0.125 um
  • support_factor (integer) – size of the support
  • normalize (boolean) – if True, use the normalization as calculated by the PSF instead of unit normalization
  • sigkf (float) – Width of wavelengths to use in the polychromatic psf, None is monochromatic. Values for kfki = kfki +- sigkf, unitless
  • nkpts (integer) – number of integration points to use for the polychromatic psf
  • cutoffval (float) – Percentage of peak psf value to cutoff using a 3-axis exp(-(r-r0)**4) function where r0 is determined by cutoffval. A good value for this is the bit depth of the camera, or possibly the SNR, so 1/2**8 or 1/50.
  • measurement_iterations (int) – number of interations used when trying to find the center of mass of the psf in a certain slice
  • k_dist (str) – Eithe [‘gaussian’, ‘gamma’] which control the wavevector distribution for the polychromatic detection psf. Default is Gaussian.
  • use_J1 (boolean) – Which hdet confocal model to use. Set to True to include the J1 term corresponding to a large-NA focusing lens, False to exclude it. Default is True
  • cutbyval (boolean) – If True, cuts the PSF based on the actual value instead of the position associated with the nearest value.
  • cutfallrate (float) – The relative value of the cutoffval over which to damp the remaining values of the psf. 0.3 seems a good fit now.
  • cutedgeval (float) – The value with which to determine the edge of the psf, typically taken around floating point, 1e-12
  • pinhole_width (Float) – The width of the line illumination, in 1/k units. Default is 1.0.
  • do_pinhole (Bool) – Whether or not to include pinhole line width in the sampling. Default is False.
  • Notes – a = ExactLineScanConfocalPSF((64,)*3) psf, (z,y,x) = a.psf_slice(1., size=51) imshow((psf*r**4)[:,:,25], cmap=’bone’)
pack_args()

Pack the parameters into the form necessary for the integration routines above. For example, packs for calculate_linescan_psf

psffunc(*args, **kwargs)

Calculates a linescan psf

class peri.comp.exactpsf.ExactPinholeConfocalPSF(shape=None, zrange=None, laser_wavelength=0.488, zslab=0.0, zscale=1.0, kfki=0.889, n2n1=0.9486166007905138, alpha=1.173, polar_angle=0.0, pxsize=0.125, support_factor=2, normalize=False, sigkf=0.0, nkpts=None, cutoffval=None, measurement_iterations=None, k_dist=’gaussian’, use_J1=True, sph6_ab=None, global_zscale=False, cutbyval=False, cutfallrate=0.25, cutedgeval=1e-12, *args, **kwargs)

PSF for a pinhole confocal microscopes that can be used with the peri framework. Calculates the spatially varying point spread function for confocal microscopes and allows them to be applied to images as a convolution.

This PSF assumes that the z extent is large compared to the image size and so calculates the local PSF for every z layer in the image.

Parameters:
  • shape (tuple) – Shape of the image in (z,y,x) pixel numbers (to be deprecated)
  • zrange (tuple) – range of z pixels over which we should calculate the psf, good pixels being zrange[0] <= z <= zrange[1]. currently must be set to the interior of the image, so [state.pad, state.image.shape[0] - state.pad]
  • laser_wavelength (float) – wavelength of light in um of the incoming laser light
  • zslab (float) – Pixel position of the optical interface where 0 is the edge of the image in the z direction
  • zscale (float) – Scaling of the z pixel size, so that each pixel is located at zscale*(z - zint), such that the interface does not move with zscale
  • kfki (float) – Ratio of outgoing to incoming light wavevectors, 2pi/lambda
  • n2n1 (float) – Ratio of the index mismatch across the optical interface. For typical glass with glycerol/water 80/20 mixture this is 1.4/1.518
  • alpha (float) – Aperture of the lens in radians, set by arcsin(n2n1)?
  • polar_angle (float) – the angle of the light polarization with respect to the scanning axis
  • pxsize (float) – the size of a xy pixel in um, defaults to cohen group size 0.125 um
  • support_factor (integer) – size of the support
  • normalize (boolean) – if True, use the normalization as calculated by the PSF instead of unit normalization
  • sigkf (float) – Width of wavelengths to use in the polychromatic psf, None is monochromatic. Values for kfki = kfki +- sigkf, unitless
  • nkpts (integer) – number of integration points to use for the polychromatic psf
  • cutoffval (float) – Percentage of peak psf value to cutoff using a 3-axis exp(-(r-r0)**4) function where r0 is determined by cutoffval. A good value for this is the bit depth of the camera, or possibly the SNR, so 1/2**8 or 1/50.
  • measurement_iterations (int) – number of interations used when trying to find the center of mass of the psf in a certain slice
  • k_dist (str) – Eithe [‘gaussian’, ‘gamma’] which control the wavevector distribution for the polychromatic detection psf. Default is Gaussian.
  • use_J1 (boolean) – Which hdet confocal model to use. Set to True to include the J1 term corresponding to a large-NA focusing lens, False to exclude it. Default is True
  • cutbyval (boolean) – If True, cuts the PSF based on the actual value instead of the position associated with the nearest value.
  • cutfallrate (float) – The relative value of the cutoffval over which to damp the remaining values of the psf. 0.3 seems a good fit now.
  • cutedgeval (float) – The value with which to determine the edge of the psf, typically taken around floating point, 1e-12
  • Notes – a = ExactLineScanConfocalPSF((64,)*3) psf, (z,y,x) = a.psf_slice(1., size=51) imshow((psf*r**4)[:,:,25], cmap=’bone’)
pack_args()

Pack the parameters into the form necessary for the integration routines above. For example, packs for calculate_linescan_psf

psffunc(x, y, z, **kwargs)

Calculates a pinhole psf

class peri.comp.exactpsf.ChebyshevPSF(cheb_degree=6, cheb_evals=8, *args, **kwargs)

Same as ExactPSF, except that the convolution is performed in the 4th dimension by employing fast Chebyshev approximates to how the PSF varies with depth into the sample. For help, see ExactPSF.

Other Parameters:
 
  • cheb_degree (integer) – degree of the Chebyshev approximant
  • cheb_evals (integer) – number of interpolation points used to create the coefficient matrix
class peri.comp.exactpsf.FixedSSChebPSF(support_size=[35, 17, 25], *args, **kwargs)

ChebyshevPSF with a fixed support size

characterize_psf()

Get support size and drift polynomial for current set of params