Source¶
-
class
simcado.source.
Source
(filename=None, lam=None, spectra=None, x=None, y=None, ref=None, weight=None, **kwargs)[source]¶ Bases:
object
Create a source object from a file or from arrays
Source class generates the arrays needed for source. It takes various inputs and converts them to an array of positions and references to spectra. It also converts spectra to photons/s/voxel. The default units for input data is ph/s/m2/bin.
The internal variables are related like so:
f(x[i], y[i]) = spectra[ref[i]] * weight[i]
- Parameters
- filenamestr
FITS file that contains either a previously saved
Source
object or a data cube with dimensions x, y, lambda. ASource
object is identified by the header keyword SIMCADO with value SOURCE.- or
- lamnp.array
[um] Wavelength bins of length (m)
- spectranp.array
[ph/s/m2/bin] A (n, m) array with n spectra, each with m spectral bins
- x, ynp.array
[arcsec] coordinates of where the emitting sources are relative to the centre of the field of view
- refnp.array
the index for .spectra which connects a position (x, y) to a spectrum f(x[i], y[i]) = spectra[ref[i]] * weight[i]
- weightnp.array
A weighting to scale the relevant spectrum for each position
Attributes Summary
Return keys of the info dict
Methods Summary
Add an EmissionCurve for the background surface brightness of the object
apply_optical_train
(self, opt_train, detector)Apply all effects along the optical path to the source photons
dump
(self, filename)Save to filename as a pickle
image_in_range
(self, psf, lam_min, lam_max, …)Find the sources that fall in the chip area and generate an image for the wavelength range [lam_min, lam_max)
load
(filename)Load :class:’.Source’ object from filename
on_grid
(self[, pix_res])Return an image with the positions of all sources.
photons_in_range
(self[, lam_min, lam_max])Number of photons between lam_min and lam_max in units of [ph/s/m2]
project_onto_chip
(self, image, chip)Re-project the photons onto the same grid as the detectors use
read
(self, filename)Read in a previously saved
Source
FITS filerotate
(self, angle[, unit, use_orig_xy])Rotates the
x
andy
coordinates byangle
[degrees]scale_spectrum
(self[, idx, mag, filter_name])Scale a certain spectrum to a certain magnitude
scale_with_distance
(self, distance_factor)Scale the source for a new distance
shift
(self[, dx, dy, use_orig_xy])Shifts the coordinates of the source by (dx, dy) in [arcsec]
write
(self, filename)Write the current Source object out to a FITS file
Attributes Documentation
-
info_keys
¶ Return keys of the info dict
Methods Documentation
-
add_background_surface_brightness
(self)[source]¶ Add an EmissionCurve for the background surface brightness of the object
-
apply_optical_train
(self, opt_train, detector, chips='all', sub_pixel=False, **kwargs)[source]¶ Apply all effects along the optical path to the source photons
- Parameters
- opt_trainsimcado.OpticalTrain
the object containing all information on what is to happen to the photons as they travel from the source to the detector
- detectorsimcado.Detector
the object representing the detector
- chipsint, str, list, optional
The IDs of the chips to be readout. “all” is also acceptable
- sub_pixelbool, optional
if sub-pixel accuracy is needed, each source is shifted individually. Default is False
- Other Parameters
- INST_DEROT_PERFORMANCEfloat
[0 .. 100] Percentage of the sky rotation that the derotator removes
- SCOPE_JITTER_FWHMfloat
[arcsec] The FWMH of the gaussian blur caused by jitter
- SCOPE_DRIFT_DISTANCEfloat
[arcsec] How far from the centre of the field of view has the telescope drifted during a DIT
Notes
Output array is in units of [ph/s/pixel] where the pixel is internal oversampled pixels - not the pixel size of the detector chips
-
image_in_range
(self, psf, lam_min, lam_max, chip, **kwargs)[source]¶ Find the sources that fall in the chip area and generate an image for the wavelength range [lam_min, lam_max)
Output is in [ph/s/pixel]
- Parameters
- psfpsf.PSF object
The PSF that the sources will be convolved with
- lam_min, lam_maxfloat
[um] the wavelength range relevant for the psf
- chipstr, detector.Chip
detector.Chip : the chip that will be seeing this image.
str : [“tiny”, “small”, “center”] -> [128, 1024, 4096] pixel chips
- Returns
- slice_arraynp.ndarray
the image of the source convolved with the PSF for the given range
-
on_grid
(self, pix_res=0.004)[source]¶ Return an image with the positions of all sources.
The pixel values correspond to the number of emitting objects in that pixel
- Parameters
- pix_resfloat
[arcsec] The grid spacing
- Returns
- im2D array
A numpy array containing an image of where the sources are
-
photons_in_range
(self, lam_min=None, lam_max=None)[source]¶ Number of photons between lam_min and lam_max in units of [ph/s/m2]
Calculate how many photons for each source exist in the wavelength range defined by lam_min and lam_max.
- Parameters
- lam_min, lam_maxfloat, optional
[um] integrate photons between these two limits. If both are
None
, limits are set at lam[0], lam[-1] for the source’s wavelength range
- Returns
- slice_photonsfloat
[ph/s/m2] The number of photons in the wavelength range
-
project_onto_chip
(self, image, chip)[source]¶ Re-project the photons onto the same grid as the detectors use
- Parameters
- imagenp.ndarray
the image to be re-projected
- chipdetector.Chip
the chip object where the image will land
-
read
(self, filename)[source]¶ Read in a previously saved
Source
FITS file- Parameters
- filenamestr
Path to the file
-
rotate
(self, angle, unit='degree', use_orig_xy=False)[source]¶ Rotates the
x
andy
coordinates byangle
[degrees]- Parameters
- anglefloat
Default is in degrees, this can set with
unit
- unitstr, astropy.Unit
Either a string with the unit name, or an
astropy.unit.Unit
object- use_orig_xybool
If the rotation should be based on the original coordinates or the current coordinates (e.g. if rotation has already been applied)
-
scale_spectrum
(self, idx=0, mag=20, filter_name='Ks')[source]¶ Scale a certain spectrum to a certain magnitude
See
simcado.source.scale_spectrum()
for examples- Parameters
- idxint
The index of the spectrum to be scaled: <Source>.spectra[idx] Default is <Source>.spectra[0]
- magfloat
[mag] new magnitude of spectrum
- filter_namestr, TransmissionCurve
Any filter name from SimCADO or a
TransmissionCurve
object (seeget_filter_set()
)
-
scale_with_distance
(self, distance_factor)[source]¶ Scale the source for a new distance
Scales the positions and brightnesses of the
Source
object according to the ratio of the new and old distancesi.e. distance_factor = new_distance / current_distance
Warning
This does not yet take into account redshift
Todo
Implement redshift
- Parameters
- distance_factorfloat
The ratio of the new distance to the current distance i.e. distance_factor = new_distance / current_distance
Examples
>>> from simcado.source import cluster >>> >>> curr_dist = 50000 # pc, i.e. LMC >>> new_dist = 770000 # pc, i.e. M31 >>> src = cluster(distance=curr_dist) >>> src.scale_with_distance( new_dist/curr_dist )
-
shift
(self, dx=0, dy=0, use_orig_xy=False)[source]¶ Shifts the coordinates of the source by (dx, dy) in [arcsec]
- Parameters
- dx, dyfloat, array
[arcsec] The offsets for each coordinate in the arrays
x
,y
. - If dx, dy are floats, the same offset is applied to all coordinates - If dx, dy are arrays, they must be the same length asx
,y
- use_orig_xybool
If the shift should be based on the original coordinates or the current coordinates (e.g. if shift has already been applied)
-
write
(self, filename)[source]¶ Write the current Source object out to a FITS file
- Parameters
- filenamestr
where to save the FITS file
Notes
Just a place holder so that I know what’s going on with the input table * The first extension [0] contains an “image” of size 4 x N where N is the number of sources. The 4 columns are x, y, ref, weight. * The second extension [1] contains an “image” with the spectra of all sources. The image is M x len(spectrum), where M is the number of unique spectra in the source list. M = max(ref) - 1