piXedfit_images

piXedfit.piXedfit_images.EBV_foreground_dust(ra, dec)

Function for estimating E(B-V) dust attenuation due to the foreground Galactic dust attenuation at a given coordinate on the sky.

Parameters:
  • ra – Right ascension coordinate in degree.

  • dec – Declination coordinate in degree.

Returns ebv:

E(B-V) value.

piXedfit.piXedfit_images.calc_pixsize(fits_image)

Function to get pixel size of an image

Parameters:

fits_image – Input image.

Returns pixsize_arcsec:

Pixel size in arc second.

piXedfit.piXedfit_images.central_brightest_pixel(flux_maps_fits, filter_id, xrange=None, yrange=None)

Function to get the central brightest pixel within a map of flux in a particular band.

Parameters:
  • flux_maps_fits – Input FITS file of the flux maps produced from the image processing.

  • filter_id – The filter of the image where the central pixel is to be found.

  • xrange – Range in x-axis for searching area.

  • yrange – Range in y-axis for searching area.

piXedfit.piXedfit_images.convert_flux_unit(wave, flux, init_unit='erg/s/cm2/A', final_unit='Jy')

Function for converting flux unit

Parameters:
  • wave – Wavelength grids.

  • flux – Flux grids.

  • init_unit – Initial unit of the input fluxes. Options are: ‘erg/s/cm2/A’, ‘erg/s/cm2’, ‘Jy’, ‘uJy’.

  • final_unit – Final unit for conversion. Options are: ‘erg/s/cm2/A’, ‘erg/s/cm2’, ‘Jy’, ‘uJy’.

piXedfit.piXedfit_images.create_psf_matching_kernel(init_PSF_name, target_PSF_name, pixscale_init_PSF, pixscale_target_PSF, window='top_hat', window_arg=1.0)

A function for creating convolution kernel for PSF matching given initial and target PSFs.

Parameters:
  • init_PSF_name – Image of input/initial PSF.

  • target_PSF_name – Image of target PSF.

  • pixscale_init_PSF – Pixel size of the initial PSF in arcsec.

  • pixscale_target_PSF – Pixel size of the target PSF in arcsec.

  • window – Options are ‘top_hat’ and ‘cosine_bell’.

  • window_arg – Coefficient value of the window function, following Photutils.

  • kernel – The data of convolution kernel in 2D array.

piXedfit.piXedfit_images.crop_stars_galregion_fits(input_fits, x_cent, y_cent, radius, name_out_fits=None)

Function for cropping foreground stars within a galaxy’s region of interst.

Parameters:
  • input_fits – The input FITS file containing the data cube that is produced using the flux_map() funciton in the images_processing class.

  • x_cent – 1D array of x coordinates of the center of the stars.

  • y_cent – 1D array of y coordinates of the center of the stars.

  • radius – 1D array of the estimated circular radii of the stars.

  • name_out_fits – Desired name for the output FITS file. If None, a generic name will be made.

piXedfit.piXedfit_images.draw_aperture_on_maps_fluxes(flux_maps_fits, ncols=6, e=[0.0], pa=[45.0], cent_x=None, cent_y=None, radius=[5.0], colors=None, lw=3, savefig=True, name_plot=None)

Function for drawing apertures on top of the multiband fluxes maps. This function can plot more than one aperture.

Parameters:
  • flux_maps_fits – input FITS file containing the maps of fluxes produced from the image processing.

  • e – Ellipticity of the apertures. Input in list format if want to make multiple apertures.

  • pa – Position angle of the elliptical apertures. Input in list format if want to make multiple apertures.

  • cent_x – The x coordinate of the central pixel. Input in list format if want to make multiple apertures.

  • cent_y – The y coordinate of the central pixel. Input in list format if want to make multiple apertures.

  • radius – The radius of the aperture. Input in list format if want to make multiple apertures.

  • colors – Colors of the apertures. Input in list format if want to make multiple apertures.

  • lw – Line width of the apertures in the plot.

  • savefig – Decide whether to save the plot or not.

  • name_plot – Name for the plot.

