cosmicray_lacosmic

ccdproc.cosmicray_lacosmic(ccd, 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=u'meanmask', fsmode=u'median', psfmodel=u'gauss', psffwhm=2.5, psfsize=7, psfk=None, psfbeta=4.765, verbose=False)[source]

Identify cosmic rays through the lacosmic technique. The lacosmic technique identifies cosmic rays by identifying pixels based on a variation of the Laplacian edge detection. The algorithm is an implementation of the code describe in van Dokkum (2001) :ref:[R1] as implemented by McCully (2014) [R2]. If you use this algorithm, please cite these two works.

Parameters:

ccd: `~ccdproc.CCDData` or `~numpy.ndarray`

Data to have cosmic ray cleaned

sigclip : float, optional

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

sigfrac : float, 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.

objlim : float, 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.

pssl : float, 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.

gain : float, optional

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

readnoise : float, optional

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

satlevel : float, 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.

niter : int, optional

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

sepmed : boolean, 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:n ‘median’: An umasked 5x5 median filtern ‘medmask’: A masked 5x5 median filtern ‘meanmask’: A masked 5x5 mean filtern ‘idw’: A masked 5x5 inverse distance weighted interpolationn Default: “meanmask”.

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

Method to build the fine structure image:n ‘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”.

psffwhm : float, optional

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

psfsize : int, optional

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

psfk : float 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.

psfbeta : float, optional

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

verbose : boolean, optional

Print to the screen or not. Default: False.

{log}

Returns:

nccd : CCDData or ndarray

An object of the same type as ccd is returned. If it is a CCDData, the mask attribute will also be updated with areas identified with cosmic rays masked.

crmask : ndarray

If an ndarray is provided as ccd, a boolean ndarray with the cosmic rays identified will also be returned.

Notes

Implementation of the cosmic ray identification L.A.Cosmic: http://www.astro.yale.edu/dokkum/lacosmic/

References

[R1](1, 2) van Dokkum, P; 2001, “Cosmic-Ray Rejection by Laplacian Edge Detection”. The Publications of the Astronomical Society of the Pacific, Volume 113, Issue 789, pp. 1420-1427. doi: 10.1086/323894
[R2](1, 2) McCully, C., 2014, “Astro-SCRAPPY”, https://github.com/astropy/astroscrappy

Examples

  1. Given an numpy.ndarray object, the syntax for running cosmicrar_lacosmic would be:

    >>> newdata, mask = cosmicray_lacosmic(data, sigclip=5)  
    

    where the error is an array that is the same shape as data but includes the pixel error. This would return a data array, newdata, with the bad pixels replaced by the local median from a box of 11 pixels; and it would return a mask indicating the bad pixels.

  2. Given an CCDData object with an uncertainty frame, the syntax for running cosmicrar_lacosmic would be:

    >>> newccd = cosmicray_lacosmic(ccd, sigclip=5)   
    

    The newccd object will have bad pixels in its data array replace and the mask of the object will be created if it did not previously exist or be updated with the detected cosmic rays.