API Reference

Covariance Construction Functions

CloudCovErr.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
CloudCovErr.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
CloudCovErr.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 Star Functions

CloudCovErr.stamp_cutterFunction
stamp_cutter(cx,cy,residimIn,star_im,maskim;Np=33) -> data_in, stars_in, kmasked2d

Cuts out local stamps around each star of the various input images to be used for per star statistics calculations.

Arguments:

  • cx: center coorindate x of the stamp
  • cy: center coorindate y of the stamp
  • residimIn: residual image with infilling from which covariance was estimated
  • star_im: input image of model of stars only. abs(modim-skyim)
  • maskim: input image of upstream masked pixels

Keywords:

  • Np: size of covariance matrix footprint around each star (default 33)

Outputs:

  • data_in: local stamp of the (non-infilled) residual image
  • stars_in: local stamp of model of stars only
  • kmasked2d: local stamp of upstream masked pixels
source
CloudCovErr.gen_pix_maskFunction
gen_pix_mask(kmasked2d,psfmodel,circmask,x_star,y_star,flux_star;Np=33,thr=20) -> psft, kstar, kpsf2d, kcond0, kcond, kpred, dnt

Assigns pixels in the local subimage around a star to either be "good", "hidden", or "ignored" based on user settings and the flux of the star. Reads in masked pixels from the quality flags on pixels coming from the community pipeline kmasked2d, a PSF model for the star, and a precomputed circular mask circmask to exclude pixels at a large radius from the stellar center since they have little impact on the regression of hidden pixels. The pixels assigned as "hidden" and to be interpolated are determined by a thr on the pixel values for flux_star times the PSF model. We use a parametric PSFs that varies with position and query the PSF at the stellar position for each star.

Arguments:

  • kmasked2d: Bool mask from upstream pixel quality flags to assign pixels as "ignored"
  • psfmodel: parametric PSF model that can be queried at different positions
  • circmask: static Bool mask assigning pixels beyond some radius of the stellar center as "ignored"
  • x_star: x-coordinate of the star (used only for flexible PSF model query)
  • y_star: y-coordinate of the star (used only for flexible PSF model query)
  • flux_star: flux of star in ADU to determine how large a region to make "hidden"

Keywords:

  • Np: size of local covariance matrix in pixels (default 33)
  • thr: threshold for psf-based masking of the residuals (larger more "hidden")

Outputs:

  • psft: static array (image) of the stellar PSF
  • kstar: Boolean indexes the NOT "good" pixels
  • kpsf2d: Boolean indexes the "hidden" pixels
  • kcond0: initial number of "good" pixels
  • kcond: final number of "good" pixels after fallbacks
  • kpred: the number of pixels "hidden"
  • dnt: quality flag bits on the solution
source
CloudCovErr.condCovEst_wdiagFunction
condCovEst_wdiag(cov_loc,μ,km,kpsf2d,data_in,stars_in,psft;Np=33,export_mean=false,n_draw=0,diag_on=true) -> out

Using a local covariance matrix estimate cov_loc, a set of masked pixels km and "hidden" pixels kpsf2d, this function computes a prediction for the mean value of the kpsf2d pixels and the covariance matrix of the kpsf2d pixels. In terms of statistics use to adjust the photometry of a star, we are only interested in the pixels masked as a result of the star (i.e. not a detector defect or cosmic ray nearby). The residual image data_in and a model of the counts above the background coming from the star stars_in for the local patch are also inputs of the function. Correction factors for the photometric flux and flux uncertainities are outputs as well as a chi2 value for the "good" 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: masked pixels (either bad or do not want to use for conditioning)
  • kpsf2d: pixels masked due to the star of interest
  • data_in: (non-infilled) residual image in local patch
  • psft: static array (image) of the stellar PSF

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)
  • diag_on: flag for adding to the pixelwise uncertainty based on the photoelectron counts of the modeled star (default true)

Outputs:

  • out[1][1]: flux uncertainity of the star
  • out[1][2]: flux uncertainity of the star assuming the covariance matrix were diagonal
  • out[1][3]: flux correction which must be added to correct the input flux estimate
  • out[1][4]: flux correction coming from the residuals (fdb_res)
  • out[1][5]: flux correction coming from the predicted background (fdb_pred)
  • out[1][6]: chi2 for the "good" pixels under cov_loc as a metric on how good our assumptions are
  • out[2]: local region (image) with "hidden" pixels replaced by the mean conditional estimate (optional output)
  • out[end]: local region (image) with "hidden" pixels replaced by the draws from the conditional distribution (optional output). Array is flattened to npix x n_draw.
source
CloudCovErr.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

Image Infill and Masking

