detect_cosmics

astroscrappy.detect_cosmics(indat, inmask=None, sigclip=4.5, sigfrac=0.3, objlim=5.0, gain=1.0, readnoise=6.5, satlevel=65536.0, pssl=0.0, niter=4, sepmed=True, cleantype='meanmask', fsmode='median', psfmodel='gauss', psffwhm=2.5, psfsize=7, psfk=None, psfbeta=4.765, verbose=False)

Detect cosmic rays in a numpy array.

If you use this code, please add this repository address in a footnote: https://github.com/astropy/astroscrappy

Please cite the original paper which can be found at: http://www.astro.yale.edu/dokkum/lacosmic/

van Dokkum 2001, PASP, 113, 789, 1420 (article : http://adsabs.harvard.edu/abs/2001PASP..113.1420V)

Parameters
indatfloat numpy array

Input data array that will be used for cosmic ray detection.

inmaskboolean numpy array, optional

Input bad pixel mask. Values of True will be ignored in the cosmic ray detection/cleaning process. Default: None.

sigclipfloat, optional

Laplacian-to-noise limit for cosmic ray detection. Lower values will flag more pixels as cosmic rays. Default: 4.5.

sigfracfloat, optional

Fractional detection limit for neighboring pixels. For cosmic ray neighbor pixels, a lapacian-to-noise detection limit of sigfrac * sigclip will be used. Default: 0.3.

objlimfloat, optional

Minimum contrast between Laplacian image and the fine structure image. Increase this value if cores of bright stars are flagged as cosmic rays. Default: 5.0.

psslfloat, optional

Previously subtracted sky level in ADU. We always need to work in electrons for cosmic ray detection, so we need to know the sky level that has been subtracted so we can add it back in. Default: 0.0.

gainfloat, optional

Gain of the image (electrons / ADU). We always need to work in electrons for cosmic ray detection. Default: 1.0

readnoisefloat, optional

Read noise of the image (electrons). Used to generate the noise model of the image. Default: 6.5.

satlevelfloat, optional

Saturation of level of the image (electrons). This value is used to detect saturated stars and pixels at or above this level are added to the mask. Default: 65536.0.

niterint, optional

Number of iterations of the LA Cosmic algorithm to perform. Default: 4.

sepmedboolean, optional

Use the separable median filter instead of the full median filter. The separable median is not identical to the full median filter, but they are approximately the same and the separable median filter is significantly faster and still detects cosmic rays well. Default: True

cleantype{‘median’, ‘medmask’, ‘meanmask’, ‘idw’}, optional

Set which clean algorithm is used:

‘median’: An umasked 5x5 median filter

‘medmask’: A masked 5x5 median filter

‘meanmask’: A masked 5x5 mean filter

‘idw’: A masked 5x5 inverse distance weighted interpolation

Default: “meanmask”.

fsmode{‘median’, ‘convolve’}, optional

Method to build the fine structure image:

‘median’: Use the median filter in the standard LA Cosmic algorithm ‘convolve’: Convolve the image with the psf kernel to calculate the fine structure image. Default: ‘median’.

psfmodel{‘gauss’, ‘gaussx’, ‘gaussy’, ‘moffat’}, optional

Model to use to generate the psf kernel if fsmode == ‘convolve’ and psfk is None. The current choices are Gaussian and Moffat profiles. ‘gauss’ and ‘moffat’ produce circular PSF kernels. The ‘gaussx’ and ‘gaussy’ produce Gaussian kernels in the x and y directions respectively. Default: “gauss”.

psffwhmfloat, optional

Full Width Half Maximum of the PSF to use to generate the kernel. Default: 2.5.

psfsizeint, optional

Size of the kernel to calculate. Returned kernel will have size psfsize x psfsize. psfsize should be odd. Default: 7.

psfkfloat numpy array, optional

PSF kernel array to use for the fine structure image if fsmode == ‘convolve’. If None and fsmode == ‘convolve’, we calculate the psf kernel using ‘psfmodel’. Default: None.

psfbetafloat, optional

Moffat beta parameter. Only used if fsmode==’convolve’ and psfmodel==’moffat’. Default: 4.765.

verboseboolean, optional

Print to the screen or not. Default: False.

Returns
crmaskboolean numpy array

The cosmic ray mask (boolean) array with values of True where there are cosmic ray detections.

cleanarrfloat numpy array

The cleaned data array.

Notes

To reproduce the most similar behavior to the original LA Cosmic (written in IRAF), set inmask = None, satlevel = np.inf, sepmed=False, cleantype=’medmask’, and fsmode=’median’.

The original IRAF version distinguishes between spectroscopic and imaging data. This version does not. After sky subtracting the spectroscopic data, this version will work well. The 1-d ‘gaussx’ and ‘gaussy’ values for psfmodel can also be used for spectroscopic data (and may even alleviate the need to do sky subtraction, but this still requires more testing).