List of Functions

Parameters

see Workflow

PointSpreadFunction generation

see Workflow

Pupils

PointSpreadFunctions.pupil_xyzFunction
pupil_xyz(sz, pp, sampling=nothing)

creates a pupil with electric field distributions in XYZ. Returns a 4D dataset with the electric field components along the 4th dimension. #Arguments

  • 'sz': size of the pupil in pixels
  • 'pp': the PSFParam structure with all the PSF parameters
  • 'sampling': the pixel sampling in the same units as the wavelength
  • 'is_proj': defines whether the pupil is to be interpreted as the projection of the McCutchen pupil or not. This yields a by 1 ./cos(Theta) modified aplanatic factor.
source
PointSpreadFunctions.field_pupilFunction
field_pupil(sz, pp, sampling)

returns the pupil polarization as a 4D array with the XY polarization components stacked along the 4th dimension.

source
PointSpreadFunctions.apply_propagatorsMethod
apply_propagators(pupil, z_planes, pp::PSFParams; sampling=nothing)

propagates a given pupil by a number of z_planes (almost) symmetrically in both directions. The result is a stack of propagated pupils. The slice-to-slice propagator is obtained via the get_propagator method.

See also:

  • get_propagator
source
PointSpreadFunctions.get_zernike_pupil_phaseMethod
get_zernike_pupil_phase(sz, pp, sampling)

calculates the phases in the pupil for a given set of aberrations as defined by J and coefficients. By default this follows the OSA nomenclature. See the help file of ZernikePolynomials.jl for more information. The pupil phase (up to the pupil border as defined by the NA in pp) is returned.

Arguments:

  • sz: size of the real-space array
  • pp: PSF parameter structure
  • sampling: pixelpitch in real space as NTuple

Example:

julia> using PointSpreadFunctions, FFTW

julia> aberr = PointSpreadFunctions.Aberrations([Zernike_Spherical, Zernike_ObliqueAstigmatism],[0.1, 0.2])
Aberrations([12, 3], [0.1, 0.2], :OSA)

julia> pp = PSFParams(580.0, 1.4, 1.518; aberrations=aberr)
PSFParams(580.0, 1.4, 1.518, Float32, ModeWidefield, PointSpreadFunctions.pol_scalar, PointSpreadFunctions.var"#42#43"(), PointSpreadFunctions.MethodPropagateIterative, nothing, Aberrations([12, 3], [0.1, 0.2], :OSA), nothing)

julia> sz = (10,10,64)
(10, 10, 64)

julia> sampling=(190,190,100)
(190, 190, 100)

julia> PointSpreadFunctions.get_zernike_pupil_phase(sz,pp,sampling)
10×10 Matrix{Float64}:
 0.0   0.0         0.0         0.0         0.0         0.0         0.0         0.0         0.0         0.0
 0.0   0.0         0.0         0.533599    0.202003   -0.0206203  -0.170662   -0.21173     0.0         0.0
 0.0   0.0         0.477274    0.186398    0.0287554  -0.104828   -0.250743   -0.3726     -0.361221    0.0
 0.0   0.533599    0.186398    0.0937363   0.0736565   0.016983   -0.112676   -0.278928   -0.3726     -0.21173
 0.0   0.202003    0.0287554   0.0736565   0.154747    0.162853    0.0615812  -0.112676   -0.250743   -0.170662
 0.0  -0.0206203  -0.104828    0.016983    0.162853    0.223607    0.162853    0.016983   -0.104828   -0.0206203
 0.0  -0.170662   -0.250743   -0.112676    0.0615812   0.162853    0.154747    0.0736565   0.0287554   0.202003
 0.0  -0.21173    -0.3726     -0.278928   -0.112676    0.016983    0.0736565   0.0937363   0.186398    0.533599
 0.0   0.0        -0.361221   -0.3726     -0.250743   -0.104828    0.0287554   0.186398    0.477274    0.0
 0.0   0.0         0.0        -0.21173    -0.170662   -0.0206203   0.202003    0.533599    0.0         0.0
source
PointSpreadFunctions.get_zernike_pupilMethod
get_zernike_pupil(sz, pp, sampling)

calculates the phases in the pupil for a given set of aberrations as defined by J and coefficients. By default this follows the OSA nomenclature. See the help file of ZernikePolynomials.jl for more information. The complex-valued pupil (up to the pupil border as defined by the NA in pp) is returned.

