Reference

This is the reference documentation for the ParticleHolography package.

Index

Types and Functions

ParticleHolography.CuLowPassFilterType
CuLowPassFilter{T<: AbstractFloat}

A struct that holds the data for the low pass filter. This can be multiplied with the Fourier transform of wavefront to get the low pass filtered wavefront after propagation.

Fields

  • data::CuArray{T,2}: The data for the low pass filter.
source
ParticleHolography.CuTransferType
CuTransfer{T<: Complex}

A struct that holds the data for the transfer function.

Fields

  • data::CuArray{T,2}: The data for the transfer function.
source
ParticleHolography.CuTransferSqrtPartType
CuTransferSqrtPart{T<: AbstractFloat}

A struct that holds the data for the square root part of the transfer function.

Fields

  • data::CuArray{T,2}: The data for the square root part of the transfer function.
source
ParticleHolography.LabonteMethod
Labonte(dict1, dict2; Dmax=50.0, α=0.005, R=50.0, Rend=0.1, β=0.9, N=10)

Implementation of the improved Labonté algorithm [14, 15]. Takes dictionaries of detected particle coordinates at two time points as input and returns a graph representing the correspondence between particles at these two time points. The graph contains UUID keys of all particles from both time points as nodes, with directed edges representing particle correspondences from particles in dict1 to particles in dict2. All correspondences can be retrieved using the enum_edge function.

Arguments

  • dict1::Dict{UUID, Vector{Float32}}: Dictionary of detected particle coordinates at the first time point.
  • dict2::Dict{UUID, Vector{Float32}}: Dictionary of detected particle coordinates at the second time point.

Optional keyword arguments

  • Dmax::Float32=50.0: Maximum allowable movement distance.
  • α::Float32=0.005: Learning rate.
  • R::Float32=50.0: Initial value of the inspection neighborhood radius.
  • Rend::Float32=0.1: Final value of the inspection neighborhood radius. The algorithm terminates when the inspection radius becomes smaller than this value through iterations.
  • β::Float32=0.9: Decay rate of the inspection neighborhood radius.
  • N::Int=10: Nodes that are not selected as winners for N consecutive times are removed. Refer to the original paper for details.

Returns

  • MetaGraph: Graph representing the correspondence between detected particles at two time points.
source
ParticleHolography.append_path!Method
append_path!(paths, g::MetaGraph)

Adds correspondences from the MetaGraph g to the trajectory array paths initialized by the enum_edge function. Specifically, if the last element of a Vector{UUID} array path in paths matches the starting point of an edge in g, the endpoint of that edge is appended to paths. If the graph is not injective, path is duplicated.

Arguments

  • paths::Vector{UUID}[]: Array of trajectories.
  • g::MetaGraph: The MetaGraph from which correspondences are to be added.

Returns

  • nothing
source
ParticleHolography.cu_apply_low_pass_filter!Method
cu_apply_low_pass_filter!(holo, lpf)

Apply a low pass filter to the wavefront holo. The low pass filter is applied by multiplying the Fourier transform of the wavefront with the Fourier transform of the low pass filter. The wavefront is then reconstructed by taking the inverse Fourier transform of the filtered Fourier transform.

Arguments

  • holo::CuWavefront: The wavefront to apply the low pass filter.
  • lpf::CuLowPassFilter: The low pass filter to apply to the wavefront.

Returns

  • nothing
source
ParticleHolography.cu_apply_low_pass_filterMethod
cu_apply_low_pass_filter(holo, lpf)

Apply a low pass filter to the wavefront holo. The low pass filter is applied by multiplying the Fourier transform of the wavefront with the Fourier transform of the low pass filter. The wavefront is then reconstructed by taking the inverse Fourier transform of the filtered Fourier transform.

Arguments

  • holo::CuWavefront: The wavefront to apply the low pass filter.
  • lpf::CuLowPassFilter: The low pass filter to apply to the wavefront.

Returns

  • CuWavefront{ComplexF32}: The wavefront after applying the low pass filter.
source
ParticleHolography.cu_connected_component_labelingMethod
cu_connected_component_labeling(input_img)

8-way connected component labeling on binary image based on the article by Playne and Hawick https://ieeexplore.ieee.org/document/8274991 and the implementation by FolkeV https://github.com/FolkeV/CUDA_CCL. It works using the CUDA.jl package and NVIDIA GPUs.