piXedfit.piXedfit_images.get_pixels_SED_fluxmap(flux_maps_fits, pix_x=None, pix_y=None, all_pixels=False)

Function to get SEDs of pixels.

Parameters:
  • flux_maps_fits – Input FITS file of the multiband fluxes.

  • pix_x – One dimensional array of x coordinates. Only relevant if all_pixels=False.

  • pix_y – One dimensional array of the y coordinates. Only relevant if all_pixels=False.

  • all_pixels – (optional) An option to get SEDs of all pixels.

Returns pix_x:

x coordinates.

Returns pix_y:

y coordinates.

Returns pix_SED_flux:

Fluxes of the pixels: (npixs,nbands)

Returns pix_SED_flux_err:

Flux uncertainties of the pixels: (npixels,nbands)

Returns photo_wave:

Central wavelength of the filters.

piXedfit.piXedfit_images.get_total_SED(flux_maps_fits)

Function to calculate total (i.e., integrated) SED from input maps of multiband fluxes.

Parameters:

flux_maps_fits – Input FITS file containing maps of fluxes produced from the image processing.

Returns tot_SED_flux:

Total fluxes in multiple bands. The format is 1D array.

Returns tot_SED_flux_err:

Total flux uncertainties in multiple bands. The format is 1D array.

Returns photo_wave:

The central wavelength of the filters.

class piXedfit.piXedfit_images.images_processing(filters, sci_img, var_img, gal_ra, gal_dec, dir_images=None, img_unit=None, img_scale=None, img_pixsizes=None, run_image_processing=True, flag_psfmatch=0, flag_reproject=0, flag_crop=0, kernels=None, gal_z=None, stamp_size=[101, 101], remove_files=True, idfil_align=None)

A Python class for processing multiband imaging data and producing a data cube containing maps of multiband fluxes that are matched in spatial resolution and sampling. The image processing basically includes the PSF matching to homogenize the spatial resolution of the multiband imaging data and spatial resampling and reprojection to homogenize the spatial sampling (i.e., pixel size) and spatil reprojection of the mulltiband imaging data. A list of imaging data sets that can be handled automatically by this class in the current version of piXedfit can be seen at List of imaging data. However, one need to download convolution kernels from this link and put those inside data/kernels/ within the piXedfit directory ($PIXEDFIT_HOME/data/kernels). These kernels are not included in the piXedfit repository because of the large file sizes. This class can also handle other imaging data. For this, one need to provide kernels for PSF matching. These kernel images should have the same pixel size as the corresponding input images.

