API Reference

Covariance Construction Functions

CloudClean.cov_avg!Function
cov_avg!(bimage, ism, bism, in_image; Np::Int=33, widx::Int=129, widy::Int=129, ftype::Int=32)

Key function for constructing the (shifted and multiplied) versions of the input image used to quickly estimate the local covariance matrix at a large number of locations. The main output is in the preallocated bism which is used as an input to build_cov!.

Arguments:

  • bimage: preallocated output array for the boxcar smoothed unshifted image
  • ism: preallocated intermediate array for the input image times itself shifted
  • bism: preallocated output array to store boxcar-smoothed image products for all shifts
  • in_image: input image the local covariance of which we want to estimate

Keywords:

  • Np::Int: size of local covariance matrix in pixels (default 33)
  • widx::Int: width of boxcar window in x which determines size of region used for samples for the local covariance estimate (default 129)
  • widy::Int: width of boxcar window in y which determines size of region used for samples for the local covariance estimate (default 129)
  • ftype::Int: determine the Float precision, 32 is Float32, otherwise Float64
source
CloudClean.cov_avg_sym!Function
cov_avg_sym!(bimage, ism, bism, in_image; Np::Int=33, widx::Int=129, widy::Int=129, ftype::Int=32)

Key function for constructing the (shifted and multiplied) versions of the input image used to quickly estimate the local covariance matrix at a large number of locations. The main output is in the preallocated bism which is used as an input to build_cov!. Version of cov_avg! using symmetric training shifts.

Arguments:

  • bimage: preallocated output array for the boxcar smoothed unshifted image
  • ism: preallocated intermediate array for the input image times itself shifted
  • bism: preallocated output array to store boxcar-smoothed image products for all shifts
  • in_image: input image the local covariance of which we want to estimate

Keywords:

  • Np::Int: size of local covariance matrix in pixels (default 33)
  • widx::Int: width of boxcar window in x which determines size of region used for samples for the local covariance estimate (default 129)
  • widy::Int: width of boxcar window in y which determines size of region used for samples for the local covariance estimate (default 129)
  • ftype::Int: determine the Float precision, 32 is Float32, otherwise Float64
source
CloudClean.boxsmooth!Function
boxsmooth!(out::AbstractArray, arr::AbstractArray, tot::Array{T,1}, widx::Int, widy::Int)

Boxcar smooths an input image (or paddedview) arr with window size widx by widy. We pass the original image size sx and sy to help handle image views.

Arguments:

  • out::AbstractArray: preallocated output array for the boxcar smoothed image
  • arr::AbstractArray: input array for which boxcar smoothing is computed (generally paddedview)
  • tot::Array{T,1}: preallocated array to hold moving sums along 1 dimension
  • widx::Int: size of boxcar smoothing window in x
  • widy::Int: size of boxcar smoothing window in y
source
CloudClean.outest_boundsFunction
outest_bounds(cx,sx) -> px0

Helper function to find maximum padding in pixels required to accomodate all query points cx outside of the image size 1:sx.

Arguments:

  • cx: list of integer star centers (in either x or y)
  • sx: image dimension along the axis indexed by cx

Outputs:

  • px0: maximum padding in pixels required to accomodate all query points
source

Per Location Functions

CloudClean.gen_pix_mask_trivialFunction
gen_pix_mask_trivial(kmasked2d; Np=33) -> kstar, kcond

Flatten a pixel mask and calculate the number of pixels used for the conditional infill.

Arguments:

  • kmasked2d: A 2D array representing the masked pixels.

Keywords:

  • Np: An optional integer specifying the number of pixels in a side (default: 33).

Returns:

  • kstar: A flattened version of the input kmasked2d array.
  • kcond: The count of unmasked pixels in the kstar array.

Examples:

julia> kmasked2d = rand(Bool, 33, 33)
julia> kstar, kcond = gen_pix_mask_trivial(kmasked2d, Np=33)
source
CloudClean.gen_pix_mask_circFunction
gen_pix_mask_circ(kmasked2d, circmask; Np=33) -> kstar, kcond

Generate a circular pixel mask and calculate the number of pixels used for the conditional infill.

Arguments:

  • kmasked2d: A 2D array representing the masked pixels.
  • circmask: A 2D array representing the circular mask.