Arguments

  • input_img::CuArray{Float32, 2}: Input binary image.

Returns

  • output_img::CuArray{UInt32, 2}: Output labeled image.
source
ParticleHolography.cu_gabor_wavefrontMethod
cu_gabor_wavefront(holo)

Create a wavefront from single hologram holo. This is for Gabor holography. The wavefront is created by taking the square root of the hologram and casting it to a complex number.

Arguments

  • holo::CuArray{Float32,2}: The hologram to create the wavefront from.

Returns

  • CuWavefront{ComplexF32}: The wavefront created from the hologram. See CuWavefront.
source
ParticleHolography.cu_get_reconst_complex_volMethod
cu_get_reconst_complex_vol(holo, transfer_front, transfer_dz, slices)

Reconstruct the observation volume from the wavefront using the transfer functions transfer_front and transfer_dz and return the complex amplitude volume. transfer_front propagates the wavefront to the front of the volume, and transfer_dz propagates the wavefront between the slices. slices is the number of slices in the volume.

Arguments

  • wavefront::CuArray{ComplexF32,2}: The wavefront to reconstruct. In Gabor's holography, this is the square root of the hologram.
  • transfer_front::CuTransfer{ComplexF32}: The transfer function to propagate the wavefront to the front of the volume. See CuTransfer.
  • transfer_dz::CuTransfer{ComplexF32}: The transfer function to propagate the wavefront between the slices.
  • slices::Int: The number of slices in the volume.

Returns

  • CuArray{ComplexF32,3}: The reconstructed complex amplitude volume.
source
ParticleHolography.cu_get_reconst_volFunction
cu_get_reconst_vol(holo, transfer_front, transfer_dz, slices, return_type)

Reconstruct the observation volume from the wavefront using the transfer functions transfer_front and transfer_dz. transfer_front propagates the wavefront to the front of the volume, and transfer_dz propagates the wavefront between the slices. slices is the number of slices in the volume.

Arguments

  • wavefront::CuArray{ComplexF32,2}: The wavefront to reconstruct. In Gabor's holography, this is the square root of the hologram.
  • transfer_front::CuTransfer{ComplexF32}: The transfer function to propagate the wavefront to the front of the volume. See CuTransfer.
  • transfer_dz::CuTransfer{ComplexF32}: The transfer function to propagate the wavefront between the slices.
  • slices::Int: The number of slices in the volume.
  • return_type::Type: The return type of the reconstructed volume. Default is N0f8.

Returns

  • CuArray{return_type,3}: The reconstructed intensity volume.
source
ParticleHolography.cu_get_reconst_vol_and_xyprojectionFunction
cu_get_reconst_vol_and_xyprojection(wavefront, transfer_front, transfer_dz, slices, return_type)

Reconstruct the observation volume from the wavefront and get the XY projection of the volume using the transfer functions transfer_front and transfer_dz. transfer_front propagates the wavefront to the front of the volume, and transfer_dz propagates the wavefront between the slices. slices is the number of slices in the volume.

Arguments

  • wavefront::CuWavefront{ComplexF32}: The wavefront to reconstruct. In Gabor's holography, this is the square root of the hologram.
  • transfer_front::CuTransfer{ComplexF32}: The transfer function to propagate the wavefront to the front of the volume. See CuTransfer.
  • transfer_dz::CuTransfer{ComplexF32}: The transfer function to propagate the wavefront between the slices.
  • slices::Int: The number of slices in the volume.
  • return_type::Type: The return type of the reconstructed volume. Default is N0f8.

Returns

  • CuArray{return_type,3}: The reconstructed volume.
  • CuArray{Float32,2}: The XY projection of the reconstructed volume.
source
ParticleHolography.cu_get_reconst_xyprojectionMethod
cu_get_reconst_xyprojectin(wavefront, transfer_front, transfer_dz, slices)

Get the XY projection of the reconstructed volume from the wavefront using the transfer functions transfer_front and transfer_dz. transfer_front propagates the wavefront to the front of the volume, and transfer_dz propagates the wavefront between the slices. slices is the number of slices in the volume.

Arguments

  • wavefront::CuWavefront{ComplexF32}: The wavefront to reconstruct. In Gabor's holography, this is the square root of the hologram.
  • transfer_front::CuTransfer{ComplexF32}: The transfer function to propagate the wavefront to the front of the volume. See CuTransfer.
  • transfer_dz::CuTransfer{ComplexF32}: The transfer function to propagate the wavefront between the slices.
  • slices::Int: The number of slices in the volume.