Arguments:

  • sz: size of the real-space array
  • pp: PSF parameter structure
  • sampling: pixelpitch in real space as NTuple

Example:

julia> using PointSpreadFunctions, FFTW

julia> aberr = PointSpreadFunctions.Aberrations([Zernike_Spherical, Zernike_ObliqueAstigmatism],[0.1, 0.2])
Aberrations([12, 3], [0.1, 0.2], :OSA)

julia> pp = PSFParams(580.0, 1.4, 1.518; aberrations=aberr)
PSFParams(580.0, 1.4, 1.518, Float32, ModeWidefield, PointSpreadFunctions.pol_scalar, PointSpreadFunctions.var"#42#43"(), PointSpreadFunctions.MethodPropagateIterative, nothing, Aberrations([12, 3], [0.1, 0.2], :OSA), nothing)

julia> sz = (10,10,64)
(10, 10, 64)

julia> sampling=(190,190,100)
(190, 190, 100)

julia> PointSpreadFunctions.get_zernike_pupil(sz,pp,sampling)
10×10 Matrix{ComplexF64}:
 1.0+0.0im        1.0+0.0im               1.0+0.0im             1.0+0.0im               1.0+0.0im       …          1.0+0.0im             1.0+0.0im               1.0+0.0im             1.0+0.0im
 1.0+0.0im        1.0+0.0im               1.0+0.0im       -0.977799-0.209546im     0.297025+0.95487im         0.478105-0.878303im   0.238146-0.971229im          1.0+0.0im             1.0+0.0im
 1.0+0.0im        1.0+0.0im         -0.989823+0.142305im   0.389073+0.921207im     0.983723+0.179694im     -0.00466964-0.999989im  -0.696362-0.717691im    -0.643318-0.765599im        1.0+0.0im
 1.0+0.0im  -0.977799-0.209546im     0.389073+0.921207im   0.831518+0.555499im     0.894807+0.446453im        0.759688-0.650288im  -0.180764-0.983527im    -0.696362-0.717691im   0.238146-0.971229im
 1.0+0.0im   0.297025+0.95487im      0.983723+0.179694im   0.894807+0.446453im     0.563395+0.826187im        0.926073+0.377344im   0.759688-0.650288im  -0.00466964-0.999989im   0.478105-0.878303im
 1.0+0.0im   0.991619-0.129199im     0.790818-0.612051im   0.994312+0.106505im     0.520607+0.853797im  …     0.520607+0.853797im   0.994312+0.106505im     0.790818-0.612051im   0.991619-0.129199im
 1.0+0.0im   0.478105-0.878303im  -0.00466964-0.999989im   0.759688-0.650288im     0.926073+0.377344im        0.563395+0.826187im   0.894807+0.446453im     0.983723+0.179694im   0.297025+0.95487im
 1.0+0.0im   0.238146-0.971229im    -0.696362-0.717691im  -0.180764-0.983527im     0.759688-0.650288im        0.894807+0.446453im   0.831518+0.555499im     0.389073+0.921207im  -0.977799-0.209546im
 1.0+0.0im        1.0+0.0im         -0.643318-0.765599im  -0.696362-0.717691im  -0.00466964-0.999989im        0.983723+0.179694im   0.389073+0.921207im    -0.989823+0.142305im        1.0+0.0im
 1.0+0.0im        1.0+0.0im               1.0+0.0im        0.238146-0.971229im     0.478105-0.878303im        0.297025+0.95487im   -0.977799-0.209546im          1.0+0.0im             1.0+0.0im
source
PointSpreadFunctions.k_pupilMethod
k_pupil(pp::PSFParams)

maxim radial (i.e. XY-) k-coordinate (1/space units) where the pupil ends

Arguments:

  • pp: PSF parameter structure
source
PointSpreadFunctions.k_dzMethod
k_dz(pp::PSFParams)

relative kz range from pupil boarder to top of Ewald sphere (in 1/space units)

Arguments:

  • pp: PSF parameter structure
source
PointSpreadFunctions.k_scaleMethod
k_scale(sz,pp,sampling)
pixelpitch (as NTuple) in k-space
  • sz: size of the real-space array
  • pp: PSF parameter structure
  • sampling: pixelpitch in real space as NTuple
source
PointSpreadFunctions.k_pupil_posMethod
k_pupil_pos(sz, pp::PSFParams, sampling)

returns the X and Y position of the pupil border in reciprocal space pixels.