Keywords:

  • Np: An optional integer specifying the number of pixels in a side (default: 33).

Returns:

  • kstar: A copy of the input kmasked2d array with circular masking applied.
  • kcond: The count of unmasked pixels in the kstar array.

Examples:

julia julia> kmasked2d = rand(Bool, 33, 33) julia> circmask = kstar_circle_mask(33,rlim=256) julia> kstar, kcond = gen_pix_mask_circ(kmasked2d, circmask, Np=33)

source
CloudClean.condCovEst_wdiagFunction
condCovEst_wdiag(cov_loc,μ,km,data_in;Np=33,export_mean=false,n_draw=0) -> out

Using a local covariance matrix estimate cov_loc and a set of known ("good") pixels km, this function computes a prediction for the mean value of masked pixels and the covariance matrix of the masked pixels. The output list can conditionally include the mean reconstruction and draws from the distribution of reconstructions.

Arguments:

  • cov_loc: local covariance matrix
  • μ: vector containing mean value for each pixel in the patch
  • km: unmasked pixels
  • data_in: input image

Keywords:

  • Np: size of local covariance matrix in pixels (default 33)
  • export_mean: when true, returns the mean conditional prediction for the "hidden" pixels (default false)
  • n_draw: when nonzero, returns that number of realizations of the conditional infilling (default 0)

Outputs:

  • out[1]: input image returned with masked pixels replaced with mean prediction
  • out[2]: input image returned with masked pixels replaced with a draw from the predicted distribution
source
CloudClean.build_cov!Function
build_cov!(cov::Array{T,2},μ::Array{T,1},cx::Int,cy::Int,bimage::Array{T,2},bism::Array{T,4},Np::Int,widx::Int,widy::Int) where T <:Union{Float32,Float64}

Constructs the local covariance matrix and mean for an image patch of size Np x Np pixels around a location of interest (cx,cy). The construction is just a lookup of pixel values from the stored boxcar-smoothed copies of the input image times itself shifted in bism. Passing the smoothed image bimage and the widths of the boxcar mean widx and widy is helpful for the mean and normalization. The covariance and mean are updated in place for speed since this operation may be performed billions of times since we construct a new covariance matrix for every detection. Math may either be performed Float32 or Float64.

Arguments:

  • cov::Array{T,2}: preallocated output array for local covariance matrix
  • μ::Array{T,1}: preallocated output vector for local mean
  • cx::Int: x-coordinate of the center of the local region
  • cy::Int: y-coordinate of the center of the local region
  • bimage::Array{T,2}: boxcar smoothed unshifted image
  • bism::Array{T,4}: boxcar-smoothed image products for all shifts
  • Np::Int: size of local covariance matrix in pixels
  • widx::Int: width of boxcar window in x which determines size of region used for samples for the local covariance estimate
  • widy::Int: width of boxcar window in y which determines size of region used for samples for the local covariance estimate
source
CloudClean.build_cov_sym!Function
build_cov_sym!(cov::Array{T,2},μ::Array{T,1},cx::Int,cy::Int,bimage::Array{T,2},bism::Array{T,4},Np::Int,widx::Int,widy::Int) where T <:Union{Float32,Float64}

Constructs the local covariance matrix and mean for an image patch of size Np x Np pixels around a location of interest (cx,cy). The construction is just a lookup of pixel values from the stored boxcar-smoothed copies of the input image times itself shifted in bism. Passing the smoothed image bimage and the widths of the boxcar mean widx and widy is helpful for the mean and normalization. The covariance and mean are updated in place for speed since this operation may be performed billions of times since we construct a new covariance matrix for every detection. Math may either be performed Float32 or Float64. Version of build_cov! that uses symmetric training shifts.

Arguments:

  • cov::Array{T,2}: preallocated output array for local covariance matrix
  • μ::Array{T,1}: preallocated output vector for local mean
  • cx::Int: x-coordinate of the center of the local region
  • cy::Int: y-coordinate of the center of the local region
  • bimage::Array{T,2}: boxcar smoothed unshifted image
  • bism::Array{T,4}: boxcar-smoothed image products for all shifts
  • Np::Int: size of local covariance matrix in pixels
  • widx::Int: width of boxcar window in x which determines size of region used for samples for the local covariance estimate
  • widy::Int: width of boxcar window in y which determines size of region used for samples for the local covariance estimate
