potpyri.primitives.image_procs
Image calibration, masking, error arrays, alignment, and stacking.
Provides WCS handling, alignment to a common grid, cosmic-ray rejection, satellite masking, master calibration combination, and stacked science output. Authors: Charlie Kilpatrick.
Functions
|
Update stack mask from per-image masks and saturation; set masked pixels to 0. |
|
Solve WCS and align reduced images to a common grid. |
|
Compute relative flux scale factors from SExtractor catalogs and RA/Dec cross-matching. |
|
Compute per-pixel error array from science image, mask, and read noise. |
|
Build a bad-pixel/saturation/cosmic-ray mask using LACosmic and optional growth. |
|
Remove row and column median from stack science data (detrend). |
|
Build a TAN WCS for a given field center and output size. |
|
Compute mean RA/Dec of image centers from a list of FITS files. |
|
Full image processing: align, mask, stack, and optionally detrend. |
|
Detect and mask satellite trails in science images using acstools.satdet. |
|
Remove PV* distortion keywords from a FITS header (modifies in place). |
|
Combine aligned CCDData list with exposure scaling (median or average). |
- potpyri.primitives.image_procs.add_stack_mask(stack, stacking_data)[source]
Update stack mask from per-image masks and saturation; set masked pixels to 0.
- Parameters:
stack (
astropy.io.fits.HDUList) – HDU list with [0]=SCI, [1]=MASK, [2]=ERROR. Modified in place.stacking_data (
listofccdproc.CCDData) – Per-image data (headers used for SATURATE).
- Returns:
The same stack reference (modified in place).
- Return type:
astropy.io.fits.HDUList
- potpyri.primitives.image_procs.align_images(reduced_files, paths, tel, binn, use_wcs=None, fieldcenter=None, out_size=None, skip_gaia=False, keep_all_astro=False, log=None)[source]
Solve WCS and align reduced images to a common grid.
Runs astrometry.net and optionally Gaia alignment, then reprojects to a common WCS.
- Parameters:
reduced_files (
listofstr) – Paths to reduced FITS files.paths (
dict) – Paths dict from options.add_paths.tel (
Instrument) – Instrument instance.binn (
str) – Binning string.use_wcs (
astropy.wcs.WCS, optional) – If provided, use this WCS instead of solving.fieldcenter (
sequence, optional) – [ra, dec] for generating WCS.out_size (
int, optional) – Output image size.skip_gaia (
bool, optional) – If True, skip Gaia alignment. Default is False.keep_all_astro (
bool, optional) – If True, keep all images regardless of astrometric dispersion.log (
ColoredLogger, optional) – Logger for progress.
- Returns:
(aligned_data, masks, errors) as lists, or None if no images aligned.
- Return type:
tupleorNone
- potpyri.primitives.image_procs.compute_relative_scales(data_images, paths, exptimes, log=None, match_radius_arcsec=2.0, min_sources=5, min_rel_err=0.01)[source]
Compute relative flux scale factors from SExtractor catalogs and RA/Dec cross-matching.
Runs Source Extractor on each aligned science image, cross-matches sources between frames in RA/Dec, and uses the inverse-variance-weighted mean of flux ratios (using FLUXERR_AUTO) to define scale factors so that combine() can stack on a common flux scale. For frames where matching fails (too few matches), the scale is set to 1.0/exptime for that frame.
- Parameters:
data_images (
listofstr) – Paths to FITS with SCI extension (e.g. *_data.fits from image_proc).paths (
dict) – Paths dict; paths.get(‘source_extractor’) used for SExtractor binary.exptimes (
array-like) – Exposure time in seconds for each image (same order as data_images).log (
ColoredLogger, optional) – Logger for progress.match_radius_arcsec (
float, optional) – Match radius in arcsec for cross-matching sources. Default 2.0.min_sources (
int, optional) – Minimum matched sources per frame to compute scale. Default 5.min_rel_err (
float, optional) – Minimum relative flux error (err/flux) used in variance of ratio to avoid zero variance. Default 0.01.
- Returns:
Relative scale factors, one per image (reference frame scale = 1.0). Frames that fail get scale = 1.0/exptime. None if catalogs cannot be built at all (caller should use exposure-time-only scaling).
- Return type:
np.ndarrayorNone
- potpyri.primitives.image_procs.create_error(science_data, mask_data, rdnoise)[source]
Compute per-pixel error array from science image, mask, and read noise.
Error = sqrt(poisson + rms^2 + rdnoise^2). Masked pixels are excluded from RMS estimation. science_data can be a path or HDU; mask_data is an HDU with mask array.
- Parameters:
science_data (
strorastropy.io.fits.HDUList) – Path to science FITS or open HDU list.mask_data (
astropy.io.fits.PrimaryHDUorImageHDU) – HDU containing boolean or integer mask (masked pixels excluded from rms).rdnoise (
float) – Read noise in electrons.
- Returns:
Error array in electrons (BUNIT=’ELECTRONS’).
- Return type:
astropy.io.fits.PrimaryHDU
- potpyri.primitives.image_procs.create_mask(science_data, saturation, rdnoise, sigclip=3.0, sigfrac=0.1, objlim=4.0, niter=6, outpath='', grow=0, cosmic_ray=True, fsmode='convolve', cleantype='medmask', log=None)[source]
Build a bad-pixel/saturation/cosmic-ray mask using LACosmic and optional growth.
- Parameters:
science_data (
str) – Path to science FITS file.saturation (
float) – Saturation level (ADU); pixels above this get flag 4.rdnoise (
float) – Read noise for LACosmic (electrons).sigclip (
float, optional) – LACosmic sigma clipping. Default is 3.0.sigfrac (
float, optional) – LACosmic sigfrac. Default is 0.10.objlim (
float, optional) – LACosmic objlim. Default is 4.0.niter (
int, optional) – LACosmic iterations. Default is 6.outpath (
str, optional) – Output path for debug files. Default is ‘’.grow (
int, optional) – Pixels to grow mask (adds flag 8). Default is 0.cosmic_ray (
bool, optional) – If True, run LACosmic. Default is True.fsmode (
str, optional) – LACosmic fsmode. Default is ‘convolve’.cleantype (
str, optional) – LACosmic cleantype. Default is ‘medmask’.log (
ColoredLogger, optional) – Logger for progress.
- Returns:
(cleaned_data_ndarray, mask_ImageHDU). Mask uses additive flags: 1=bad, 2=cosmic ray, 4=saturated, 8=neighbor.
- Return type:
tuple
- potpyri.primitives.image_procs.detrend_stack(stack)[source]
Remove row and column median from stack science data (detrend).
- Parameters:
stack (
astropy.io.fits.HDUList) – HDU list with [0]=SCI, [1]=MASK. Modified in place.- Returns:
The same stack reference (modified in place).
- Return type:
astropy.io.fits.HDUList
- potpyri.primitives.image_procs.generate_wcs(tel, binn, fieldcenter, out_size)[source]
Build a TAN WCS for a given field center and output size.
- Parameters:
tel (
Instrument) – Instrument instance (for pixel scale).binn (
str) – Binning string (e.g. ‘22’).fieldcenter (
sequence) – [ra_deg, dec_deg] or parseable coordinate.out_size (
int) – Output image size (out_size x out_size pixels).
- Returns:
WCS instance.
- Return type:
astropy.wcs.WCS
- potpyri.primitives.image_procs.get_fieldcenter(images)[source]
Compute mean RA/Dec of image centers from a list of FITS files.
- Parameters:
images (
listofstr) – Paths to FITS files with WCS in primary header.- Returns:
[mean_ra_deg, mean_dec_deg].
- Return type:
list
- potpyri.primitives.image_procs.image_proc(image_data, tel, paths, skip_skysub=False, fieldcenter=None, out_size=None, satellites=True, cosmic_ray=True, skip_gaia=False, keep_all_astro=False, relative_calibration=False, log=None)[source]
Full image processing: align, mask, stack, and optionally detrend.
Orchestrates WCS solving, alignment, satellite/cosmic-ray masking, stacking, and optional sky subtraction. Optionally calibrates frames to each other using SExtractor catalogs and RA/Dec cross-matching before stacking.
- Parameters:
image_data (
astropy.table.Table) – File table subset (e.g. one TargType) with ‘File’ column.tel (
Instrument) – Instrument instance.paths (
dict) – Paths dict from options.add_paths.skip_skysub (
bool, optional) – If True, skip sky subtraction. Default is False.fieldcenter (
sequence, optional) – [ra, dec] for alignment.out_size (
int, optional) – Output image size.satellites (
bool, optional) – If True, mask satellite trails. Default is True.cosmic_ray (
bool, optional) – If True, run cosmic-ray rejection. Default is True.skip_gaia (
bool, optional) – If True, skip Gaia alignment. Default is False.keep_all_astro (
bool, optional) – If True, keep all images regardless of astrometric dispersion.relative_calibration (
bool, optional) – If True, run SExtractor on each aligned frame, cross-match sources in RA/Dec, and use relative fluxes to set scale factors for the stack combine. Default is False.log (
ColoredLogger, optional) – Logger for progress.
- Returns:
Path to stacked output FITS, or None if processing failed.
- Return type:
strorNone
- potpyri.primitives.image_procs.mask_satellites(images, filenames, log=None)[source]
Detect and mask satellite trails in science images using acstools.satdet.
- Parameters:
images (
listofccdproc.CCDData) – Science images; data is modified in place (trail pixels set to nan).filenames (
listofstr) – Paths corresponding to images (used for temp files).log (
ColoredLogger, optional) – Logger for progress.
- Return type:
None
- potpyri.primitives.image_procs.remove_pv_distortion(header)[source]
Remove PV* distortion keywords from a FITS header (modifies in place).
- Parameters:
header (
astropy.io.fits.Header) – Header to modify.- Returns:
The same header reference (modified in place).
- Return type:
astropy.io.fits.Header
- potpyri.primitives.image_procs.stack_data(stacking_data, tel, masks, errors, mem_limit=8000000000.0, log=None, relative_scales=None)[source]
Combine aligned CCDData list with exposure scaling (median or average).
Masks are applied (masked pixels set to nan) before combination. Optionally uses relative flux scales from source cross-matching so that frames are combined on a common flux scale.
- Parameters:
stacking_data (
listofccdproc.CCDData) – Aligned science images.tel (
Instrument) – Instrument instance (for stack_method and get_exptime).masks (
listofastropy.io.fits.ImageHDU) – Per-image masks.errors (
listofastropy.io.fits.ImageHDU) – Per-image error arrays.mem_limit (
float, optional) – Memory limit in bytes for ccdproc.combine. Default is 8e9.log (
ColoredLogger, optional) – Logger for progress.relative_scales (
np.ndarray, optional) – Relative flux scale factors (one per image; reference = 1.0). If provided, scale = relative_scales (exposure-time scaling is not applied).
- Returns:
HDU list [science, mask, error] (mask/error from first image).
- Return type:
astropy.io.fits.HDUList