Arguments:

  • sz: size of the real-space array
  • pp: PSF parameter structure
  • sampling: pixelpitch in real space as NTuple
source
PointSpreadFunctions.k_0_posMethod
k_0_pos(sz, pp::PSFParams, sampling)

returns the X, Y and Z position of the Ewald-sphere border in reciprocal space pixels.

Arguments:

  • sz: size of the real-space array
  • pp: PSF parameter structure
  • sampling: pixelpitch in real space as NTuple
source

Polarization

These functions can be conveniently supplied to PSFParams() via the named argument polarization.

PointSpreadFunctions.pol_scalarFunction
pol_scalar(T, xypos)

ignores polarization aspects in the calculation but calculates only based on (high-NA) scalar theory. This is a lot faster but not as accurate.

source
PointSpreadFunctions.pol_scalar_spiralFunction
pol_scalar_spiral(T, xypos)

ignores polarization aspects in the calculation but calculates only based on (high-NA) scalar theory. This version still includes a (scalar) phase spiral. This is a lot faster but not as accurate.

source
PointSpreadFunctions.pol_xFunction
pol_x(T, xypos)

assumes x-polarization in illumination or an x-oriented polarizer in detection. In a high-NA objective this is converted into XYZ electric fields at the focus.

source
PointSpreadFunctions.pol_yFunction
pol_y(T, xypos)

assumes y-polarization in illumination or an x-oriented polarizer in detection. In a high-NA objective this is converted into XYZ electric fields at the focus.

source
PointSpreadFunctions.pol_circFunction
pol_circ(T, xypos)

assumes circular polarization in illumination or an x-oriented polarizer in detection. In a high-NA objective this is converted into XYZ electric fields at the focus.

source
PointSpreadFunctions.pol_circ_spiralFunction
pol_circ_spiral(T, xypos)

assumes circular polarization in illumination or an x-oriented polarizer in detection. This version includes phase spiral defining the local (xypos-dependent) phase of both x and y polarization.

source
PointSpreadFunctions.pol_circ_tophatFunction
pol_circ_tophat(T, xypos)

assumes circular polarization in illumination. This version includes phase 0/π tophat defining the local (xypos-dependent) phase of both x and y polarization.

source
PointSpreadFunctions.pol_radialFunction
pol_radial(T, xypos)

assumes radial polarization in (illumination/detection) of at the pupil.

Example:

julia> using PointSpreadFunctions, View5D
julia> pp_em = PSFParams(0.532, 1.3, 1.52; mode=ModeWidefield, pol=pol_radial);
julia> h_p = apsf(MethodPropagate, sz, pp_em, sampling=samp);
julia> @vv real.(h_p[:,:,1,3])
source
PointSpreadFunctions.pol_radial_annulusFunction
pol_radial_annulus(T, xypos; r0= 0.8, σ=0.05) # e.g. for STED microscopy

returns a pupil function with radial polarization and an annulus at relative pupil radius r0 of exponential apodization width σ.

source
PointSpreadFunctions.pupil_annulusFunction
pupil_annulus(T, xypos; r0= 0.8, σ=0.05)

returns a pupil function with an annulus at relative pupil radius r0 of exponential apodization width σ. Note that this pupil function should not be used by itself, but needs to be combined with a polarization function.

source
PointSpreadFunctions.pupil_apodize_hannFunction
pupil_apodize_hann(T, xypos; r0= 0.9) #

returns a pupil function with a hanning-type apodization starting at radius r0, extending to the pupil limit. Both, cos-window and Hann window apodization should yield a finite standard deviation of the PSF. Note that this pupil function should not be used by itself, but needs to be combined with a polarization function.

source
PointSpreadFunctions.pupil_apodize_cosFunction
pupil_apodize_cos(T, xypos; r0= 0.9) #

returns a pupil function with a cos-type apodization starting at radius r0, extending to the pupil limit. The cos-type apodization reaches zero at the border of the pupil, but with a slope. This cos-filter goes back to Gabor. Both, cos-window and Hann window apodization should yield a finite standard deviation of the PSF. Note that this pupil function should not be used by itself, but needs to be combined with a polarization function.

source

Aplanatic factors

PointSpreadFunctions.aplanatic_illuminationFunction
aplanatic_illumination = (θ) -> sqrt.(max.(0,cos.(θ)))

This is the aplanatic factor typically used in illumination of (randomly oriented) fluorophores. Note that it is identical to detection.