Returns

  • CuArray{Float32,2}: The XY projection of the reconstructed volume.
source
ParticleHolography.cu_phase_retrieval_holoMethod
cu_phase_retrieval_holo(holo1, holo2, transfer, invtransfer, priter, datlen)

Perform the Gerchberg-Saxton algorithm-based phase retrieving on two holograms and return the retrieved wavefront at the z-coordinate point of holo1. The algorithm is repeated priter times. holo1 and holo2 are the holograms (I = |phi|^2) of the object at two different z-coordinates. transfer and invtransfer are the transfer functions for the propagation from holo1 to holo2 and vice versa. datlen is the size of the holograms.

Arguments

  • holo1::CuArray{Float32,2}: The hologram at the z-cordinate of closer to the object.
  • holo2::CuArray{Float32,2}: The hologram at the z-coordinate of further from the object.
  • transfer::CuTransfer{ComplexF32}: The transfer function from holo1 to holo2.
  • invtransfer::CuTransfer{ComplexF32}: The transfer function from holo2 to holo1.
  • priter::Int: The number of iterations to perform the algorithm.
  • datlen::Int: The size of the holograms.

Returns

  • CuWavefront{ComplexF32}: The retrieved wavefront at the z-coordinate of holo1. See CuWavefront.
source
ParticleHolography.cu_rectangle_filterMethod
cu_rectangle_filter(prop_dist::AbstractFloat, wavlen::AbstractFloat, imglen::Int, pixel_picth::AbstractFloat)