Parameters:
  • filters – List of photometric filters names. To check the filters currently available and manage them (e.g., adding new ones) please see this page. For this input, it is not mandatory to make the input filters in a wavelength order.

  • sci_img – Dictionary containing names of the science images. An example of input: sci_img={‘sdss_u’:’img1.fits’, ‘sdss_g’:’img2.fits’}

  • var_img – Dictionary containing names of the variance images. It has a similar format to sci_img.

  • gal_ra – Right Ascension (RA) coordinate of the target galaxy. This should be in degree.

  • gal_dec – Declination (DEC) coordinate of the target galaxy. This should be in degree.

  • dir_images – Directory where images are stored.

  • img_unit – (optional) Unit of pixel value in the multiband images. The acceptable format of this input is a python dictionary, similar to that of sci_img. This input is optional. This input will only be considered (and required) if the input images are not among the default list of recognized imaging data in piXedfit (i.e. GALEX, SDSS, 2MASS, WISE, Spitzer, and Herschel). The allowed units are: (1)”erg/s/cm2/A”, (2) “Jy”, and (3) “MJy/sr”.

  • img_scale – (optional) Scale of the pixel value with respect to the unit in img_unit. For instance, if image is in unit of MJy, the img_unit can be set to be “Jy” and img_scale is set to be 1e+6. This input is only relevant if the input images are not among the default list of recognized images in piXedfit. The format of this input should be in python dictionary, similar to sci_img.

  • img_pixsizes – (optional) Pixel sizes (in arcsecond) of the input imaging data. This input should be in dictionary format, similar to sci_img. If not provided, pixel size will be calculated based on the WCS information in the header of the FITS file.

  • flag_psfmatch – (optional) Flag stating whether the multiband imaging data have been PSF-matched or not. The options are: (1) 0 means hasn’t been PSF-matched, and (2)1 means has been PSF-matched.

  • flag_reproject – (optional) Flag stating whether the multiband imaging data have been spatially-resampled and matched in the projection. The options are: (1)0 means not yet, and (2)1 means has been carried out.

  • flag_crop – (optional) Flag stating whether the multiband imaging data have been cropped around the target galaxy. The options are: (1)0 means not yet, and (2)1 means has been cropped. If flag_crop=0, cropping will be done according to the input stamp_size. If flag_crop=1, cropping will not be done.

  • kernels – (optional) Dictionary containing names of FITS files for the kernels to be used for the PSF matching process. If None, the code will first look for kernels inside piXedfit/data/kernels directory. If an appropriate kernel is not found, the kernel input remain None and PSF matching will not be done for the corresponding image. If external kerenels are avaiable, the kernel images should have the same pixel size as the corresponding input images and the this input should be in a dictionary format, which is similar to the input sci_img and var_img.

  • gal_z – Galaxy’s redshift. This information will not be used in the image processing and only intended to be saved in the heder of a produced FITS file.

  • stamp_size – Desired size for the reduced maps of multiband fluxes. This is a list data type with 2 elements. Accepted struture is: [dim_y,dim_x]. Only relevant if flag_crop=0.

  • remove_files – If True, the unnecessary image files produced during the image processing will be removed. This can save disk space. If False, those files will not be removed.

  • idfil_align – Index of the filter of which the image will be used as the reference in the spatial reprojection and sampling processes.

flux_map(gal_region, output_stamps=None, Gal_EBV=None, scale_unit=1e-17, mag_zp_2mass=None, unit_spire='Jy_per_beam', name_out_fits=None)

Function for calculating maps of multiband fluxes from the stamp images produced by the reduced_stamps().

Parameters:
  • gal_region – A 2D array containing the galaxy’s region of interest. It is preferably the output of the gal_region(), but one can also make this input region. The 2D array should has the same size as that of the output stamps and the pixel value is 1 for the galaxy’s region and 0 otherwise.

  • output_stamps – (optional) Supply output_stamps dictionary input. This input is optional. If run_image_processing=False, this input is mandatory.

  • Gal_EBV – The E(B-V) dust attenuation due to the foreground Galactic dust. This is optional parameter. If None, this value will be retrive from the IRSA data server through the astroquery package.

  • scale_unit – Normalized unit for the fluxes in the output fits file. The unit is flux density in erg/s/cm^2/Ang.

  • mag_zp_2mass – Magnitude zero-points of 2MASS images. Shoud be in 1D array with three elements: [magzp-j,magzp-h,magzp-k]. This is optional parameter. If not given (i.e. None), the values will be taken from the header of the FITS files.

  • unit_spire – Unit of SPIRE images, in case Herschel/SPIRE image is included in the analysis. Therefore, this input is only relevant if Herschel/SPIRE image is among the images that are analyzed. Options are: [‘Jy_per_beam’, ‘MJy_per_sr’, ‘Jy_per_pixel’]

  • name_out_fits – Desired name for the output FITS file. If None, a default name will be adopted.

galaxy_region(segm_maps_ids=None, use_ellipse=False, x_cent=None, y_cent=None, ell=0, pa=45.0, radius_sma=30.0)

Define galaxy’s region of interest for further analysis.