source

Sampling

PointSpreadFunctions.get_Abbe_limitMethod
get_Abbe_limit(pp::PSFParams)

returns the Abbe limit of incoherent imaging in real space as a Tuple with 3 entries. Note that the coherent limit needs only a factor of two less sampling, as long as no intensity is calculated. This allows upsampling right before calculating intensities.

See also:

  • getrequiredamp_sampling()

Example:

julia> using PointSpreadFunctions

julia> pp = PSFParams(580.0, 1.4, 1.518)
PSFParams(580.0, 1.4, 1.518, Float32, ModeWidefield, PointSpreadFunctions.pol_scalar, PointSpreadFunctions.var"#42#43"(), PointSpreadFunctions.MethodPropagateIterative, nothing, Aberrations(Any[], Any[], :OSA), nothing)

julia> PointSpreadFunctions.get_Abbe_limit(pp)
(191.04085f0, 191.04085f0, 622.8464f0)
source
PointSpreadFunctions.get_required_amp_samplingMethod
get_required_amp_sampling(sz::NTuple, pp::PSFParams)

returns the necessary sampling (in real space) for sampling amplitudes. This is almost corresponding to the Abbe limit. Factor of two less because of amplitude but twice because of the Nyquist theorem. However, this is slight modified by requiring slightly higher sampling (one empty pixel on each side of Fourier space) to stay clear of ambiguities.

See also:

  • getAbbelimit()

Example:

julia> using PointSpreadFunctions

julia> pp = PSFParams(580.0, 1.4, 1.518)
PSFParams(580.0, 1.4, 1.518, Float32, ModeWidefield, PointSpreadFunctions.pol_scalar, PointSpreadFunctions.var"#42#43"(), PointSpreadFunctions.MethodPropagateIterative, nothing, Aberrations(Any[], Any[], :OSA), nothing)

julia> sz = (256,256,64)
(256, 256, 64)

julia> PointSpreadFunctions.get_required_amp_sampling(sz,pp)
(188.05583f0, 188.05583f0, 219.73679f0)
source
PointSpreadFunctions.check_amp_sampling_xyMethod
check_amp_sampling_xy(sz, pp,sampling)

issues a warning if the amplitude sampling along X and Y is not within the required limits.

See also:

  • getAbbelimit()
  • getrequiredamp_sampling()
  • checkampsampling()
  • checkampsampling_z()

Arguments:

  • sz: size of the real-space array
  • pp: PSF parameter structure
  • sampling: pixelpitch in real space as NTuple
source
PointSpreadFunctions.check_amp_sampling_zMethod
check_amp_sampling_z(sz, pp,sampling)

issues a warning if the amplitude sampling along Z is not within the required limits.

See also:

  • getAbbelimit()
  • getrequiredamp_sampling()
  • checkampsampling_xy()
  • checkampsampling()

Arguments:

  • sz: size of the real-space array
  • pp: PSF parameter structure
  • sampling: pixelpitch in real space as NTuple
source
PointSpreadFunctions.check_amp_samplingMethod
check_amp_sampling(sz, pp,sampling)

issues a warning if the amplitude sampling (X,Y and Z) is not within the required limits.

See also:

  • getAbbelimit()
  • getrequiredamp_sampling()
  • checkampsampling_xy()
  • checkampsampling_z()

Arguments:

  • sz: size of the real-space array
  • pp: PSF parameter structure
  • sampling: pixelpitch in real space as NTuple
source
PointSpreadFunctions.check_amp_sampling_sincrMethod
check_amp_sampling_sincr(sz, pp,sampling)

issues a warning if the amplitude sampling (X,Y and Z) is not within the required limits for the SincR method, as this requires sampling the Ewald-sphere according to Nyquist.

See also:

  • getAbbelimit()
  • getrequiredamp_sampling()
  • checkampsampling()
  • getEwaldsampling()

Arguments:

  • sz: size of the real-space array
  • pp: PSF parameter structure
  • sampling: pixelpitch in real space as NTuple
source
PointSpreadFunctions.get_Nyquist_limitFunction
get_Nyquist_limit(pp::PSFParams)

returns the Nyquist sampling limit of incoherent imaging in real space as a Tuple with 3 entries. Note that the coherent limit needs only a factor of two less sampling, as long as no intensity is calculated. This allows upsampling right before calculating intensities.

See also:

  • getAbbelimit()

Example:

julia> using PointSpreadFunctions

julia> pp = PSFParams(580.0, 1.4, 1.518)
PSFParams(580.0, 1.4, 1.518, Float32, ModeWidefield, PointSpreadFunctions.pol_scalar, PointSpreadFunctions.var"#42#43"(), PointSpreadFunctions.MethodPropagateIterative, nothing, Aberrations(Any[], Any[], :OSA), nothing)

julia> PointSpreadFunctions.get_Nyquist_limit(pp)
(95.520424f0, 95.520424f0, 311.4232f0)
source
PointSpreadFunctions.limit_θFunction
limit_θ(theta, α_max)

limits the maximum angle to the theoretical α_max. This is particularly important for the 1/sqrt(cos θ) aplanatic factor as this factor diverges and wins over the roundoff errors of the smooth pupil.

source
PointSpreadFunctions.limit_kzMethod
limit_k_z(ft_shell, pp::PSFParams, sampling)

limits the k_z range of the ewald-sphere. returns the extracted region and the new sampling

source
PointSpreadFunctions.AU_per_pixelFunction
AU_in_pixels(pp, sampling)

calculates the size of the Airy unites in pixels. #arguments

  • pp: PSF parameter structure
  • sampling: sampling of the image pixels in object coordinates
source

Utilities

PointSpreadFunctions.amp_to_intMethod
amp_to_int(field, pp)

converts a complex-valued amplitude field to intensity via abs2. and summing over the 4th dimension and dropping the 4th dimension. If pp.transition_dipole is provided, the dipole transition probability is calculated, otherwise a directionally independed effect is assumed.

#Arguments

  • field: The 4D field to extract the intensity from
  • pp: Is checked for whether pp.transition_dipole is not nothing.
source
PointSpreadFunctions.calc_with_resamplingFunction
calc_with_resampling(fct, sz, sampling, args)
calculates a function `fct` on a twice downsampled grid and upsamples the result.

#Parameters +fct: Function that performs the calculation. Its first argument needs to be the N-dimensional data size sz and the sampling. +sz: size of the array to calculate +sampling: pixel sizes of the final array to calculate

  • norm_amp: decides whether amplitude or intensity is normalized to account for the size change
source
PointSpreadFunctions.normalize_amp_to_planeFunction
normalize_amp_to_plane(apsf, plane=nothing, mydim=3)

normalizes an amplitude PSF such that the intensity sum over the central plane position is one. Dimension mydim defines the coordinate referring the the slicing.

source
PointSpreadFunctions.has_z_symmetryMethod
has_z_symmetry(pp::PSFParams)

checks whether the point spread function is expected to be symmetric along the z-direction. Currently this is defined by the aberration list being empty.

source
PointSpreadFunctions.get_McCutchen_kz_centerMethod
get_McCutchen_kz_center(ft_shell, pp::PSFParams, sampling)

calculates the (rounded) pixels position half way between both, the kz-borders of the McCutchen pupil to extract from the full sized Ewald sphere. The pixel z position is returned together with the corresponding kz position.

source
PointSpreadFunctions.disc_pinhole_ftFunction
disc_pinhole_ft(sz, pp, pinhole)

calculates the Fourier transform of a disc-shaped pinhole

  • pp: PSF parameter structure
  • sampling: sampling of the image pixels in object coordinates
  • pinhole: diameter of the pinhole in pixels
source
PointSpreadFunctions.my_discMethod
my_disc(sz, pp)

creates a disc such that there is no overalp upon wrap-around, when convolving it with itself on this grid. This is the radius where there is no overlap. However, in the cornes of the calculation, the results are not 100% accurate as something is missing.

source
PointSpreadFunctions.box_pinhole_ftFunction
disc_pinhole_ft(sz, pp, pinhole)

calculates the Fourier transform of a box-shaped pinhole

  • pp: PSF parameter structure
  • sampling: sampling of the image pixels in object coordinates
  • pinhole: tuple of side lengths of the pinhole in pixels
source
PointSpreadFunctions.sinc_rMethod
sinc_r(sz::NTuple, pp::PSFParams; sampling=nothing)

calculates the 3-dimensional sinc(abs(position)) scaled such that the Fourier transformation yields the k-sphere. Note that for this to work the sampling needs to be sufficient, which may be problematic especially along the z-direction. This can be checked with the help of the get_Ewald_sampling() function.

See also:

  • getEwaldsampling()
  • jincr2d