Creates a low pass filter with a rectangular window. This can be multiplied with the Fourier transform of wavefront to get the low pass filtered wavefront. See Eq. 13 and 14 in (Fugal, 2009, https://doi.org/10.1088/0957-0233/20/7/075501)

Arguments

  • prop_dist::AbstractFloat: The maximum propagation distance of recorded objects.
  • wavlen::AbstractFloat: The wavelength of the light.
  • imglen::Int: The length of the image.
  • pixel_picth::AbstractFloat: The pixel pitch of the image.

Returns

  • CuLowPassFilter: The low pass filter as a CuLowPassFilter object.
source
ParticleHolography.cu_super_gaussian_filterMethod
cu_super_gaussian_filter(prop_dist::AbstractFloat, wavlen::AbstractFloat, imglen::Int, pixel_picth::AbstractFloat)

Creates a low pass filter with a super Gaussian window. This can be multiplied with the Fourier transform of wavefront to get the low pass filtered wavefront. See Eq. 15 in (Fugal, 2009, https://doi.org/10.1088/0957-0233/20/7/075501)

Arguments

  • prop_dist::AbstractFloat: The maximum propagation distance of recorded objects.
  • wavlen::AbstractFloat: The wavelength of the light.
  • imglen::Int: The length of the image.
  • pixel_picth::AbstractFloat: The pixel pitch of the image.

Returns

  • CuLowPassFilter: The low pass filter as a CuLowPassFilter object.
source
ParticleHolography.cu_transferMethod
cu_transfer(z0, datLen, wavLen, d_sqr)

Create a CuArray of size datLen x datLen with the values of the transfer function for a given propagated distance z0. d_sqr can be obtained with cutransfersqrtarr(datlen, wavlen, dx).

Arguments

  • z0::AbstractFloat: The distance to propagate the wave.
  • datLen::Int: The size of the CuArray.
  • wavLen::AbstractFloat: The wavelength of the light.
  • d_sqr::::CuTransferSqrtPart{Float32}: The square of the distance from the center of the hologram, obtained with cutransfersqrtarr(datlen, wavlen, dx).

Returns

  • CuTransfer{Float32}: The transfer function for the propagation. See CuTransfer.
source
ParticleHolography.cu_transfer_sqrt_arrMethod
cu_transfer_sqrt_arr(datlen, wavlen, dx)

Create a CuArray of size datlen x datlen with the values of the square-root part of the transfer function.

Arguments

  • datlen::Int: The size of the CuArray.
  • wavlen::AbstractFloat: The wavelength of the light.
  • dx::AbstractFloat: The pixel size of the hologram.

Returns

  • CuTransferSqrtPart{Float32}: The square-root part of the transfer function. See CuTransferSqrtPart.
source
ParticleHolography.dictloadMethod
dictload(filename)

Load a particle dictionary from a file in JSON format. The dictionary should have UUID keys and values as Vector{Float32}, which includes the coordinates (and diameters) of the particles.

Arguments

  • filename::String: The path to the file.

Returns

  • Dict: The loaded dictionary.
source
ParticleHolography.dictsaveMethod
dictsave(filename, dict)

Save a particle dictionary to a file in JSON format. The dictionary should have UUID keys and values as Vector{Float32}, which includes the coordinates (and diameters) of the particles.

Arguments

  • filename::String: The path to the file.
  • dict::Dict: The dictionary to save.

Returns

  • nothing
source
ParticleHolography.enum_edgeMethod
enum_edge(g::MetaGraph)

Enumerates all edges in the given MetaGraph g and returns them as an array of arrays of UUIDs. Each inner array contains the labels of the source and destination nodes of an edge.

Arguments

  • g::MetaGraph: The MetaGraph to enumerate edges from.

Returns

  • Vector{UUID}[]: An array of arrays of UUIDs, where each inner array contains the labels of the source and destination nodes of an edge.
source
ParticleHolography.finalize_particle_neighborhoods!Method
finalize_particle_neighborhoods!(particle_neighborhoods)

Finalizes the particle neighborhoods by removing duplicates and particles that are too small or too elongated in the x-y plane. Detailed criteria are as follows:

  • Duplicate bounding boxes
  • Bounding boxes with a length-to-width (x-y) ratio greater than 3 or less than 1/3
  • Bounding boxes with an area less than $\sqrt{10}$ pixels
  • Bounding boxes with a depth of 1

Arguments

  • particle_neighborhoods: The particle neighborhoods to be finalized.

Returns

  • nothing
source
ParticleHolography.find_external_contoursMethod
find_external_contours(image)

Finds non-hole contours in binary images. This function is excuted on the CPU. Equivalent to CVRETREXTERNAL and CVCHAINAPPROX_NONE modes of the findContours() function provided in OpenCV.

Arguments

  • image: The binary image. the image should be a 2D array of 0 and 1.

Returns

  • Vector{Vector{CartesianIndex}}: A vector of contours. Each contour is a vector of CartesianIndex.
source
ParticleHolography.gen_fulldictMethod
gen_fulldict(dicts::Vector{Dict{UUID, Vector{Float32}}})

Generates a dictionary containing the UUIDs of particles for all conceivable frames as keys and their coordinates as values. The input is a vector of dictionaries with UUID keys and values as Vector{Float32}, which includes the coordinates (and diameters) of the particles.

Arguments

  • dicts::Vector{Dict{UUID, Vector{Float32}}}: Vector of dictionaries with UUID keys and values as Vector{Float32}.

Returns

  • Dict{UUID, Vector{Float32}}: Dictionary containing the UUIDs of particles for all conceivable frames as keys and their coordinates as values.
source
ParticleHolography.gen_fulldictMethod
gen_fulldict(filepaths::Vector{String})

Generates a dictionary containing the UUIDs of particles for all conceivable frames as keys and their coordinates as values. The input is a vector of file paths containing JSON files with particle coordinates.

Arguments

  • filepaths::Vector{String}: Vector of file paths containing JSON files with particle coordinates.

Returns

  • Dict{UUID, Vector{Float32}}: Dictionary containing the UUIDs of particles for all conceivable frames as keys and their coordinates as values.
source
ParticleHolography.getErrorVecFunction
getErrorVec(vecMap, coefa, gridSize = 128, imgSize = 1024)

Return the vector `errorVec` obtained by calculating the difference between `targetX` `targetY` and `procX` `procY` obtained by converting `targetX` `targetY` from the `gridx` `gridy` set in the first image to the second image. Also returns `gridx` `gridy` for Jacobian acquisition.

Arguments

  • vecMap: The PIV map.
  • coefa: The coefficients for the quadratic distortion correction.
  • gridSize: The size of the grid.
  • imgSize: The size of the image.
source
ParticleHolography.getYacobianFunction
getYacobian(imgSize = 1024, gridSize = 128)

Get the Jacobian matrix for the quadratic distortion correction. The default image size is 1024, and the grid size is 128.

Arguments

  • imgSize::Int: The size of the image.
  • gridSize::Int: The size of the grid.

Returns

  • Array{Float64,2}: The Jacobian matrix.
source
ParticleHolography.get_distortion_coefficientsMethod
get_distortion_coeficients(img1, img2, gridSize = 128, intrSize = 128, srchSize = 256)

Get the distortion coefficients from the two images img1 and img2. The default grid size is 128, the default search size is 256, and the default image size is 128.

Arguments

  • img1::Array{<:AbstractFloat,2}: The first image.
  • img2::Array{<:AbstractFloat,2}: The second image.
  • gridSize::Int: The size of the grid.
  • intrSize::Int: The size of the search.
  • srchSize::Int: The size of the search.

Returns

  • Vector{<:AbstractFloat}: The distortion coefficients.
source
ParticleHolography.load_gray2floatMethod
load_gray2float(path)

Load a grayscale image from a file and return it as a Array{Float32, 2} array.

Arguments

  • path::String: The path to the image file.

Returns

  • Array{Float32, 2}: The image as a Float32 array.
source
ParticleHolography.make_backgroundMethod
make_background(pathlist; mode=:mode)

Make a background image from a list of image paths. The background image is calculated by taking the mean or mode of the images in the list. The default mode is :mode.

Arguments

  • pathlist::Vector{String}: A list of image paths. glob() can be used to generate this list.
  • mode::Symbol: The mode to use for calculating the background. Options are :mean or :mode. Default is :mode.

Returns

  • Array{Float64, 2}: The background image.
source
ParticleHolography.modified_Cholesky_decompositionMethod
modified_Cholesky_decomposition(A)

This function decomposes the target matrix A to A = LDL^T. The diagonal components of the return matrix are the inverse of D, and the lower left components correspond to L. 対象行列Aを A = LDL^T に変換します。戻り値行列の対角成分は D の逆数、左下成分は L に一致します。

Arguments

  • A::Array{<:AbstractFloat,2}: The target matrix.

Returns

  • Array{AbstractFloat,2}: The decomposed matrix.
source
ParticleHolography.particle_bounding_boxesMethod
particle_bounding_boxes(d_bin_vol)

Detects the particles in a binary volume and returns the bounding boxes of the particles. This function performs three-dimensional element connected labeling on binary volumes. However, please note that it does not perform strict adjacent connectivity in the optical axis (z) direction. Strict adjacent connectivity in the optical axis direction may result in artifacts not being connected, potentially leading to the detection of many ghost particles. This function assumes that two particles never overlap at exactly the same X-Y coordinates. All elements that overlap in X-Y coordinates are considered connected. Therefore, this method is not suitable for accurate position detection of particles that overlap in X-Y coordinates. Additionally, this processing may have some effects, such as slightly elongating the bounding box of particles in the optical axis direction. However, this has minimal impact on the accuracy of particle position detection.

Arguments

  • d_bin_vol::CuArray{Bool, 3}: The binary volume of reconstructed holographic volume. Binarization can be done by thresholding the reconstructed volume. true values represent the particles and neighboring voxels and vice versa.

Returns

  • Dict{UUID, Vector{Int}}: The bounding boxes of the particles.
source
ParticleHolography.particle_coor_diamsFunction
particle_coor_diams(particle_bbs, d_vol, d_lpf_vol = nothing; depth_metrics = tamura, profile_smoothing_kernel = Kernel.gaussian(5,), diameter_metrics = equivalent_diameter)

Calculates the coordinates and diameters of the particles in the reconstructed volume with the bounding boxe dictionary. The depth of the particles is the maximum of the profile that is calculated using the depth_metrics function at each slice of the bounding box. The profile is then smoothed using the profile_smoothing_kernel. The x and y coordinates are calculated by finding the center of mass of the slice with the detected depth. The low pass filtered volume would be better for coordinate detection. The diameter of the particles is calculated using the diameter_metrics function. If the low pass filtered volume is provided, the coordinate is calculated using the low pass filtered volume.

Arguments

  • particle_bbs::Dict{UUID, Vector{Int}}: The bounding boxes of the particles.
  • d_vol::CuArray{N0f8, 3}: The reconstructed volume.

Optional arguments

  • d_lpf_vol::Union{CuArray{N0f8,3}, Nothing} = nothing: The low pass filtered volume.

Optional keyword arguments

  • depth_metrics::Function = tamura: The function that calculates the depth profile of the particles.
  • profile_smoothing_kernel = Kernel.gaussian((5,)): The kernel used for smoothing the depth profile.
  • diameter_metrics::Function = equivalent_diameter: The function that calculates the diameter of the particles.

Returns

  • Dict{UUID, Vector{Float32}}: The coordinates and diameters of the particles.
source
ParticleHolography.particle_coordinatesMethod
particle_coordinates(particle_bbs, d_vol; depth_metrics = tamura, profile_smoothing_kernel = Kernel.gaussian(5,))

Calculates the coordinates of the particles in the reconstructed volume with the bounding boxe dictionary. The depth of the particles is the maximum of the profile that is calculated using the depth_metrics function at each slice of the bounding box. The profile is then smoothed using the profile_smoothing_kernel. The x and y coordinates are calculated by finding the center of mass of the slice with the detected depth. The low pass filtered volume would be better for coordinate detection.

Arguments

  • particle_bbs::Dict{UUID, Vector{Int}}: The bounding boxes of the particles.
  • d_vol::CuArray{N0f8, 3}: The reconstructed volume.
  • depth_metrics::Function = tamura: The function that calculates the depth profile of the particles.
  • profile_smoothing_kernel = Kernel.gaussian((5,)): The kernel used for smoothing the depth profile.

Returns

  • Dict{UUID, Vector{Float32}}: The coordinates of the particles.
source
ParticleHolography.particleplot!Method
particleplot(data::Dict{UUID, Vector{Float32}}; scaling=(1.0, 1.0, -1.0), shift=(0.0, 0.0, 0.0), kwargs...)
particleplot!(data::Dict{UUID, Vector{Float32}}; scaling=(1.0, 1.0, -1.0), shift=(0.0, 0.0, 0.0), kwargs...)

Plots the particles in the 3D space. The data should be a dictionary with UUID keys and values as Vector{Float32}, which includes the coordinates of the particles. The scaling and shift parameters are used to transform the coordinates to the observed space.

In this plot, we display the x-y plane in the foreground by mapping the data's z-axis to the y-axis of Plots.scatter(). To align the data's y-axis with the hologram image coordinates, we set zflip = true by default. This results in a left-handed coordinate system, so we additionally set yflip = true to correct for this. The z-coordinate of the data (slice number) is opposite to the optical axis direction shown in the Introduction, so we set the z-axis data scaling to -1.0 for correct display.

Arguments

  • data::Dict{UUID, Vector{Float32}}: The coordinates of the particles.

Optional keyword arguments

  • scaling::Tuple{Float32, Float32, Float32} = (1.0, 1.0, -1.0): The scaling of the coordinates.
  • shift::Tuple{Float32, Float32, Float32} = (0.0, 0.0, 0.0): The shift of the coordinates.
  • kwargs...: Additional keyword arguments passed to Plots.scatter().
source
ParticleHolography.particleplot!Method
particleplot(data::Dict{UUID, Vector{Float32}}; scaling=(1.0, 1.0, -1.0), shift=(0.0, 0.0, 0.0), kwargs...)
particleplot!(data::Dict{UUID, Vector{Float32}}; scaling=(1.0, 1.0, -1.0), shift=(0.0, 0.0, 0.0), kwargs...)

Plots the particles in the 3D space. The data should be a dictionary with UUID keys and values as Vector{Float32}, which includes the coordinates of the particles. The scaling and shift parameters are used to transform the coordinates to the observed space.

In this plot, we display the x-y plane in the foreground by mapping the data's z-axis to the y-axis of Plots.scatter(). To align the data's y-axis with the hologram image coordinates, we set zflip = true by default. This results in a left-handed coordinate system, so we additionally set yflip = true to correct for this. The z-coordinate of the data (slice number) is opposite to the optical axis direction shown in the Introduction, so we set the z-axis data scaling to -1.0 for correct display.

Arguments

  • data::Dict{UUID, Vector{Float32}}: The coordinates of the particles.

Optional keyword arguments

  • scaling::Tuple{Float32, Float32, Float32} = (1.0, 1.0, -1.0): The scaling of the coordinates.
  • shift::Tuple{Float32, Float32, Float32} = (0.0, 0.0, 0.0): The shift of the coordinates.
  • kwargs...: Additional keyword arguments passed to Plots.scatter().
source
ParticleHolography.particleplotMethod
particleplot(data::Dict{UUID, Vector{Float32}}; scaling=(1.0, 1.0, -1.0), shift=(0.0, 0.0, 0.0), kwargs...)
particleplot!(data::Dict{UUID, Vector{Float32}}; scaling=(1.0, 1.0, -1.0), shift=(0.0, 0.0, 0.0), kwargs...)

Plots the particles in the 3D space. The data should be a dictionary with UUID keys and values as Vector{Float32}, which includes the coordinates of the particles. The scaling and shift parameters are used to transform the coordinates to the observed space.

In this plot, we display the x-y plane in the foreground by mapping the data's z-axis to the y-axis of Plots.scatter(). To align the data's y-axis with the hologram image coordinates, we set zflip = true by default. This results in a left-handed coordinate system, so we additionally set yflip = true to correct for this. The z-coordinate of the data (slice number) is opposite to the optical axis direction shown in the Introduction, so we set the z-axis data scaling to -1.0 for correct display.

Arguments

  • data::Dict{UUID, Vector{Float32}}: The coordinates of the particles.

Optional keyword arguments

  • scaling::Tuple{Float32, Float32, Float32} = (1.0, 1.0, -1.0): The scaling of the coordinates.
  • shift::Tuple{Float32, Float32, Float32} = (0.0, 0.0, 0.0): The shift of the coordinates.
  • kwargs...: Additional keyword arguments passed to Plots.scatter().
source
ParticleHolography.quadratic_distortion_correctionMethod
quadratic_distortion_correction(img, coefa)

Correct the quadratic distortion in the grayscale image img using the coefficients coefa. img have to be square, and coefa have to be a vector of 12 coefficients.

Arguments

  • img::Array{<:AbstractFloat,2}: The image to correct.
  • coefa::Vector{<:AbstractFloat}: The coefficients to correct the distortion.

Returns

  • Array{AbstractFloat,2}: The corrected image.
source
ParticleHolography.simultanious_equation_solverMethod
simultanious_equation_solver(chlskyMat,yacob,errorArray)

Solve the equation chlskyMat x = - yacob errorArray using the Cholesky decomposition matrix chlskyMat, Jacobian yacob, and error vector errorArray.

Arguments

  • chlskyMat::Array{<:AbstractFloat,2}: The Cholesky decomposition matrix.
  • yacob::Array{<:AbstractFloat,2}: The Jacobian matrix.
  • errorArray::Array{<:AbstractFloat,1}: The error vector.

Returns

  • Array{AbstractFloat,1}: The solution vector.
source
ParticleHolography.tamuraMethod
tamura(arr)

Calculates the Tamura coefficient of an array. The Tamura coefficient is defined as the standard deviation divided by the mean of the array. Please refer to the use in digital holography https://doi.org/10.1364/OL.36.001945

Arguments

  • arr::Array{Float32, 2}: The array for which the Tamura coefficient is calculated.

Returns

  • Float32: The Tamura coefficient of the array.
source
ParticleHolography.trajectoryplot!Method
trajectoryplot(paths::Vector{Vector{UUID}}, fulldict::Dict{UUID, Vector{Float32}}; colors=palette(:tab10), framerange=(0, 1024), scaling=(1.0, 1.0, -1.0), shift=(0.0, 0.0, 0.0), kwargs...)
trajectoryplot!(paths::Vector{Vector{UUID}}, fulldict::Dict{UUID, Vector{Float32}}; colors=palette(:tab10), framerange=(0, 1024), scaling=(1.0, 1.0, -1.0), shift=(0.0, 0.0, 0.0), kwargs...)

Plots particle trajectories with different colors in 3D space. paths is a Vector{Vector{UUID}}, where each element represents the trajectory of particles considered identical. fulldict is a dictionary containing the UUIDs of particles for all conceivable frames as keys and their coordinates as values, generated using gen_fulldict. colors specifies the palette for trajectory colors, defaulting to palette(:tab10). For scaling and shift, refer to particleplot. framerange specifies the range of frames to plot. For example, when generating fulldict for N consecutive frames, setting framerange=(N-10+1, N) will plot only the trajectories of the last 10 frames. When creating animations with the @anim macro, setting framerange from (1,1) to (1,N) will create an N-frame animation.

Arguments

  • paths::Vector{Vector{UUID}}: The trajectories of the particles.
  • fulldict::Dict{UUID, Vector{Float32}}: The dictionary of all particles.

Optional keyword arguments

  • colors::Vector{Colorant} = palette(:tab10): The colors of the trajectories.
  • framerange::Tuple{Int, Int} = (0, 1024): The range of frames to plot.
  • scaling::Tuple{Float32, Float32, Float32} = (1.0, 1.0, -1.0): The scaling of the coordinates.
  • shift::Tuple{Float32, Float32, Float32} = (0.0, 0.0, 0.0): The shift of the coordinates.
  • kwargs...: Additional keyword arguments passed to Plots.@series.
source
ParticleHolography.trajectoryplot!Method
trajectoryplot(paths::Vector{Vector{UUID}}, fulldict::Dict{UUID, Vector{Float32}}; colors=palette(:tab10), framerange=(0, 1024), scaling=(1.0, 1.0, -1.0), shift=(0.0, 0.0, 0.0), kwargs...)
trajectoryplot!(paths::Vector{Vector{UUID}}, fulldict::Dict{UUID, Vector{Float32}}; colors=palette(:tab10), framerange=(0, 1024), scaling=(1.0, 1.0, -1.0), shift=(0.0, 0.0, 0.0), kwargs...)

Plots particle trajectories with different colors in 3D space. paths is a Vector{Vector{UUID}}, where each element represents the trajectory of particles considered identical. fulldict is a dictionary containing the UUIDs of particles for all conceivable frames as keys and their coordinates as values, generated using gen_fulldict. colors specifies the palette for trajectory colors, defaulting to palette(:tab10). For scaling and shift, refer to particleplot. framerange specifies the range of frames to plot. For example, when generating fulldict for N consecutive frames, setting framerange=(N-10+1, N) will plot only the trajectories of the last 10 frames. When creating animations with the @anim macro, setting framerange from (1,1) to (1,N) will create an N-frame animation.

Arguments

  • paths::Vector{Vector{UUID}}: The trajectories of the particles.
  • fulldict::Dict{UUID, Vector{Float32}}: The dictionary of all particles.

Optional keyword arguments

  • colors::Vector{Colorant} = palette(:tab10): The colors of the trajectories.
  • framerange::Tuple{Int, Int} = (0, 1024): The range of frames to plot.
  • scaling::Tuple{Float32, Float32, Float32} = (1.0, 1.0, -1.0): The scaling of the coordinates.
  • shift::Tuple{Float32, Float32, Float32} = (0.0, 0.0, 0.0): The shift of the coordinates.
  • kwargs...: Additional keyword arguments passed to Plots.@series.
source
ParticleHolography.trajectoryplotMethod
trajectoryplot(paths::Vector{Vector{UUID}}, fulldict::Dict{UUID, Vector{Float32}}; colors=palette(:tab10), framerange=(0, 1024), scaling=(1.0, 1.0, -1.0), shift=(0.0, 0.0, 0.0), kwargs...)
trajectoryplot!(paths::Vector{Vector{UUID}}, fulldict::Dict{UUID, Vector{Float32}}; colors=palette(:tab10), framerange=(0, 1024), scaling=(1.0, 1.0, -1.0), shift=(0.0, 0.0, 0.0), kwargs...)

Plots particle trajectories with different colors in 3D space. paths is a Vector{Vector{UUID}}, where each element represents the trajectory of particles considered identical. fulldict is a dictionary containing the UUIDs of particles for all conceivable frames as keys and their coordinates as values, generated using gen_fulldict. colors specifies the palette for trajectory colors, defaulting to palette(:tab10). For scaling and shift, refer to particleplot. framerange specifies the range of frames to plot. For example, when generating fulldict for N consecutive frames, setting framerange=(N-10+1, N) will plot only the trajectories of the last 10 frames. When creating animations with the @anim macro, setting framerange from (1,1) to (1,N) will create an N-frame animation.

Arguments

  • paths::Vector{Vector{UUID}}: The trajectories of the particles.
  • fulldict::Dict{UUID, Vector{Float32}}: The dictionary of all particles.

Optional keyword arguments

  • colors::Vector{Colorant} = palette(:tab10): The colors of the trajectories.
  • framerange::Tuple{Int, Int} = (0, 1024): The range of frames to plot.
  • scaling::Tuple{Float32, Float32, Float32} = (1.0, 1.0, -1.0): The scaling of the coordinates.
  • shift::Tuple{Float32, Float32, Float32} = (0.0, 0.0, 0.0): The shift of the coordinates.
  • kwargs...: Additional keyword arguments passed to Plots.@series.
source