source

Image Preprocessing

CloudClean.prelim_infill!Function
prelim_infill!(testim,bmaskim,bimage,bimageI,testim2,bmaskim2,goodpix;widx=19,widy=19,ftype::Int=32,widmult=1.4)

This intial infill replaces masked pixels with a guess based on a smoothed boxcar. For large masked regions, the smoothing scale is increased. If this iteration takes too long/requires too strong of masking, the masked pixels are replaced with the median of the image.

We use 3 copies of the input image and mask image. The first is a view (with reflective boundary condition padding) with the pixels to be infilled replaced with zeros, the second is allocated to hold various smoothings of the image, and the third holds the output image which contains our best infill guess. A final Bool array of size corresponding to the image is used to keep track of pixels that have safe infill values.

Arguments:

  • testim: input image which requires infilling
  • bmaskim: input mask indicating which pixels require infilling
  • bimage: preallocated array for smoothed version of input image
  • bimageI: preallocated array for smoothed mask counting the samples for each estimate
  • testim2: inplace modified ouptut array for infilled version of image input
  • bmaskim2: inplace modified mask to keep track of which pixels still need infilling
  • goodpix: preallocated array for Bool indexing pixels with good infill

Keywords:

  • widx: initial size of boxcar smoothing window in x (default 19)
  • widy: initial size of boxcar smoothing window in y (default 19)
  • ftype::Int: determine the Float precision, 32 is Float32, otherwise Float64
  • widmult: multiplicative factor for increasing the smoothing scale at each iteration step
source
CloudClean.sig_iqrFunction
sig_iqr(x)

Calculate the normalized interquartile range (IQR) of an array as a robust measure of standard deviation.

Arguments

  • x: A 1D array or iterable.

Returns

  • The normalized IQR, computed as the IQR divided by 1.34896.

Examples

julia> x = [1, 2, 3, 4, 5, 6, 7, 8, 9, 10]
julia> result = sig_iqr(x)
source
CloudClean.im_subrngFunction
im_subrng(jx,jy,cx,cy,sx,sy,px0,py0,stepx,stepy,padx,pady,tilex,tiley) -> xrng, yrng, star_ind

Computes the flux a star must have so that the PSF-based masking using thr would require a larger stamp area. Used for computational savings.

Arguments:

  • jx: tile index along x
  • jy: tile index along y
  • cx: list of stellar x-coordinates
  • cy: list of stellar y-coordinates
  • sx: size of original image in x
  • sy: size of original image in y
  • px0: maximal padding in x to account for stars outside image
  • py0: maximal padding in y to account for stars outside image
  • stepx: tiling step size in x
  • stepy: tiling step size in y
  • padx: tile padding in x required to account for local stamp size, sample size, and pixels outside the image
  • pady: tile padding in x required to account for local stamp size, sample size, and pixels outside the image
  • tilex: total number of tile divisions along x
  • tiley: total number of tile divisions along y

Outputs:

  • xrng: slicing range of the tile in x
  • yrng: slicing range of the tile in y
  • star_ind: Bool mask of all stars falling within the tile (subimage)
source
CloudClean.add_sky_noise!Function
add_sky_noise!(testim2, maskim, sig_iqr; seed=2021)

Add sky noise to pixels in an image specified by a given mask.

Arguments:

  • testim2: A mutable 2D array representing the image.
  • maskim: A 2D array representing the mask.
  • sig_iqr: The standard deviation of the noise distribution, generally calculated as the normalized IQR.

Keywords:

  • seed: An optional integer specifying the random number generator seed (default: 2021).

Returns:

  • Modifies testim2 in place by adding sky noise to the masked pixels.

Examples:

julia> testim2 = rand(100, 100)
julia> maskim = rand(Bool, 100, 100)
julia> sig_iqr = 0.5
julia> add_sky_noise!(testim2, maskim, sig_iqr, seed=2021)
source
CloudClean.add_noise!Function
add_noise!(testim2,gain;seed=2021)

Adds noise to an image that matches the Poisson noise of the pixel counts. A random seed to set a local random generator is provided for reproducible unit testing.