source
PointSpreadFunctions.sinc_r_2dFunction
sinc_r_2d(sz::NTuple, side_length=(1.0,1.0), dtype=Float32)

calculates a sinc(abs(position)) function in 2D such that its inverse Fourier transform corresponds to the box-shaped pupil of diameter which is a Tuple of both long axes diameters. dtype specifies the destination type. Note that the result is normalized such that the disk as obtained by an inverse Fourier transformation is filled with a value of 1.0.

See also:

  • sinc_r()

Example

julia> using PointSpreadFunctions, View5D

julia> sz = (100,100);d=20.0; w=jinc_r_2d(sz,d); 

julia> q = sinc_r_2d(sz, d; r_func=PointSpreadFunctions.rr_rfft); # a version in RFFT space

julia> @vt rr(sz) .< d/2.0 real.(ift(w)) fftshift(irfft(q,sz[1]))
source
PointSpreadFunctions.jinc_r_2dMethod
jinc_r_2d(sz::NTuple, pp::PSFParams; sampling=nothing)

calculates a jinc(abs(position)) function in 2D such that its Fourier transform corresponds to the disk-shaped pupil (indcluding the effect ot the numerical aperture). If sampling is not provided, sampling is assumed such that the Nyquist sampling would correspond to sampling the wavelength (i.e. the Ewals sphere) correctly for this amplitude.

See also:

  • sinc_r()
source
PointSpreadFunctions.jinc_r_2dFunction
jinc_r_2d(sz::NTuple, pp::PSFParams; sampling=nothing)

calculates a jinc(abs(position)) function in 2D such that its Fourier transform corresponds to the disk-shaped pupil (indcluding the effect ot the numerical aperture). If sampling is not provided, sampling is assumed such that the Nyquist sampling would correspond to sampling the wavelength (i.e. the Ewals sphere) correctly for this amplitude.

See also:

  • sinc_r()
source
jinc_r_2d(sz::NTuple, diameter=(1.0,1.0), dtype=Float32)

calculates a jinc(abs(position)) function in 2D such that its inverse Fourier transform corresponds to the disk-shaped pupil of diameter which is a Tuple of both long axes diameters. dtype specifies the destination type. Note that the result is normalized such that the disk as obtained by an inverse Fourier transformation is filled with a value of 1.0.

See also:

  • sinc_r()

Example

julia> using PointSpreadFunctions, View5D

julia> sz = (100,100);d=20.0; w=jinc_r_2d(sz,d); 

julia> q = jinc_r_2d(sz, d; r_func=PointSpreadFunctions.rr_rfft); # a version in RFFT space

julia> @vt rr(sz) .< d/2.0 real.(ift(w)) fftshift(irfft(q,sz[1]))
source
PointSpreadFunctions.theta_zMethod
theta_z(sz)

a θ function along the z direction, being one for z positon > 0 and zero elsewhere.

Arguments:

  • sz: size of the array
source
PointSpreadFunctions.xx_rfftFunction
xx_rfft(dtype, sz; scale=1.0)
generates the positions in rfft space based on the distance to the zero coordinate. `sz` corresponds to a size of the full (inverse Fourier-transformed) array.
source
PointSpreadFunctions.yy_rfftFunction
yy_rfft(dtype, sz; scale=1.0)
generates the positions in rfft space based on the distance to the zero coordinate. `sz` corresponds to a size of the full (inverse Fourier-transformed) array.
source
PointSpreadFunctions.rr_rfftFunction
rr_rfft(dtype, sz; scale=1.0)
generates the positions in rfft space based on the distance to the zero coordinate. `sz` corresponds to a size of the full (inverse Fourier-transformed) array.
source
PointSpreadFunctions.confocal_intFunction
confocal_int(psf_ex, psf_em, pinhole=nothing, pinhole_ft=disc_pinhole_ft, sampling=nothing, use_resampling=true, return_amp=nothing, pinhole_positions=[(0.0,0.0)], ex_modifier=modify_ident)

helper function to calculate the confocal detection, once the excitation and emission PSFs and pinhole parameters are defined. Two photon excitation and ISM detection are also included.

source
PointSpreadFunctions.kz_mid_posFunction
kz_mid_pos(sz, pp::PSFParams, sampling)

reciprocal space distance from the center, where the mid-point between top and bottom of the McCutchen pupil is.

Arguments:

  • sz: Size of the dataset
  • pp: PSF parameter structure
source