Imphot source-code documentation

The following documentation comes from the doc-strings of the imphot module and its public functions and classes.

imphot

The imphot module provides functions and scripts that can be used to determine the photometric parameters of individual MUSE images. These parameters include the width of the point-spread function, the scaling and zero offset of the flux in the image and the pointing error of the observation. The photometric parameters of a particular MUSE image are estimated by comparing that image to an HST image of the same region of sky. Two algorithms are provided:

  1. The fit_image_photometry() function uses the HST image as a good estimate of the flux distribution on the sky, and sees what observing conditions would be needed to reduce the resolution of the HST image, and change its scaling and offsets to match the MUSE observation.

    To efficiently process many images of a single area, a multi-process iterator is also provided that hands the job of calling fit_image_photometry() on multiple images to multiple processes. This is called FitImagePhotometryMP().

  2. The fit_star_photometry() function fits a Moffat PSF to a given star or a bright point source in both the MUSE image and the HST image. The PSF parameters are determined directly from the star fit in the MUSE image. The flux scaling, flux offset and pointing errors are determined by comparing the total fluxes, flux offsets and positions of the star fits in the two images.

    To efficiently process many images of a single area, a multi-process iterator is also provided that hands the job of calling fit_star_photometry() on multiple images to multiple processes. This is called FitStarPhotometryMP().

In regions that have one or more bright stars, the two methods can be used to cross check each other. For this purpose, the fit_image_photometry() function can be asked to restrict the area of the image that it fits to be the same area that is used to fit the star in fit_star_photometry(). This is requested by passing an optional star=[ra,dec,radius] argument. When both fit_image_photometry() and fit_star_photometry() are passed the same values for this star argument, they operate on exactly the same pixels of the images. A multiprocessing iterator is also provided for this case, called FitImageAndStarPhotometryMP(). For each MUSE image this returns the results from both methods.

Beware that some bright stars suffer from an effect that acts like saturation. It isn’t clear whether this is in the HST observations or the MUSE observations. The result is that the flux scale-factor that is found on these stars is higher than elsewhere in the images. Note that this is seen regardless of whether the fitting is performed by fit_image_photometry() or fit_star_photometry(). This implies that the problem is not an artefact of the fitting method. If fit_image_photometry() is applied to an image that contains a bright star like this, a poor fit is obtained for both the star and other objects in the image. The residual image shows that the star ends up being under-subtracted, while the dimmer objects in the image become over-subtracted. Again, this indicates that the calibration factor for the star is different from other sources in the image. In such cases it may be better to exclude the star from the fit. This can be done by using the optional regions=<filename> argument of fit_image_photometry() to specify a ds9 region file that excludes a small region around the star.

In fit_image_photometry(), ds9 region files can be used to mask out unwanted areas of an image or to selectively indicate which areas of the image to include in the fit. Note that in ds9, each region can be marked as inclusive or exclusive, using the Property menu of the dialog that appears when one double clicks on a region.

The main() functions of the fit_photometry, regrid_hst_to_muse and make_wideband_image scripts are also provided. For example, the fit_photometry script is comprised of the following few lines that just invoke the main function within the imphot module:

import sys
import imphot
imphot.fit_photometry_main(sys.argv)
class imphot.FittedValue(param=None, value=None, stdev=None, fixed=None)[source]

A parameter value returned by a least-squares fitting function.

Parameters:
paramlmfit.parameter.Parameter

The best-fit value of a parameter.

Attributes:
valuefloat

The best-fit value.

stdevfloat

The estimated 1-sigma uncertainty in the best-fit value.

fixedbool

True if the value of this parameter was fixed during the fit. False if the parameter was free to vary.

exception imphot.UserError(message)[source]

An exception that is raised when an error occurs that is due to user-specified inputs, rather than a programming error.

with_traceback()

Exception.with_traceback(tb) – set self.__traceback__ to tb and return self.

exception imphot.WorkerError(message, traceback)[source]

An exception class for reporting exceptions that occur in separate worker processes.

Parameters:
messagestr

A short message about the exception.

tracebacklist of str

The traceback of the exception, in the form of a list of text lines.

with_traceback()

Exception.with_traceback(tb) – set self.__traceback__ to tb and return self.

class imphot.ImphotArgumentParser(functions, *args, **kwargs)[source]

A parser of command-line arguments, that recognizes options that are used by fit_image_photometry() and/or fit_star_photometry().

For example, to create a parser that recognizes the optional arguments of both fit_image_photometry() and fit_star_photometry(), one would do the following:

parser = imphot.ImphotArgumentParser(imphot.FIT_IMAGE +
                               imphot.FIT_STAR, prog=argv[0])

To use this parser to process the command-line arguments of a program, one would then use the parse_args() function of argparse.ArgumentParser as follows:

options = parser.parse_args(argv[1:])

This returns an object that contains initialized variables for each of the arguments recognized by imphot.ImphotArgumentParser. From this one can obtain separate dictionaries of the keyword/value pairs of the arguments of the fit_image_photometry() and fit_star_photometry() functions by using the extract_function_args() function as follows:

imfit_kwargs = imphot.extract_function_args(options,
                           imphot.fit_image_photometry)
stfit_kwargs = imphot.extract_function_args(options,
                           imphot.fit_star_photometry))

These dictionaries can then be passed as keyword/value arguments to fit_image_photometry() and fit_star_photometry() as follows:

imfit = imphot.fit_image_photometry(hst, muse, **imfit_kwargs)
stfit = imphot.fit_star_photometry(hst, muse, **stfit_kwargs)
Parameters:
functionsint

The fitting functions that we want to obtain arguments for, expressed as a sum of one or both of imphot.FIT_IMAGE and imphot.FIT_STAR.

add_argument(dest, ..., name=value, ...) add_argument(option_string, option_string, ..., name=value, ...)
add_argument(option_string, option_string, ..., name=value, ...) None
error(message: string)

Prints a usage message incorporating the message to stderr and exits.

If you override this in a subclass, it should not return – it should either exit or raise an exception.

imphot.extract_function_args(options, function)[source]

Given a dictionary of key/value pairs or a Namespace object returned by argparse.ArgumentParser.parse_args(), extract the subset of argument values that are recognized by a specified function and return them in a dictionary of keyword/value pairs. The caller can then pass this dictionary to the specified function as keyword/value arguments.