Parameters:
  • segm_maps_ids – Array of IDs of selected segmentation maps (i.e., filters) to be used (merged) for constructing the galaxy’s region. If None, all segmentation maps are merged together.

  • use_ellipse – Alternative of defining galaxy’s region using elliptical aperture centered at the target galaxy. Set use_ellipse=True if you want to use this option.

  • x_cent – x coordinate of the ellipse center. If x_cent=None, the ellipse center is assumed to be the same as the image center.

  • y_cent – y coordinate of the ellipse center. If y_cent=None, the ellipse center is assumed to be the same as the image center.

  • ell – Ellipticity of the elliptical aperture.

  • pa – Position angle of the elliptical aperture.

  • radius_sma – Radal distance along the semi-major axis of the elliptical aperture. This radius is in pixel unit.

Returns gal_region:

Output galaxy’s region of interest.

get_output_stamps()

Get the names of output stamp images in a dictionary format.

plot_gal_region(gal_region, output_stamps=None, ncols=6, savefig=True, name_plot=None)

Plot the defined galaxy’s region.

Parameters:
  • gal_region – The defined galaxy’s region.

  • output_stamps – (optional) Supply output_stamps dictionary input. This input is optional. If run_image_processing=False, this input is mandatory.

  • ncols – Number of columns

  • savefig – (default=True) Decide to save the plot or not.

  • name_plot – (optional) Name for the output plot.

plot_image_stamps(output_stamps=None, ncols=6, savefig=True, name_plot_sci=None, name_plot_var=None)

Plotting resulted image stamps from the image processing.

Parameters:
  • output_stamps – (optional) Supply output_stamps dictionary input. This input is optional. If run_image_processing=False, this input is mandatory.

  • ncols – Number of columns in the multipanel plots.

  • savefig – (default=True) Decide to save the plot or not.

  • name_plot_sci – (optional) Name for the output plot of science image stamps.

  • name_plot_var – (optional) Name for the output plot of variance image stamps.

plot_segm_maps(ncols=6, savefig=True, name_plot=None)

Plotting segmentation maps.

Parameters:
  • ncols – Number of columns in the multipanel plots.

  • savefig – (default=True) Decide to save the plot or not.

  • name_plot – (optional) Name for the output plot.

rectangular_regions(output_stamps=None, x=None, y=None, ra=None, dec=None, make_plot=True, ncols=6, savefig=True, name_plot=None)

Define rectangular aperture and get pixels within it

Parameters:
  • output_stamps – (optional) Supply output_stamps dictionary input. This input is optional. If run_image_processing=False, this input is mandatory.

  • x – x coordinates of the rectangular corners. The shape should be (n,4) with n is the number of rectangular apertures.

  • y – y coordinates of the rectangular corners. The shape should be (n,4) with n is the number of rectangular apertures.

  • ra – RA coordinates of the rectangular corners. The shape should be (n,4) with n is the number of rectangular apertures. If x input is given, this input will be ignored.

  • dec – DEC coordinates of the rectangular corners. The shape should be (n,4) with n is the number of rectangular apertures. If y input is given, this inpur will be ignored.

  • pa – Position angle of the rectangular aperture.

  • make_plot – Decide to make a plot or not.

  • ncols – Number of columns

  • savefig – (default=True) Decide to save the plot or not.

  • name_plot – (optional) Name for the output plot.

Returns rect_region:

Rectangular aperture region.

reduced_stamps()

Run the image processing that includes PSF matching, spatial resampling and reprojection, and cropping around the target galaxy.

segmentation_sep(output_stamps=None, thresh=1.5, minarea=30, deblend_nthresh=32, deblend_cont=0.005)

Get segmentation maps of a galaxy in multiple bands using the SEP (a Python version of the SExtractor).

Parameters:
  • output_stamps – (optional) Supply output_stamps dictionary input. This input is optional. If run_image_processing=False, this input is mandatory.

  • thresh – Detection threshold for the source detection and segmentation.

  • minarea – Minimum number of pixels (above threshold) required for an object to be detected.

  • deblend_nthresh – Number of deblending sub-thresholds. Default is 32.

  • deblend_cont – Minimum contrast ratio used for object deblending. Default is 0.005. To entirely disable deblending, set to 1.0.