CloudCovErr.prelim_infill!Function
prelim_infill!(testim,bmaskim,bimage,bimageI,testim2,bmaskim2,goodpix,ccd;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
  • ccd: string name of FITS extension for verbose cmdline printing

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
CloudCovErr.gen_mask_staticPSF!Function
gen_mask_staticPSF!(maskd, psfstamp, x_stars, y_stars, flux_stars, thr=20)

Generate a mask for an input image (which is usually an image of model residuals) that excludes the cores of stars (which are often mismodeled). In this function, we use a fixed PSF psfstamp for all sources, and adjust the masking fraction based on the stellar flux and a threshold thr. A more general position dependent PSF model could be used with a slight generalization of this function, but is likely overkill for the problem of making a mask.

Arguments:

  • maskd: bool image to which mask will be added (bitwise or)
  • psfstamp: simple 2D array of a single PSF to be used for the whole image
  • x_stars: list of source x positions
  • y_stars: list of source y positions
  • flux_stars: list of source fluxes

Keywords:

  • thr: threshold used for flux-dependent masking (default 20)
source
CloudCovErr.gen_mask_staticPSF2!Function
gen_mask_staticPSF2!(maskd, psfstamp, psfstamp1, x_stars, y_stars, flux_stars, thr=20)

Generate a mask for an input image (which is usually an image of model residuals) that excludes the cores of stars (which are often mismodeled). In this function, we use a small fixed PSF psfstamp1 for all faint sources, and adjust the masking fraction based on the stellar flux and a threshold thr. Only for source bright enough to need a larger PSF stamp do we use psfstamp, which saves some computational cost.

Arguments:

  • maskd: bool image to which mask will be added (bitwise or)
  • psfstamp: simple 2D array of a single PSF to be used for bright stars in the whole image
  • psfstamp1: simple 2D array of a single PSF to be used for faint stars in the whole image
  • x_stars: list of source x positions
  • y_stars: list of source y positions
  • flux_stars: list of source fluxes

Keywords:

  • thr: threshold used for flux-dependent masking (default 20)
source
CloudCovErr.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
CloudCovErr.add_sky_noise!Function
add_sky_noise!(testim2,maskim,skyim,gain;seed=2021)

Adds noise to the infill that matches the Poisson noise of a rough estimate for the sky background. A random seed to set a local random generator is provided for reproducible unit testing.

Arguments:

  • testim2: input image which had infilling
  • maskim: mask of pixels which were infilled
  • skyim: rough estimate of sky background counts
  • gain: gain of detector to convert from photon count noise to detector noise

Keywords:

  • seed: random seed for random generator
source
CloudCovErr.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
CloudCovErr.findmaxpsfFunction
findmaxpsf(psfstamp1;thr=20) -> flim

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:

  • psfstamp1: a small image of a representative PSF

Keywords:

  • thr: threshold used to determine which pixels are bright enough to be "hidden"

Outputs:

  • flim: maximum flux that can be masked by thr without exceeding the PSF stamp footprint
source
CloudCovErr.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

DECam Specific Functions

CloudCovErr.decam.read_decamFunction
read_decam(base,date,filt,vers,ccd;corrects7=true) -> ref_im, d_im

Read in raw image files associated with exposures obtain on the DarkEnergyCamera. Returns the image and a quality flag mask image. See NOAO handbook for more details on what is contained in each file and how they are obtained.

Arguments:

  • base: parent directory and file name prefix for exposure files
  • date: date_time of the exposure
  • filt: optical filter used to take the exposure
  • vers: NOAO community processing version number
  • ccd: which ccd we are pulling the image for

Keywords:

  • corrects7: use crowdsource load to read ccd "S7" to correct for floating amplifier on half of the chip (default true)

Output:

  • ref_im: image of photoelectron counts from observation on DECam
  • d_im: quality flag mask image from NOAO community pipeline

Example

ref_im, d_im = read_decam("/n/fink2/decaps/c4d_","170420_040428","g","v1","N14")
source
CloudCovErr.decam.read_crowdsourceFunction
read_crowdsource(basecat,date,filt,vers,ccd) -> x_stars, y_stars, flux_stars, decapsid, gain, mod_im, sky_im, wcol, w

Read in outputs of crowdsource, a photometric pipeline. To pair with an arbitrary photometric pipeline, an analogous read in function should be created. The relevant outputs are the model image (including the sources) so that we can produce the residual image, the sky/background model (no sources), and the coordinates of the stars. The survey id number is also readout of the pipeline solution file to help cross-validate matching of the CloudCovErr outputs and the original sources. The empirical gain is read out of the header (for other photometric pipelines which don't perform this estiamte, the gain from DECam is likely sufficient). All columns from the photometric catalogue are also read in at this point to be rexported with the CloudCovErr outputs.

Arguments:

  • basecat: parent directory of the cat directory holding all of the single-epoch crowdsource catalogue files
  • date: date_time of the exposure
  • filt: optical filter used to take the exposure
  • vers: NOAO community processing version number
  • ccd: which ccd we are pulling the image for

Keywords:

  • corrects7: use crowdsource load to read ccd "S7" to correct for floating amplifier on half of the chip (default true)

Output:

  • x_stars: list of source x-coordinates (accounting for indexing order and start point)
  • y_stars: list of source y-coordinates (accounting for indexing order and start point)
  • flux_stars: list of stellar fluxes in ADU
  • decapsid: list of survey id number for each detection
  • gain: gain of detector to convert from photon count noise to detector noise
  • mod_im: model image (including the sources) from photometric pipeline
  • sky_im: sky image (background estimate) from photometric pipeline
  • wcol: list of all column names in photometric catalogue
  • w: list of all column values in photometric catalogue
source
CloudCovErr.decam.inject_renameFunction
inject_rename(fname) -> ifname

Convenience renaming file paths to read data from injection tests which are stored separately from data obtained on DECam for data provenance purposes.

Arguments:

  • fname: file name for exposure data from DECam

Output:

  • ifname: corresponding file name for injection tests into that exposure
source
CloudCovErr.decam.load_psfmodel_csFunction
load_psfmodel_cs(base,date,filt,vers,ccd) -> psfmodel

Julia wrapper function for the PyCall that reads the position dependent psfmodel produced by crowdsource from the catalogue file for a given exposure and ccd. The returned psfmodel takes an x- and y-position for the source location and the size of the desired psfstamp (the stamps are square and required to be odd).

Arguments:

  • base: parent directory and file name prefix for catalogue files
  • date: date_time of the exposure
  • filt: optical filter used to take the exposure
  • vers: NOAO community processing version number
  • ccd: which ccd we are pulling the image for

Output:

  • psfmodel: function that returns PSF stamp from parametric PSF model that is a function of position
source
CloudCovErr.decam.save_fxnFunction
save_fxn(wcol,w,basecat,date,filt,vers,ccd)

Saves CloudCovErr.jl outputs and initial photometric catalogue outputs to a new single-epoch catalogue. Massages types of columns to reduce data storage size. Converts the native CloudCovErr.jl output of the bias offset value into a cflux corrected flux column for the ease of catalogue users.

Arguments:

  • wcol: list of all column names in photometric catalogue
  • w: list of all column values in photometric catalogue
  • basecat: parent directory of the cat directory holding all of the single-epoch crowdsource catalogue files
  • date: date_time of the exposure
  • filt: optical filter used to take the exposure
  • vers: NOAO community processing version number
  • ccd: which ccd we are pulling the image for
source
CloudCovErr.decam.get_catnamesFunction
get_catnames(f) -> extnames

Reads list of extension names from an open FITS file to determine which CCDs have completed photometric catalogues and are eligible for CloudCovErr.jl.

Arguments:

  • f: an open FITS file handle containing crowdsource catalogues

Outputs:

  • extnames: list of CCDs that have photometric catalogues in f
source
CloudCovErr.decam.proc_ccdFunction
proc_ccd(base,date,filt,vers,basecat,ccd;thr=20,outthr=20000,Np=33,corrects7=true,widx=129,widy=widx,tilex=1,tiley=tilex,ftype::Int=32)

Primary run function for a given CCD image of a larger exposure.

Arguments:

  • base: parent directory and file name prefix for exposure files
  • date: date_time of the exposure
  • filt: optical filter used to take the exposure
  • vers: NOAO community processing version number
  • basecat: parent directory of the cat directory holding all of the single-epoch crowdsource catalogue files
  • ccd: which ccd we are pulling the image for

Keywords:

  • thr: threshold used for flux-dependent masking (default 20)
  • outthr: threshold for residual-based masking (default 20000)
  • Np: size of local covariance matrix in pixels (default 33)
  • corrects7: use crowdsource load to read ccd "S7" to correct for floating amplifier on half of the chip (default true)
  • widx: width of boxcar window in x which determines size of region used for samples for the local covariance estimate (default 129)
  • widy: width of boxcar window in y which determines size of region used for samples for the local covariance estimate (default 129)
  • tilex: total number of tile divisions along x (default 1)
  • tiley: total number of tile divisions along y (default tilex)
  • ftype::Int: determine the Float precision, 32 is Float32, otherwise Float64
source
CloudCovErr.decam.proc_allFunction
proc_all(base,date,filt,vers,basecat;ccdlist=String[],resume=false,corrects7=true,thr=20,outthr=20000,Np=33,widx=129,widy=widx,tilex=1,tiley=tilex,ftype::Int=32)

Exposure level run function the manages which ccds to run and calls proc_ccd serially.

Arguments:

  • base: parent directory and file name prefix for exposure files
  • date: date_time of the exposure
  • filt: optical filter used to take the exposure
  • vers: NOAO community processing version number
  • basecat: parent directory of the cat directory holding all of the single-epoch crowdsource catalogue files

Keywords:

  • ccdlist: run only ccds in this list
  • resume: if the exposure is partially complete, resume running from where it left off (default false)
  • corrects7: use crowdsource load to read ccd "S7" to correct for floating amplifier on half of the chip (default true)
  • thr: threshold used for flux-dependent masking (default 20)
  • outthr: threshold for residual-based masking (default 20000)
  • Np: size of local covariance matrix in pixels (default 33)
  • widx: width of boxcar window in x which determines size of region used for samples for the local covariance estimate (default 129)
  • widy: width of boxcar window in y which determines size of region used for samples for the local covariance estimate (default 129)
  • tilex: total number of tile divisions along x (default 1)
  • tiley: total number of tile divisions along y (default tilex)
  • ftype::Int: determine the Float precision, 32 is Float32, otherwise Float64
source