Parameters:
optionsargparse.Namespace or dict

The return value of a call to argparse.ArgumentParser.parse_args().

functionfunction

The function whose arguments are needed.

Returns:
outdict

A dictionary of keyword/value pairs that can be passed to the specified function.

class imphot.FittedPhotometry(method, muse, fit_report, scale, bg, dx, dy, fwhm, beta, rms_error)[source]

The superclass of FittedImagePhotometry and FittedStarPhotometry which contains the fitted photometry parameters of a MUSE image, the identity of the original MUSE image, and a report from the least-squares fit.

Parameters:
methodstr

A short word that identifies the method used to fit the photometry.

musempdaf.obj.Image

The MUSE image that the fit was performed on.

fit_reportstr

A multi-line string that contains a verbose report about the final least-squares fit or fits.

scaleFittedValue

The best-fit value and error of the calibration scale factor, MUSE/HST.

bgFittedValue

The best-fit value and error of the calibration offset, MUSE-HST.

dxFittedValue

The best-fit value and error of the x-axis pointing offset, MUSE.x-HST.x.

dyFittedValue

The best-fit value and error of the y-axis pointing offset, MUSE.y-HST.y.

fwhmFittedValue

The best-fit value and error of the FWHM of the Moffat PSF.

betaFittedValue

The best-fit value and error of the beta parameter of the Moffat PSF.

rms_errorfloat

The root-mean square of the pixel residuals of the fit in the MUSE image, in the same units as the pixels of the original MUSE image.

Attributes:
namestr

The basename of the MUSE FITS without the .fits extension.

fit_reportstr

A printable report on the fit from the least-squares fitting function.

scaleFittedValue

The best-fit value and error of the calibration scale factor, MUSE/HST.

bgFittedValue

The best-fit value and error of the calibration offset, MUSE-HST.

dxFittedValue

The best-fit value and error of the x-axis pointing offset, MUSE.x-HST.x (arcsec). This is the distance that features in the HST image had to be moved to the right, in the direction of increasing x-axis pixel index, to line them up with the same features in the MUSE image. One way to correct the pointing error of the MUSE observation, is to divide ‘dx’ by the pixel size along the x-axis (usually 0.2 arcsec), then add the resulting pixel offset to the CRPIX1 header parameter.

dyFittedValue

The best-fit value and error of the y-axis pointing offset, MUSE.y-HST.y (arcsec). This is the distance that features in the HST image had to be moved upwards, in the direction of increasing y-axis pixel index, to line them up with the same features in the MUSE image. One way to correct the pointing error of the MUSE observation, is to divide ‘dy’ by the pixel size along the y-axis (usually 0.2 arcsec), then add the resulting pixel offset to the CRPIX2 header parameter.

draFittedValue

The right-ascension error (arcsec) that corresponds to the pointing error dx,dy. This is the angular distance that features in the HST image had to be moved towards increased right-ascension, to line them up with the same feaures in the MUSE image. One way to correct the pointing error of the MUSE observation is to subtract ‘dra’ from the CRVAL1 header value of the MUSE observation.

ddecFittedValue

The declination error (arcsec) that corresponds to the pointing error dx,dy. This is the angular distance that features in the HST image had to be moved towards increased declination, to line them up with the same feaures in the MUSE image. One way to correct the pointing error of the MUSE observation is to subtract ‘ddec’ from the CRVAL2 header value of the MUSE observation.

fwhmFittedValue

The best-fit value and error of the FWHM of the Moffat PSF.

betaFittedValue

The best-fit value and error of the beta parameter of the Moffat PSF.

rms_errorfloat

The root-mean square of the pixel residuals of the fit in the MUSE image, in the same units as the pixels of the original MUSE image.

summary(header=True)[source]

Return a string that summarizes the results, optionally listed under a 3-line heading of column descriptions.

Parameters:
headingbool

True to include 3 header lines above the values (default=True)

class imphot.HstFilterInfo(hst)[source]

An object that contains the filter characteristics of an HST image.

Parameters:
hstmpdaf.obj.Image or str

The HST image to be characterized, or the name of an HST filter, such as “F606W.

Attributes:
filter_namestr

The name of the filter.

abmag_zerofloat