piXedfit.piXedfit_images.kpc_per_pixel(z, arcsec_per_pix, cosmo='flat_LCDM', H0=70.0, Om0=0.3)

Function for calculating a physical scale (in kpc) corresponding to a single pixel in an image.

Parameters:
  • z – Redshift of the galaxy.

  • arsec_per_pix – Pixel size in arcsecond.

  • cosmo – Choices for the cosmology. The options are: (a)’flat_LCDM’ or 0, (b)’WMAP5’ or 1, (c)’WMAP7’ or 2, (d)’WMAP9’ or 3, (e)’Planck13’ or 4, (f)’Planck15’ or 5, (g)’Planck18’ or 6. These are similar to the choices available in the Astropy Cosmology package.

  • H0 – The Hubble constant at z=0. Only relevant if cosmo=’flat_LCDM’.

  • Om0 – The Omega matter at z=0.0. Only relevant if cosmo=’flat_LCDM’.

Returns kpc_per_pix:

Output physical scale in kpc.

piXedfit.piXedfit_images.open_fluxmap_fits(flux_maps_fits)

Function to extract data from the FITS file of fluxes maps produced from the image processing.

Parameters:

flux_maps_fits – Input FITS file containing the flux maps.

Returns filters:

The set of filters.

Returns gal_region:

The galaxy’s spatial region of interest.

Returns flux_map:

The data cube of flux maps: (nbands,ny,nx).

Returns flux_err_map:

The data cube of flux uncertainties maps: (nbands,ny,nx).

Returns unit_flux:

The flux unit. In this case, the scale factor to erg/s/cm^2/Angstrom.

piXedfit.piXedfit_images.photometry_within_aperture(flux_maps_fits, e=0.0, pa=45.0, cent_x=None, cent_y=None, radius=None)

Function to get photometry within a given aperture (elliptical/circular) from input flux maps

Parameters:
  • flux_maps_fits – input FITS file containing the maps of fluxes produced from the image processing.

  • e – Ellipticity of the apertures.

  • pa – Position angle of the elliptical apertures.

  • cent_x – The x coordinate of the central pixel.

  • cent_y – The y coordinate of the central pixel.

  • radius – The radius of the aperture.

Returns tot_fluxes:

Output total fluxes.

Returns tot_flux_errors:

Output total flux uncertainties.

piXedfit.piXedfit_images.plot_SED_pixels(flux_maps_fits, pix_x=None, pix_y=None, all_pixels=False, logscale_y=True, logscale_x=True, wunit='angstrom', yrange=None, xrange=None, savefig=True, name_plot=None)

Function for plotting SEDs of pixels.

Parameters:
  • flux_maps_fits – Input FITS file of the multiband fluxes.

  • pix_x – One dimensional array of x coordinates. Only relevant if all_pixels=False.

  • pix_y – One dimensional array of the y coordinates. Only relevant if all_pixels=False.

  • all_pixels – (optional) An option to get SEDs of all pixels.

  • logscale_y – (optional) Option to set the y axis in logarithmic scale.

  • logscale_x – (optional) Option to set the x axis in logarithmic scale.

  • wunit – (optional) Wavelength unit. Options are ‘angstrom’ and ‘micron’.

  • yrange – (optional) Range in y axis.

  • xrange – (optional) Range in x axis.

  • savefig – (optional) Option to save the plot.

  • name_plot – (optional) Name for the output plot.

piXedfit.piXedfit_images.plot_SNR_radial_profile(flux_maps_fits, e=0.0, pa=45.0, cent_x=None, cent_y=None, yrange=None, xrange=None, savefig=True, name_plot=None)

Function for plotting S/N ratios of pixels.

