Reconstruction
We assume that the necessary preprocessing has been completed and an array of type ParticleHolography.CuWaveFront
, d_wavefront
, has been obtained. Using this array, we perform a 3D reconstruction of the observed volume. In Gabor holography, we have already described the functions for obtaining the projection images in the optical axis direction of the observed volume. In addition to this, depending on the required information, you can select the following reconstruction functions.
ParticleHolography.cu_get_reconst_xyprojection
— Functioncu_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. SeeCuTransfer
.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.
ParticleHolography.cu_get_reconst_vol
— Functioncu_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. SeeCuTransfer
.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 isN0f8
.
Returns
CuArray{return_type,3}
: The reconstructed intensity volume.
ParticleHolography.cu_get_reconst_complex_vol
— Functioncu_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. SeeCuTransfer
.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.
ParticleHolography.cu_get_reconst_vol_and_xyprojection
— Functioncu_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. SeeCuTransfer
.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 isN0f8
.
Returns
CuArray{return_type,3}
: The reconstructed volume.CuArray{Float32,2}
: The XY projection of the reconstructed volume.