source_from_image

simcado.source.source_from_image(images, lam, spectra, plate_scale, oversample=1, units='ph/s/m2', flux_threshold=0, center_offset=(0, 0), conserve_flux=True, **kwargs)[source]

Create a Source object from an image or a list of images.

Note

plate_scale is the original plate scale of the images. If this is not the same as the plate scale of the Detector then you will need to specify oversample to interpolate between the two scales. I.e. oversample = Image plate scale / Detector plate scale

Parameters
imagesnp.ndarray, list

A single or list of np.ndarrays describing where the flux is coming from The spectrum for each pixel in the image is weighted by the pixel value.

lamnp.ndarray

An array contains the centres of the wavelength bins for the spectra

spectranp.ndarray

A (n,m) array with n spectra, each with m bins

plate_scalefloat

[arcsec] The plate scale of the images in arcseconds (e.g. 0.004”/pixel)

oversampleint

The factor with which to oversample the image. Each image pixel is split into (oversample)^2 individual point sources.

unitsstr, optional

The energy units of the spectra. Default is [ph/s/m2]

flux_thresholdfloat, optional

If there is noise in the image, set threshold to the noise limit so that only real photon sources are extracted. Default is 0.

center_offset(float, float)

[arcsec] If the centre of the image is offset, add this offset to (x,y) coordinates.

conserve_fluxbool, optional

If True, when the image is rescaled, flux is conserved i.e. np.sum(image) remains constant If False, the maximum value of the image stays constant after rescaling i.e. np.max(image) remains constant

Returns
srcSource object

Notes

Currently only one object per image is supported.

Examples

To create a Source object we need an image that describes the spatial distribution of the object of interest and spectrum. For the sake of ease we will assign a generic elliptical galaxy spectrum to the image.:

>>> from astropy.io import fits
>>> from simcado.source import SED, source_from_image
>>>
>>> im = fits.getdata("galaxy.fits")
>>> lam, spec = SED("elliptical")
>>> src = source_from_image(im, lam, spec,
                            plate_scale=0.004)

Note Here we have assumed that the plate scale of the image is the same as the MICADO wide-field mode, i.e. 0.004 arcseconds. If the image is from a real observation, or it was generated with a different pixel scale, we will need to tell SimCADO about this:

>>> src = source_from_image(im, lam, spec,
                            plate_scale=0.01,
                            oversample=2.5)

If the image is from real observations, chances are good that the background flux is higher than zero. We can set a threshold in order to tell SimCADO to ignore all pixel with values below the background level:

>>> src = source_from_image(im, lam, spec,
                            plate_scale=0.01,
                            oversample=2.5,
                            flux_threshold=0.2)

Finally, if the image centre is not the centre of the observation, we can shift the image relative to the MICADO field of view. The units for the offset are [arcsec]:

>>> src = source_from_image(im, lam, spec,
                            plate_scale=0.01,
                            oversample=2.5,
                            flux_threshold=0.2,
                            center_offset=(10,-15))