Parameters:
  • flux_maps_fits – Input FITS file produced from the image processing.

  • e – Ellipticity of the elliptical apertures that will be used for deriving the S/N ratios.

  • pa – The position angle of the elliptical apertures.

  • cent_x – The x coordinate of the central pixel.

  • cent_y – The y coordinate of the central pixel.

  • yrange – Range in y-axis for the plot.

  • xrange – Range in x-axis for the plot.

  • savefig – Decide whether to save the plot or not.

  • name_plot – Name for the output plot.

piXedfit.piXedfit_images.plot_maps_fluxes(flux_maps_fits, ncols=5, savefig=True, name_plot_mapflux=None, name_plot_mapfluxerr=None)

Function for plotting maps of multiband fluxes.

Parameters:
  • flux_maps_fits – Input FITS file of multiband flux maps.

  • ncols – Number of columns in the plots.

  • savefig – Decide whether to save the plot or not.

  • name_plot_mapflux – Name of the output plot for the maps of multiband fluxes.

  • name_plot_mapfluxerr – Name of the output plot for the maps of multiband flux uncertainties.

piXedfit.piXedfit_images.radial_profile_psf(psf, pixsize, e=0.0, pa=45.0, dr_arcsec=None)

Function to get the radial profile of PSF.

Parameters:
  • psf – Input PSF, either in a FITS file or the data array (2D).

  • e – (default = 0) Ellipticity of the apertures. The default is zero (i.e., circular).

  • pa – The position angle of the elliptical apertures.

  • dr_arsec – Radial increment in units of arc second.

Pixsize:

Pixel size in arc second.

Returns curve_rad:

The radius array of the PSF radial profile.

Returns curve_val:

The normalized fluxes of the PSF radial profile.

piXedfit.piXedfit_images.subtract_background(fits_image, hdu_idx=0, sigma=3.0, box_size=None, mask_region=None, mask_sources=True, var_image=None, thresh=1.5, minarea=5, deblend_nthresh=32, deblend_cont=0.005)

Function for estimating 2D background and subtracting it from the input image. This function also produce RMS image. This function adopts the Background2D function from the photutils package. To estimate 2D background, the input image is gridded and sigma clipping is done to each bin (grid). Then 2D interpolation is performed to construct the 2D background image. A more information can be seen at photutils website.

Parameters:
  • fits_image – Input image.

  • hdu_idx – The extension (HDU) where the image is stored. Default is 0 (HDU0).

  • sigma – Sigma clipping threshold value.

  • box_size – The box size for the image gridding. The format is: [ny, nx]. If None, both axes will be divided into 10 grids.

  • mask_region – Region within the image that are going to be excluded. mask_region should be 2D array with the same size as the input image.

  • mask_sources – If True, source detection and segmentation will be performed with SEP (Pyhton version of SExtractor) and the regions associated with the detected sources will be excluded. This help reducing contamination by astronomical sources.

  • var_image – Variance image (in FITS file format) to be used in the sources detection process. This input argument is only relevant if mask_sources=True.

  • thresh – Detection threshold for the sources detection. If variance image is supplied, the threshold value for a given pixel is interpreted as a multiplicative factor of the uncertainty (i.e. square root of the variance) on that pixel. If var_image=None, the threshold is taken to be 2.5 percentile of the pixel values in the image.

  • minarea – Minimum number of pixels above threshold triggering detection.

  • deblend_nthresh – Number of thresholds used for object deblending.

  • deblend_cont – Minimum contrast ratio used for object deblending. To entirely disable deblending, set to 1.0.

piXedfit.piXedfit_images.test_psfmatching_kernel(init_PSF_name, target_PSF_name, kernel, pixscale_init_PSF, pixscale_target_PSF, dr_arcsec=None, xrange_arcsec=[-0.5, 0.5], savefig=False, name_plot=None)

Function for testing the convolution kernel. This includes convolving the initial PSF with the kernel, comparing the radial profiles of the convolved initial PSF and the target PSF along with the original initial PSF.