Arguments:

  • testim2: input image which had infilling
  • gain: gain of detector to convert from photon count noise to detector noise

Keywords:

  • seed: random seed for random generator
source
CloudClean.kstar_circle_maskFunction
kstar_circle_mask(Np;rlim=256) -> circmask

Generates a Bool mask for pixels beyond a given (squared) radius of the center of an image.

Arguments:

  • Np: size of image stamp

Keywords:

  • rlim: squared radius (in pixels^2) beyond which pixels should be masked (default 256)

Outputs:

  • circmask: static Bool mask used for assigning pixels beyond some radius of the stellar center as "ignored"
source

Full Image Processing Functions

CloudClean.proc_discreteFunction
proc_discrete(x_locs, y_locs, raw_image, mask_image; Np=33, widx=129, widy=widx, tilex=1, tiley=tilex, seed=2021, ftype::Int=32, ndraw=0) -> out_mean, out_draw

Process an image with a mask, replacing masked pixels with either a mean or draw from a distribution resembling the local pixel-pixel covariance structure in the image.

Arguments:

  • x_locs: A 1D array representing the location centers (in the x coordinate) for infilling.
  • y_locs: A 1D array representing the location centers (in the y coordinate) for infilling.
  • raw_image: A 2D array representing the input image.
  • mask_image: A 2D array representing the mask.

Keywords:

  • Np: An optional integer specifying the number of pixels in a side (default: 33).
  • widx: An optional integer specifying the width of the region used for training the local covariance in the x-direction (default: 129).
  • widy: An optional integer specifying the width of the region used for training the local covariance in the y-direction (default: widx).
  • tilex: An optional integer specifying the number of tiles in the x-direction for subdividing the image (default: 1).
  • tiley: An optional integer specifying the number of tiles in the y-direction for subdividing the image (default: tilex).
  • seed: An optional integer specifying the random number generator seed (default: 2021).
  • ftype: An optional integer specifying the floating-point precision type (32 or 64) (default: 32).
  • rlim: Radius limit for the radial mask beyond which pixels are not used for conditioning (units are pixels^2). (default: Inf)
  • ndraw: An optional integer specifying the number of draws of samples from the statistical distribution of possible masked pixel values (default: 0).

Returns

  • If ndraw is 0, returns the debiased image as a 2D array.
  • If ndraw is greater than 0, returns the debiased image as a 2D array and an array of ndraw draws.

Examples

julia> raw_image = rand(100, 100)
julia> mask_image = kstar_circle_mask(100,rlim=256)
julia> result = proc_continuous([50],[50],raw_image, mask_image, Np=33, widx=129, seed=2021)
source
CloudClean.proc_continuousFunction
proc_continuous(raw_image, mask_image; Np=33, widx=129, widy=widx, tilex=1, tiley=tilex, seed=2021, ftype::Int=32, ndraw=0) -> out_mean, out_draw

Process an image with a mask, replacing masked pixels with either a mean or draw from a distribution resembling the local pixel-pixel covariance structure in the image.

Arguments:

  • raw_image: A 2D array representing the input image.
  • mask_image: A 2D array representing the mask.

Keywords:

  • Np: An optional integer specifying the number of pixels in a side (default: 33).
  • widx: An optional integer specifying the width of the region used for training the local covariance in the x-direction (default: 129).
  • widy: An optional integer specifying the width of the region used for training the local covariance in the y-direction (default: widx).
  • tilex: An optional integer specifying the number of tiles in the x-direction for subdividing the image (default: 1).
  • tiley: An optional integer specifying the number of tiles in the y-direction for subdividing the image (default: tilex).
  • seed: An optional integer specifying the random number generator seed (default: 2021).
  • ftype: An optional integer specifying the floating-point precision type (32 or 64) (default: 32).
  • ndraw: An optional integer specifying the number of draws of samples from the statistical distribution of possible masked pixel values (default: 0).

Returns

  • If ndraw is 0, returns the debiased image as a 2D array.
  • If ndraw is greater than 0, returns the debiased image as a 2D array and an array of ndraw draws.

Examples

julia> raw_image = rand(100, 100)
julia> mask_image = rand(Bool, 100, 100)
julia> result = proc_continuous(raw_image, mask_image, Np=33, widx=129, seed=2021)
source