The AB magnitude that corresponds to the zero electrons/s from the camera (see https://archive.stsci.edu/prepds/xdf).

photflamfloat

The calibration factor to convert pixel values in electrons/s to erg cm-2 s-1 Angstrom-1. Calculated using:

photflam = 10**(-abmag_zero/2.5 - 2*log10(photplam) - 0.9632)

which is a rearranged form of the following equation from http://www.stsci.edu/hst/acs/analysis/zeropoints:

abmag_zero = -2.5 log10(photflam)-5 log10(photplam)-2.408
photplamfloat

The pivot wavelength of the filter (Angstrom) (from http://www.stsci.edu/hst/acs/analysis/bandwidths)

photbwfloat

The effective bandwidth of the filter (Angstrom) (from http://www.stsci.edu/hst/acs/analysis/bandwidths)

imphot.rescale_hst_like_muse(hst, muse, inplace=True)[source]

Rescale an HST image to have the same flux units as a given MUSE image.

Parameters:
hstmpdaf.obj.Image

The HST image to be resampled.

musempdaf.obj.Image or mpdaf.obj.Cube

A MUSE image or cube with the target flux units.

inplacebool

(This defaults to True, because HST images tend to be large) If True, replace the contents of the input HST image object with the rescaled image. If False, return a new Image object that contains the rescaled image.

Returns:
outmpdaf.obj.Image

The rescaled HST image.

imphot.regrid_hst_like_muse(hst, muse, inplace=True)[source]

Resample an HST image onto the spatial coordinate grid of a given MUSE image or MUSE cube.

Parameters:
hstmpdaf.obj.Image

The HST image to be resampled.

musempdaf.obj.Image of mpdaf.obj.Cube

The MUSE image or cube to use as the template for the HST image.

inplacebool

(This defaults to True, because HST images tend to be large) If True, replace the contents of the input HST image object with the resampled image. If False, return a new Image object that contains the resampled image.

Returns:
outmpdaf.obj.Image

The resampled HST image.

imphot.image_grids_aligned(im1, im2, tolerance=0.01)[source]

Return True if two images sample the same array of sky positions to within a specified tolerance.

Parameters:
im1mpdaf.obj.Image

The first of the images to be compared.

im2mpdaf.obj.Image

The second of the images to be compared.

tolerancefloat

The tolerance within which positions must match, in arcseconds.

Returns:
outbool

True if the images have matching coordinate grids; False if they don’t.

imphot.apply_corrections(im, imfit, corrections='xy,scale,zero', resample=False, inplace=False)[source]

Correct an image for pointing-errors, calibration errors and/or zero offsets, using the fitted errors returned by a preceding call to imphot.fit_image_photometry() or imphot.fit_star_photometry().

Note that this function is called on the user’s behalf by imphot.fit_image_photometry() and imphot.fit_star_photometry() if those functions are given the argument, apply=True.

Parameters:
immpdaf.obj.Image

The image to be corrected.

imfitimphot.FittedPhotometry

Fitted image properties returned by imphot.fit_image_photometry() or imphot.fit_star_photometry().

correctionsstr

A comma-separated list of the corrections that are to be applied, indicated by the following words.

xy - The x-axis and y-axis pointing error. scale - The flux scale-factor. zero - The flux zero-offset.

The default string is “x,y,scale,zero”, which applies all of the corrections.

resamplebool

This parameter controls how pointing errors are corrected. They can be corrected either by changing the coordinates of the pixels (resample=False), or by resampling the pixels to shift the image within the original coordinate grid (resample=True).

inplacebool

By default a new image is returned, and the input image is not changed. However by setting inplace=True, the image is corrected within the input container.

Returns:
outmpdaf.obj.Image

The corrected image.

imphot.fit_image_photometry(hst, muse, regions=None, fix_scale=None, fix_bg=None, fix_dx=None, fix_dy=None, fix_fwhm=None, fix_beta=None, min_scale=0, hst_fwhm=0.085, hst_beta=1.6, margin=2.0, segment=False, display=False, nowait=False, hardcopy=None, title=None, star=None, save=False, fig=None, taper=9, apply=False, resample=False, extramask=None, init_dx=None, init_dy=None)[source]

Given a MUSE image and an HST image that has been regridded and aligned onto the same coordinate grid as the MUSE image, use the HST image as a calibrated reference to fit for the flux-scale, the FWHM and beta parameter of a Moffat PSF, the position offset, and the background offset of the MUSE image relative to the HST image.

More specifically, it uses the HST image to reproduce the MUSE image by convolving the HST image with a PSF that best reproduces the resolution of the MUSE image, scales and offsets the HST image with numbers that best reproduce the fluxes of the MUSE image, and shifts the HST image by an amount that best lines up features in the two images. The PSF parameters, calibration factors and pointing offsets that best reproduce the MUSE image, are the outputs of this algorithm.

Optionally, if apply=True is passed to this function, pointing corrections and calibration corrections are derived from the fitted parameters, and a corrected version of the MUSE image is written to a new FITS file. See also the resample option, which controls how pointing errors are corrected.

Parameters:
hstmpdaf.obj.Image or filename

An HST image of the same area as the MUSE image, and that has been gridded onto the same pixel grid as the MUSE image. This can be given as an MPDAF Image object, or by the filename of the FITS file.

musempdaf.obj.Image or filename

The MUSE image to be characterized. This can be given as an MPDAF Image object, or by the filename of the FITS file.

regionsNone, str, or iterable

This can be None, “none”, or “”, to indicate that no regions are needed, “star” to restrict the fit to pixels within any region defined by the star=(ra,dec,radius) argument, “notstar” to restrict the fit to pixels outside any region defined the by star=(ra,dec,radius) argument, the name of a filename of the ds9 region file, or an iterable that returns successive lines of a ds9 region file.

The regions can be used to exclude problematic areas of an image or sources that would degrade the global PSF fit, such as saturated stars, stars with significant proper motion, and variable sources. Alternatively it can be used to restrict the fit to one or more objects by masking everything except small regions around these objects.

Only ds9 circle, ellipse and box regions are supported. Other region types and most of the configuration parameters found in ds9 region files, are simply ignored. If the file is created within ds9, be careful to tell ds9 whether you want the region to be included or excluded. Also be careful to request either fk5 or physical pixel coordinates, because other coordinate systems are not supported.

fix_scalefloat or None

The calibration scale factor, (MUSE_flux / HST_flux) is fixed to the specified value while fitting, unless the value is None.

fix_bgfloat or None

The calibration zero-offset, (MUSE_flux - HST_flux) is fixed to the specified value while fitting, unless the value is None.

fix_dxfloat or None

The x-axis pointing offset, (MUSE_x - HST_x) is fixed to the specified value while fitting, unless the value is None.

fix_dyfloat or None

The y-axis pointing offset, (MUSE_y - HST_y) is fixed to the specified value while fitting, unless the value is None.

fix_fwhmfloat or None

The FWHM of the Moffat PSF is fixed to the specified value while fitting, unless the value is None.

fix_betafloat or None

The beta exponent of the Moffat PSF is fixed to the specified value while fitting, unless the value is None.

hst_fwhmfloat

The FWHM of a Moffat model of the effective PSF of the HST. The default value that is used if this parameter is not specified, came from Moffat fits to stars within HST UDF images. To obtain the closest estimate to the dithered instrumental PSF, these fits were made to images with the smallest available pixel size (30mas).

hst_betafloat

The beta parameter of a Moffat model of the effective PSF of the HST. The default value came from Moffat fits to stars within HST UDF images, as described above for the hst_fwhm parameter.

marginfloat

The width (arcsec) of a margin of zeros to add around the image before processing. A margin is needed because most of the processing is performed using discrete Fourier transforms, which are periodic in the width of the image. Without a margin, features at one edge of the image would spill over to the opposite edge of the image when a position shift was applied, or when features were widened by convolving them with a larger PSF. The margin width should be the maximum of the largest expected position error between the two input images, and the largest expected PSF width.

segmentbool

If True, ignore areas that don’t contain significant objects by ignoring pixels that are below the median value in a morphologically opened version of the HST image.

displaybool

If True (the default), display the plot.

hardcopystr or None

Unless this is an empty string or the word, “none”, then it should contain a graphics file suffix supported by matplotlib, such as “pdf”, “jpg”, “png” or “eps”. Plots of the image fit will be written to a filename that starts with the name of the MUSE input file, after removing any .fits suffix, followed by “_image_fit.<suffix>”.

nowaitbool

When this argument is False, wait for the user to dismiss the plot before returning. This allows the user to interact with the plot. When this argument is True, the plot will dissapear as soon as another plot is drawn and the user will not be able to interact with it.

titlestr or None

A specific plot title, None or “none” to request the default title, or “” to indicate that no title is wanted.

star(float, float, float)

This option can be used to restrict the fitting procedure to a single star within the MUSE image. Its arguments are the right ascension and declination of the star in decimal degrees, and the radius of the area over which to perform the fit, in arcseconds. Beware that this option overrides the regions argument.

savebool

If True, save the result images of each input image to FITS files.

figmatplotlib.figure.Figure

When display==True, a matplot figure is automatically generated for the plots if fig is None. Alternatively a specific matplotlib figure can be specified via this argument.

taperint

This argument controls how transitions between unmasked and masked regions are softened. Because the fitting algorithm replaces masked pixels with zeros, bright sources that are truncted by masked regions cause sharp changes in brightness that look like real features and bias the fitted position error. To reduce this effect, pixels close to the boundary of a masked region are smoothly tapered towards zero over a distance specified by the taper argument. The method used to smooth the transition requires that the taper argument be an odd number of pixels, so if a floating point or even-valued integer is specified, this is quietly rounded up to the next highest odd number. Alternatively, the softening algorithm can be disabled by passing 0 (or any value below 2).

applybool

If True, derive corrections from the fitted position errors and calibration errors, apply these to the MUSE image, and write the resulting image to a FITS file. The name of the output file is based on the name of the input file, by replacing its “.fits” extension with “_aligned.fits”. If the input muse image was not read from a file, then a file called “muse_aligned.fits” is written in the current directory.

resamplebool

When apply==True, this argument determines how position errors are corrected. If resample=False, then the coordinate reference pixel (CRPIX1, CRPIX2) is adjusted to change the coordinates of the pixels without changing any pixel values. If resample=True, then the pixel values are resampled to shift the image without changing the coordinates of the pixels.

extramaskstr or mpdaf.obj.Image or numpy.ndarray

An optional masking image to apply to the MUSE image (and the shifted HST image). This can be a numpy bool array of the same shape as the MUSE image, with False elements denoting unmasked pixels and True elements denoting masked pixels. Alternatively it can be an MPDAF Image object of the same shape as the MUSE image, which contains an image of integer valued pixels, where 0 denotes unmasked pixels, and 1 denotes masked pixels. Finally, it can be a FITS file that contains an IMAGE extension called ‘DATA’, which must contain an image that has the same shape as the MUSE image, and should contain integers that are 0 for unmasked pixels, and 1 for masked pixels. Beware that if the mask is specified as an MPDAF Image or a FITS image, then the WCS information must match the MUSE image.

init_dxfloat or None

An initial guess for the x-axis pointing offset, (MUSE_x - HST_x), or None to request the default (currently 0.0). This argument is ignored if the fix_dx argument has been given a value.

init_dyfloat or None

An initial guess for the y-axis pointing offset, (MUSE_y - HST_y), or None to request the default (currently 0.0). This argument is ignored if the fix_dy argument has been given a value.

Returns
——-
outFittedImagePhotometry

An object that contains the fitted parameters and a textual report about the fit.

class imphot.FittedImagePhotometry(muse, results, rms_error)[source]

The class of the object that fit_image_photometry() returns.

Note that both this class and FittedStarPhotometry are derived from the FittedPhotometry class, which contains the essential features of the fitted parameters, without features that are specific to the fitting method.

Parameters:
musempdaf.obj.Image

The MUSE image that the fit was performed on.

resultslmfit.ModelResult

The model fitting results.

rms_errorfloat

The root-mean square of the residual image pixels, in the same units as the pixels of the original MUSE image.

Attributes:
namestr

The basename of the MUSE FITS without the .fits extension.

fit_reportstr

A printable report on the fit from the least-squares fitting function.

scaleFittedValue

The best-fit value and error of the calibration scale factor, (MUSE.flux / HST.flux).

bgFittedValue

The best-fit value and error of the calibration offset, (MUSE.flux - HST.flux).

dxFittedValue

The best-fit value and error of the x-axis pointing offset, MUSE.x-HST.x (arcsec). This is the distance that features in the HST image had to be moved to the right, in the direction of increasing x-axis pixel index, to line them up with the same features in the MUSE image. One way to correct the pointing error of the MUSE observation, is to divide ‘dx’ by the pixel size along the x-axis (usually 0.2 arcsec), then add the resulting pixel offset to the CRPIX1 header parameter.

dyFittedValue

The best-fit value and error of the y-axis pointing offset, MUSE.y-HST.y (arcsec). This is the distance that features in the HST image had to be moved upwards, in the direction of increasing y-axis pixel index, to line them up with the same features in the MUSE image. One way to correct the pointing error of the MUSE observation, is to divide ‘dy’ by the pixel size along the y-axis (usually 0.2 arcsec), then add the resulting pixel offset to the CRPIX2 header parameter.

draFittedValue

The right-ascension error (arcsec) that corresponds to the pointing error dx,dy. This is the angular distance that features in the HST image had to be moved towards increased right-ascension, to line them up with the same feaures in the MUSE image. One way to correct the pointing error of the MUSE observation is to subtract ‘dra’ from the CRVAL1 header value of the MUSE observation.

ddecFittedValue

The declination error (arcsec) that corresponds to the pointing error dx,dy. This is the angular distance that features in the HST image had to be moved towards increased declination, to line them up with the same feaures in the MUSE image. One way to correct the pointing error of the MUSE observation is to subtract ‘ddec’ from the CRVAL2 header value of the MUSE observation.

fwhmFittedValue

The best-fit value and error of the FWHM of the Moffat PSF.

betaFittedValue

The best-fit value and error of the beta parameter of the Moffat PSF

rchifloat

The reduced chi-squared value of the fit. Beware that the absolute value of this number is not very useful, because the fit is performed to the complex pixels of zero-padded FFTs, which are not independent observables, rather than to image pixels.

rms_errorfloat

The root-mean square of the residual image pixels, in the same units as the pixels of the original MUSE image.

summary(header=True)

Return a string that summarizes the results, optionally listed under a 3-line heading of column descriptions.

Parameters:
headingbool

True to include 3 header lines above the values (default=True)

class imphot.FitImagePhotometryMP(hst_filename, muse_filenames, kwargs={}, nworker=0)[source]

A multiprocessing iterator that creates a pool or worker processes to repeatedly call fit_image_photometry() for each of a list of MUSE files, returning the results via the iterator as they become available.

Parameters:
hst_filenamestr

The name of a FITS fie that contains an HST has the same coordinate grid as the MUSE image FITS files that are to be processed.

muse_filenameslist of str

A list of filenames of the FITS files of the MUSE images that are to be processed. These must all be of the same field, with the same image coordinate grid as the HST file.

kwargsdict

An optional dictionary of keyword/value arguments to be passed to fit_image_photometry(). The default is the empty dictionary, {}.

nworkerint

The number of worker processes to use. The default is 0, which creates multiprocessing.cpu_count() processes. Alternatively, if a negative number, -n, is specified, then max(multiprocessing.cpu_count()-n,1) processes are created.

next()[source]

Return the results from the next image in the list of input MUSE images.

Returns:
outFittedImagePhotometry

The fitting results from the next image.

imphot.fit_star_photometry(hst, muse, star, fix_fwhm=None, fix_beta=None, display=False, nowait=False, hardcopy=None, title=None, fig=None, apply=False, resample=False)[source]

Given a MUSE image and an HST image that have been regridded and aligned onto the same coordinate grid as the MUSE image, fit moffat profiles to a specified star within both of these images, then compare the fits to determine the X,Y pointing offset, pixel zero-offset and pixel scaling error of the MUSE image. Also use the fitted FWHM and beta parameters to indicate the characteristics of the MUSE PSF.

Optionally, if apply=True is passed to this function, pointing corrections and calibration corrections are derived from the fitted parameters, and a corrected version of the MUSE image is written to a new FITS file. See also the resample option, which controls how pointing errors are corrected.

Parameters:
hstmpdaf.obj.Image

An HST image of the same area as the MUSE image, and that has been gridded onto the same pixel grid as the MUSE image.

musempdaf.obj.Image

The MUSE image to be characterized.

starfloat, float, float

The central right-ascension, declination and radius of the area within which to fit the stellar brightness profile. The right ascension and declination are in degrees. The radius is in arcseconds.

fix_fwhmfloat or None

The FWHM of the Moffat PSF in the MUSE image is fixed to the specified value while fitting, unless the value is None. This does not affect the PSF that is fitted to the star in the HST image.

fix_betafloat or None

The beta exponent of the Moffat PSF in the MUSE image is fixed to the specified value while fitting, unless the value is None. This does not affect the PSF that is fitted to the star in the HST image.

displaybool

If True (the default), display the image.

hardcopystr or None

Unless this is an empty string or the word, “none”, then it should contain a graphics file suffix supported by matplotlib, such as “pdf”, “jpg”, “png” or “eps”. Plots of the star fits will be written to a filename that starts with the name of the MUSE input file, after removing any .fits suffix, followed by “_star_fit.<suffix>”.

nowaitbool

When this argument is False, wait for the user to dismiss the plot before returning. This allows the user to interact with the plot. When this argument is True, the plot will dissapear as soon as another plot is drawn and the user will not be able to interact with it.

titlestr or None

A specific plot title, None or “none” to request the default title, or “” if no title is wanted.

applybool

If True, derive corrections from the fitted position errors and calibration errors, apply these to the MUSE image, and write the resulting image to a FITS file. The name of the output file is based on the name of the input file, by replacing its “.fits” extension with “_aligned.fits”. If the input muse image was not read from a file, then a file called “muse_aligned.fits” is written in the current directory.

resamplebool

When apply==True, this argument determines how position errors are corrected. If resample=False, then the coordinate reference pixel (CRPIX1, CRPIX2) is adjusted to change the coordinates of the pixels without changing any pixel values. If resample=True, then the pixel values are resampled to shift the image without changing the coordinates of the pixels.

Returns:
outFittedStarPhotometry

The results are returned in an object.

class imphot.FittedStarPhotometry(muse, muse_results, hst_results, rms_error, ra, dec)[source]

A class for returning the results of fit_star_photometry().

Note that both this class and FittedImagePhotometry are derived from the FittedPhotometry class, which contains the essential features of the fitted parameters, without features that are specific to the fitting method.

Parameters:
musempdaf.obj.Image

The MUSE image that the fit was performed on.

muse_resultslmfit.ModelResult

The results of fitting a 2D Moffat to a star in the MUSE image.

hst_resultslmfit.ModelResult

The results of fitting a 2D Moffat to a star in the HST image.

rafloat

The right-ascension used for the origin of the fit (degrees)

decfloat

The declination used for the origin of the fit (degrees)

rms_errorfloat

The root-mean square of the pixel residuals of the fit to the star in the MUSE image, in the same units as the pixels of the original MUSE image.

Attributes:
namestr

The basename of the MUSE FITS without the .fits extension.

fit_reportstr

A multi-line string listing verbose reports on both the MUSE and the HST star fits.

scaleFittedValue

The best-fit value and error of the calibration scale factor, MUSE/HST. This is derived from muse_flux/hst_flux.

bgFittedValue

The best-fit value and error of the calibration offset, MUSE-HST. This is derived from muse_bg - hst_bg.

dxFittedValue

The best-fit value and error of the x-axis pointing offset, muse_dx-hst_dx (arcsec). This is the distance that features in the HST image had to be moved to the right, in the direction of increasing x-axis pixel index, to line them up with the same features in the MUSE image. One way to correct the pointing error of the MUSE observation, is to divide ‘dx’ by the pixel size along the x-axis (usually 0.2 arcsec), then add the resulting pixel offset to the CRPIX1 header parameter.

dyFittedValue

The best-fit value and error of the y-axis pointing offset, muse_y-hst_y (arcsec). This is the distance that features in the HST image had to be moved upwards, in the direction of increasing y-axis pixel index, to line them up with the same features in the MUSE image. One way to correct the pointing error of the MUSE observation, is to divide ‘dy’ by the pixel size along the y-axis (usually 0.2 arcsec), then add the resulting pixel offset to the CRPIX2 header parameter.

draFittedValue

The right-ascension error (arcsec) that corresponds to the pointing error dx,dy. This is the angular distance that features in the HST image had to be moved towards increased right-ascension, to line them up with the same feaures in the MUSE image. One way to correct the pointing error of the MUSE observation is to subtract ‘dra’ from the CRVAL1 header value of the MUSE observation.

ddecFittedValue

The declination error (arcsec) that corresponds to the pointing error dx,dy. This is the angular distance that features in the HST image had to be moved towards increased declination, to line them up with the same feaures in the MUSE image. One way to correct the pointing error of the MUSE observation is to subtract ‘ddec’ from the CRVAL2 header value of the MUSE observation.

fwhmFittedValue

The best-fit value and error of the FWHM of the Moffat PSF (arcsec). This is the same as muse_fwhm.

betaFittedValue

The best-fit value and error of the beta parameter of the Moffat PSF. This is the same as muse_beta.

rms_errorfloat

The root-mean square of the pixel residuals of the fit to the star in the MUSE image, in the same units as the pixels of the original MUSE image.

rafloat

The right-ascension used for the origin of the fit (degrees)

decfloat

The declination used for the origin of the fit (degrees)

muse_fwhmFittedValue

The best-fit value and error of the FWHM of the Moffat PSF fitted to the star in the MUSE image (arcsec)

muse_betaFittedValue

The best-fit value and error of the beta parameter of the Moffat PSF fitted to the star in the MUSE image.

muse_dxFittedValue

The x-axis image offset of the center of the star fitted in the MUSE image to the coordinate self.ra,self.dec (arcsec).

muse_dyFittedValue

The y-axis image offset of the center of the star fitted in the MUSE image to the coordinate self.ra,self.dec (arcsec).

muse_bgFittedValue

The fitted flux zero-offset under the fitted Moffat PSF in the MUSE image of the star.

muse_fluxFittedValue

The fitted total flux of the Moffat PSF that was fitted to the star in the MUSE image, using the flux units of the MUSE image.

muse.rchifloat

The reduced chi-squared of the fit of the Moffat PSF to the star in the MUSE image.

hst_fwhmFittedValue

The best-fit value and error of the FWHM of the Moffat PSF fitted to the star in the HST image (arcsec)

hst_betaFittedValue

The best-fit value and error of the beta parameter of the Moffat PSF fitted to the star in the HST image.

hst_dxFittedValue

The x-axis image offset of the center of the star fitted in the HST image to the coordinate self.ra,self.dec (arcsec).

hst_dyFittedValue

The y-axis image offset of the center of the star fitted in the HST image to the coordinate self.ra,self.dec (arcsec).

hst_bgFittedValue

The fitted flux zero-offset under the fitted Moffat PSF in the HST image of the star.

hst_fluxFittedValue

The fitted total flux of the Moffat PSF fitted to the star in the HST image, using the flux units of the MUSE image.

hst_rchifloat

The reduced chi-squared of the fit of the Moffat PSF to the star in the HST image.

summary(header=True)

Return a string that summarizes the results, optionally listed under a 3-line heading of column descriptions.

Parameters:
headingbool

True to include 3 header lines above the values (default=True)

class imphot.FitStarPhotometryMP(hst_filename, muse_filenames, kwargs, nworker=0)[source]

A multiprocessing iterator that creates a pool or worker processes to repeatedly call fit_star_photometry() for each of a list of MUSE files, returning the results via the iterator as they become available.

Parameters:
hst_filenamestr

The name of a FITS fie that contains an HST has the same coordinate grid as the MUSE image FITS files that are to be processed.

muse_filenameslist of str

A list of filenames of the FITS files of the MUSE images that are to be processed. These must all be of the same field, with the same image coordinate grid as the HST file.

nworkerint

The number of worker processes to use. The default is 0, which creates multiprocessing.cpu_count() processes. Alternatively, if a negative number, -n, is specified, then max(multiprocessing.cpu_count()-n,1) processes are created.

kwargsdict

A dictionary of keyword/value arguments to be passed to fit_star_photometry(). This should always include a value for the star=(ra,dec,radius) argument.

next()[source]

Return the results from the next image in the list of input MUSE images.

Returns:
outFittedStarPhotometry

The fitting results from the next image.

imphot.fit_image_and_star_photometry(hst, muse, star, regions='star', fix_scale=None, fix_bg=None, fix_dx=None, fix_dy=None, fix_fwhm=None, fix_beta=None, hst_fwhm=0.085, hst_beta=1.6, margin=2.0, segment=False, display=False, nowait=False, hardcopy=None, title=None, save=False, fig=None)[source]

Given a MUSE image and an HST image that has been regridded and aligned onto the same coordinate grid as the MUSE image, use the HST image as a calibrated reference to fit for the flux-scale, the FWHM and beta parameter of a Moffat PSF, the position offset, and the background offset of the MUSE image relative to the HST image. The fit is performed two ways; first by calling fit_image_photometry(), then again by calling fit_star_photometry(). Both sets of results are returned.

Parameters:
hstmpdaf.obj.Image or filename

An HST image of the same area as the MUSE image, and that has been gridded onto the same pixel grid as the MUSE image. This can be given as an MPDAF Image object, or by the filename of the FITS file.

musempdaf.obj.Image or filename

The MUSE image to be characterized. This can be given as an MPDAF Image object, or by the filename of the FITS file.

star(float, float, float)

This argument indicates the position of a single star within the MUSE image and the radius of a circular area of pixels around it that should be used in the fitting process by fit_star_photometry().

If the optional regions argument is set to the special string, “star”, then the star’s position and the associated radius are also passed to fit_image_photometry() to indicate that the image fit should also be restricted to the same circular area around the star.

The value of the star argument is a tuple of the right ascension and declination of the star in decimal degrees, and the radius of the area over which to perform the fit, in arcseconds.

regionsstr or iterable or None

This can be None, “none”, or “”, to indicate that no regions are needed in the global image fit, “star” to restrict the image fit to pixels within any region defined by the star=(ra,dec,radius) argument, “notstar” to restrict the image fit to pixels outside any region defined the by star=(ra,dec,radius) argument, the name of a filename of the ds9 region file, or an iterable that returns successive lines of a ds9 region file.

These regions are passed to fit_image_photometry(). They can be used to exclude problematic areas of an image, or to exclude sources that would degrade the global PSF fit, such as saturated stars, stars with significant proper motion, and variable sources. Alternatively it can be used to restrict the fit to one or more objects by masking everything except small regions around these objects. In particular, if the value of the regions argument is “star”, then the fit is restricted to the circular area specified by star argument.

Only ds9 circle, ellipse and box regions are supported. Other region types and most of the configuration parameters found in ds9 region files, are simply ignored. For each region, be careful to tell ds9 whether you want the region to be included or excluded. Also be careful to request either fk5 or physical pixel coordinates, because other coordinate systems are not supported.

fix_scalefloat or None

The calibration scale factor, (MUSE_flux / HST_flux) is fixed to the specified value while fitting, unless the value is None. This parameter affects the image fit but not the star-profile fit.

fix_bgfloat or None

The calibration zero-offset, (MUSE_flux - HST_flux) is fixed to the specified value while fitting, unless the value is None. This parameter affects the image fit but not the star-profile fit.

fix_dxfloat or None

The x-axis pointing offset, (MUSE_x - HST_x) is fixed to the specified value while fitting, unless the value is None. This parameter affects the image fit but not the star-profile fit.

fix_dyfloat or None

The y-axis pointing offset, (MUSE_y - HST_y) is fixed to the specified value while fitting, unless the value is None. This parameter affects the image fit but not the star-profile fit.

fix_fwhmfloat or None

The FWHM of the Moffat PSF is fixed to the specified value while fitting, unless the value is None. This parameter affects both the image fit and the FWHM of the star profile that is fit to the MUSE image.

fix_betafloat or None

The beta exponent of the Moffat PSF is fixed to the specified value while fitting, unless the value is None. This parameter affects both the image fit and the beta parameter of the star profile that is fit to the MUSE image.

hst_fwhmfloat

The FWHM of a Moffat model of the effective PSF of the HST. The default value that is used if this parameter is not specified, came from Moffat fits to stars within HST UDF images. To obtain the closest estimate to the dithered instrumental PSF, these fits were made to images with the smallest available pixel size (30mas).

hst_betafloat

The beta parameter of a Moffat model of the effective PSF of the HST. The default value that is used if this parameter is not specified, came from Moffat fits to stars within HST UDF images, as described above for the hst_fwhm parameter.

marginfloat

The width (arcsec) of a margin of zeros to add around the image before processing. A margin is needed because most of the processing is performed using discrete Fourier transforms, which are periodic in the width of the image. Without a margin, features at one edge of the image would spill over to the opposite edge of the image when a position shift was applied, or when features were widened by convolving them with a larger PSF. The margin width should be the maximum of the largest expected position error between the two input images, and the largest expected PSF width.

segmentbool

If True, ignore areas that don’t contain significant objects by ignoring pixels that are below the median value in a morphologically opened version of the HST image.

displaybool

If True (the default), display the plot.

hardcopystr or None

If this is a non-empty string, then it should contain a graphics file suffix supported by matplotlib, such as “pdf”, “jpg”, “png” or “eps”. Plots of the fits will be written to filenames that start with the name of the MUSE input file, after removing any .fits suffix, followed by “_image_fit.<suffix>” for the plot of the image fit, and “_star_fit.<suffix>” for the plots of the star fits.

nowaitbool

When this argument is False, wait for the user to dismiss the plot before returning. This allows the user to interact with the plot. When this argument is True, the plot will dissapear as soon as another plot is drawn and the user will not be able to interact with it.

titlestr or None

A specific plot title, or None to request the default title. Specify “” if no title is wanted.

savebool

If True, save the result images of each input image to FITS files.

Returns:
out(FittedImagePhotometry, FittedStarPhotometry)

An object that contains the fitted parameters from both fit_image_photometry() and fit_star_photometry().

class imphot.FitImageAndStarPhotometryMP(hst_filename, muse_filenames, kwargs, nworker=0)[source]

A multiprocessing iterator that creates a pool or worker processes to repeatedly call fit_image_and_star_photometry() for each of a list of MUSE files, returning successive results from fit_image_and_star_photometry() via the iterator, as they become available.

Parameters:
hst_filenamestr

The name of a FITS fie that contains an HST has the same coordinate grid as the MUSE image FITS files that are to be processed.

muse_filenameslist of str

A list of filenames of the FITS files of the MUSE images that are to be processed. These must all be of the same field, with the same image coordinate grid as the HST file.

kwargsdict

A dictionary of keyword/value arguments to be passed to fit_image_and_star_photometry().

nworkerint

The number of worker processes to use. The default is 0, which creates multiprocessing.cpu_count() processes. Alternatively, if a negative number, -n, is specified, then max(multiprocessing.cpu_count()-n,1) processes are created.

next()[source]

Return the results from the next image in the list of input MUSE images.

Returns:
outFittedStarPhotometry

The fitting results from the next image.

imphot.bandpass_image(cube, wavelengths, sensitivities, unit_wave=Unit('Angstrom'), interpolation='linear', truncation_warning=True, nprocess=0)[source]

Given a cube of images versus wavelength and the bandpass filter-curve of a wide-band monochromatic instrument, extract an image from the cube that has the spectral response of the monochromatic instrument.

For example, this can be used to create a MUSE image that has the same spectral characteristics as an HST image. The MUSE image can then be compared to the HST image without having to worry about any differences caused by different spectral sensitivities.

For each spectral plane of the cube, the filter-curve is integrated over the width of that spectral plane to obtain a weight, w[n]. The output image is then given by the following weighted mean:

output_image = sum(w[n] * cube_image[n]) / sum(w[n])

In practice, to accomodate masked pixels, the w[n] array is expanded into a cube w[n,y,x], and the weights of individual masked pixels in the cube are zeroed before the above equation is applied.

If the wavelength axis of the cube only partly overlaps the bandpass of the filter-curve, the filter curve is truncated to fit within the bounds of the wavelength axis. A warning is printed to stderr if this occurs, because this results in an image that lacks flux from some of the wavelengths of the requested bandpass.

This function is equivalent to Cube.bandpass_image() in MPDAF, but it is much faster, because of its use of multiple processes.

Parameters:
cubempdaf.obj.Cube

The cube to be processed.

wavelengthsnumpy.ndarray

An array of the wavelengths of the filter curve, listed in ascending order of wavelength. Outside the listed wavelengths the filter-curve is assumed to be zero.

sensitivitiesnumpy.ndarray

The relative flux sensitivities at the wavelengths in the wavelengths array. These sensititivies will be normalized, so only their relative values are important.

unit_waveastropy.units.Unit

The units used in the array of wavelengths. The default is angstroms. To specify pixel units, pass None.

interpolationstr

The form of interpolation to use to integrate over the filter curve. This should be one of:

"linear"     : Linear interpolation
"cubic"      : Cubic spline interpolation (very slow)

The default is linear interpolation. If the filter curve is well sampled and its sampling interval is narrower than the wavelength pixels of the cube, then this should be sufficient. Alternatively, if the sampling interval is significantly wider than the wavelength pixels of the cube, then cubic interpolation should be used instead. Beware that cubic interpolation is much slower than linear interpolation.

truncation_warningbool

By default, a warning is reported if the filter bandpass exceeds the wavelength range of the cube. The warning can be disabled by setting this argument to False.

nprocessint

The number of worker processes to use. The default is 0, which creates multiprocessing.cpu_count() processes. Alternatively, if a negative number, -n, is specified, then max(multiprocessing.cpu_count()-n,1) processes are created.

Returns:
outImage

An image formed from the filter-weighted mean of spectral planes in the cube that overlap the bandpass of the filter curve.

imphot.fit_photometry_main()[source]

This is the main function of the fit_photometry script. It attempts to determine the photometry characteristics of a MUSE image via comparisons with an HST image, according to arguments passed by the user.

imphot.regrid_hst_to_muse_main()[source]

This is the main function of the regrid_hst_to_muse script.

Given a MUSE image or MUSE cube and one or more HST images with 30mas pixels, it resamples the HST images onto the pixel grid of the MUSE image, and scales their flux units from electrons/s to the flux units of the MUSE image (usually 1e-20 erg/cm^2/s/Angstrom).

imphot.make_wideband_image_main()[source]

This is the main function of the make_wideband_image script.

Given a MUSE cube and the wavelength sensitivity curve of a monochromatic camera, extract an image from the cube that has the same wavelength sensitivity curve as the camera.

imphot.ds9regions source-code documentation

class imphot.ds9regions.Ds9Regions(regions)[source]

Parse circle, ellipse and box regions from a ds9 region file.

Only circle, ellipse, and box regions are recognized. Other types of ds9 region and configuration parameters are simply ignored. The center coordinates of each region are assumed to be in equatorial coordinates and sizes are assumed to be angular sizes, not pixel sizes.

The Ds9Regions object can be treated like an indexable tuple. For example:

regions = ds9regions.Ds9Regions("myfile.reg")

# Indexing and len() are supported.

print("First region = ", regions[0])
print("Last region = ", regions[-1])

# Repeatable iteration like a tuple is supported:

for region in regions:
   print(region)
Parameters:
regionsfilename or iterable

Either a string that contains the name of a ds9 region file, or an iterable that returns successive lines of a ds9 file.

class imphot.ds9regions.Ds9Region(shape, exclude, system, x, y)[source]

A generic Ds9Region.

Parameters:
shapestr

The name of the shape (“circle”, “ellipse”, or “box”).

excludebool

True if the region is to be excluded.

systemstr

The name of the coordinate system.

xfloat

The x coordinate of the center of the shape (degrees).

yfloat

The y coordinate of the center of the shape (degrees).

Attributes:
shapestr

The name of the shape (“circle”, “ellipse”, or “box”).

excludebool

True if the region is to be excluded.

systemstr

The name of the coordinate system.

xfloat

The x coordinate of the center of the shape (degrees).

yfloat

The y coordinate of the center of the shape (degrees).

class imphot.ds9regions.Ds9Circle(exclude, system, x, y, radius)[source]

A circular ds9 region.

Parameters:
excludebool

True if the region is to be excluded.

systemstr

The name of the coordinate system.

xfloat

The x-axis coordinate of the center of the shape (degrees).

yfloat

The declination of the center of the shape (degrees).

radiusfloat

The radius of the circle (degrees).

Attributes:
shapestr

The name of the shape (“circle”, “ellipse”, or “box”).

excludebool

True if the region is to be excluded.

systemstr

The name of the coordinate system.

xfloat

The x coordinate of the center of the shape (degrees).

yfloat

The y coordinate of the center of the shape (degrees).

radiusfloat

The radius of the circle (degrees).

class imphot.ds9regions.Ds9Ellipse(exclude, system, x, y, width, height, pa)[source]

An elliptical ds9 region.

Parameters:
excludebool

True if the region is to be excluded.

systemstr

The name of the coordinate system.

xfloat

The x-axis coordinate of the center of the shape (degrees).

yfloat

The declination of the center of the shape (degrees).

widthfloat

The width (degrees) of the ellipse when the position angle is 0.

heightfloat

The height (degrees) of the ellipse when the position angle is 0.

pafloat

The position angle of the rectangle, east of celestial north for sky coordinate systems.

Attributes:
shapestr

The name of the shape (“circle”, “ellipse”, or “box”).

excludebool

True if the region is to be excluded.

systemstr

The name of the coordinate system.

xfloat

The x coordinate of the center of the shape (degrees).

yfloat

The y coordinate of the center of the shape (degrees).

widthfloat

The width (degrees) of the ellipse when the position angle is 0.

heightfloat

The height (degrees) of the ellipse when the position angle is 0.

pafloat

The position angle of the rectangle, east of celestial north for sky coordinate systems.

class imphot.ds9regions.Ds9Box(exclude, system, x, y, width, height, pa)[source]

A rectangular ds9 region.

Parameters:
excludebool

True if the region is to be excluded.

systemstr

The name of the coordinate system.

xfloat

The x-axis coordinate of the center of the shape (degrees).

yfloat

The declination of the center of the shape (degrees).

widthfloat

The width (degrees) of the rectangle when the position angle is 0.

heightfloat

The height (degrees) of the rectangle when the position angle is 0.

pafloat

The position angle of the rectangle, east of celestial north for sky coordinate systems.

Attributes:
shapestr

The name of the shape (“circle”, “ellipse”, or “box”).

excludebool

True if the region is to be excluded.

systemstr

The name of the coordinate system.

xfloat

The x coordinate of the center of the shape (degrees).

yfloat

The y coordinate of the center of the shape (degrees).

widthfloat

The width (degrees) of the rectangle when the position angle is 0.

heightfloat

The height (degrees) of the rectangle when the position angle is 0.

pafloat

The position angle of the rectangle, east of celestial north for sky coordinate systems.