Parameters:
  • init_PSF_name – Image of the initial PSF.

  • target_PSF_name – Image of the target PSF.

  • kernel – The convolution kernel data (2D array).

  • pixscale_init_PSF – Pixel size of the initial PSF in arcsec.

  • pixscale_target_PSF – Pixel size of the target PSF in arcsec.

  • dr_arsec – Radial increment for the radial profile in arcsec.

  • xrange_arcsec – Range of x-axis in arsec.

  • savefig – Decide whether to save the plot or not.

  • name_plot – Name for the output plot.

piXedfit.piXedfit_images.var_img_2MASS(sci_img, skyrms_img=None, skyrms_value=None, name_out_fits=None)

Function for constructing a variance image from an input 2MASS image. The estimation of flux uncertainty follows the information provided on the 2MASS website.

Parameters:
  • sci_img – Input background-subtracted science image.

  • skyrms_img – FITS file of the RMS background image. If skyrms_value=None, this parameter should be provided. The background subtraction and calculation of RMS image can be performed using the subtract_background() function.

  • skyrms_value – Scalar value of median or mean of the RMS background image. This input will be considered when skyrms_img=None.

  • name_out_fits – Desired name for the output FITS file.

piXedfit.piXedfit_images.var_img_GALEX(sci_img, skybg_img, filter_name, name_out_fits=None)

Function for calculating variance image from an input GALEX image

Parameters:
  • sci_img – Input GALEX science image (i.e., background subtracted). This type of image is provided in the GALEX website as indicated with “-intbgsub” (i.e., background subtracted intensity map).

  • skybg_img – Input sky background image .

  • filter_name – A string of filter name. Options are: ‘galex_fuv’ and ‘galex_nuv’.

  • name_out_fits – Desired name for the output variance image. If None, a generic name will be used.

piXedfit.piXedfit_images.var_img_WISE(sci_img, unc_img, filter_name, skyrms_img, name_out_fits=None)

Function for constructing variance image from an input WISE image. The uncertainty estimation follows the information from the WISE website

Parameters:
  • sci_img – Input background-subtracted science image.

  • unc_img – Input uncertainty image. This type of WISE image is provided in the IRSA website and indicated with ‘-unc-’ keyword.

  • filter_name – The filter name. Options are: ‘wise_w1’, ‘wise_w2’, ‘wise_w3’, and ‘wise_w4’

  • skyrms_img – Input RMS background image. This image is produced in the background subtraction with the subtract_background() function.

  • name_out_fits – (optional, default: None) Desired name for the output FITS file. If None, a generic name will be made.

piXedfit.piXedfit_images.var_img_from_unc_img(unc_image, name_out_fits=None)

Function for constructing a variance image from an input of uncertainty image. This function simply takes square of the uncertainty image and store it into a new FITS file while retaining the header information.

Parameters:
  • unc_img – Input uncertainty image.

  • name_out_fits – (optional, default: None) Name of output FITS file. If None, a generic name will be generated.

piXedfit.piXedfit_images.var_img_from_weight_img(wht_image, name_out_fits=None)

Function for constructing a variance image from an input weight (i.e., inverse variance) image. This funciton simply takes inverse of the weight image and store it into a new FITS file while retaining the header information.

Parameters:

wht_image – Input of weight image (i.e., inverse variance).

Returns name_out_fits:

(optional, default: None) Name of output FITS file. If None, a generic name will be made.

piXedfit.piXedfit_images.var_img_sdss(fits_image, filter_name, name_out_fits=None)

A function for constructing a variance image from an input SDSS image

Parameters:
  • fits_image – Input SDSS image (corrected frame type).

  • filter_name – A string of filter name. Options are: ‘sdss_u’, ‘sdss_g’, ‘sdss_r’, ‘sdss_i’, and ‘sdss_z’.

Returns name_out_fits:

Name of output FITS file. If None, a generic name will be used.