Instrument

class poppy.Instrument(name='', *args, **kwargs)[source]

Bases: object

A generic astronomical instrument, composed of
  1. an optical system implemented using POPPY, optionally with several configurations such as selectable image plane or pupil plane stops, and
  2. some defined spectral bandpass(es) such as selectable filters, implemented using pysynphot.

This provides the capability to model both the optical and spectral responses of a given system. PSFs may be calculated for given source spectral energy distributions and output as FITS files, with substantial flexibility.

It also provides capabilities for modeling some PSF effects not due to wavefront aberrations, for instance blurring caused by pointing jitter.

This is a base class for Instrument functionality - you cannot easily use this directly, but rather should subclass it for your particular instrument of interest. Some of the complexity of this class is due to splitting up functionality into many separate routines to allow users to subclass just the relevant portions for a given task. There’s a fair amount of functionality here but the learning curve is steeper than elsewhere in POPPY.

You will at a minimum want to override the following class methods:

  • _get_optical_system
  • _get_filter_list
  • _get_default_nlambda
  • _get_default_fov
  • _get_fits_header

For more complicated systems you may also want to override:

  • _validate_config
  • _get_synphot_bandpass
  • _apply_jitter

Attributes Summary

filter Currently selected filter name (e.g.
filter_list List of available filter names for this instrument
name
options A dictionary capable of storing other arbitrary options, for extensibility.
pixelscale Detector pixel scale, in arcseconds/pixel (default: 0.025)
pupil Aperture for this optical system.
pupilopd Pupil OPD for this optical system.

Methods Summary

calc_datacube(wavelengths, *args, **kwargs) Calculate a spectral datacube of PSFs
calc_psf([outfile, source, nlambda, …]) Compute a PSF.
display() Display the currently configured optical system on screen
get_optical_system(*args, **kwargs) Return an OpticalSystem instance corresponding to the instrument as currently configured.

Attributes Documentation

filter

Currently selected filter name (e.g. F200W)

filter_list = None

List of available filter names for this instrument

name = 'Instrument'
options = {}

A dictionary capable of storing other arbitrary options, for extensibility. The following are all optional, and may or may not be meaningful depending on which instrument is selected.

source_offset_r : float
Radial offset of the target from the center, in arcseconds
source_offset_theta : float
Position angle for that offset
pupil_shift_x, pupil_shift_y : float
Relative shift of a coronagraphic pupil in X and Y, expressed as a decimal between 0.0-1.0 Note that shifting an array too much will wrap around to the other side unphysically, but for reasonable values of shift this is a non-issue.
jitter : string “gaussian” or None
Type of jitter model to apply. Currently only convolution with a Gaussian kernel of specified width jitter_sigma is implemented. (default: None)
jitter_sigma : float
Width of the jitter kernel in arcseconds (default: 0.007 arcsec)
parity : string “even” or “odd”
You may wish to ensure that the output PSF grid has either an odd or even number of pixels. Setting this option will force that to be the case by increasing npix by one if necessary.
pixelscale = 0.025

Detector pixel scale, in arcseconds/pixel (default: 0.025)

pupil = None

Aperture for this optical system. May be a FITS filename, FITS HDUList object, or poppy.OpticalElement

pupilopd = None

Pupil OPD for this optical system. May be a FITS filename, or FITS HDUList. If the file contains a datacube, you may set this to a tuple (filename, slice) to select a given slice, or else the first slice will be used.

Methods Documentation

calc_datacube(wavelengths, *args, **kwargs)[source]

Calculate a spectral datacube of PSFs

wavelengths : iterable of floats
List or ndarray or tuple of floating point wavelengths in meters, such as you would supply in a call to calc_psf via the “monochromatic” option
calc_psf(outfile=None, source=None, nlambda=None, monochromatic=None, fov_arcsec=None, fov_pixels=None, oversample=None, detector_oversample=None, fft_oversample=None, overwrite=True, display=False, save_intermediates=False, return_intermediates=False, normalize='first')[source]

Compute a PSF. The result can either be written to disk (set outfile=”filename”) or else will be returned as a FITS HDUlist object.

Output sampling may be specified in one of two ways:

  1. Set oversample=. This will use that oversampling factor beyond detector pixels for output images, and beyond Nyquist sampling for any FFTs to prior optical planes.
  2. set detector_oversample= and fft_oversample=. This syntax lets you specify distinct oversampling factors for intermediate and final planes.

By default, both oversampling factors are set equal to 2.

More advanced PSF computation options (pupil shifts, source positions, jitter, …) may be set by configuring the options dictionary attribute of this class.

source : pysynphot.SourceSpectrum or dict
specification of source input spectrum. Default is a 5700 K sunlike star.
nlambda : int
How many wavelengths to model for broadband? The default depends on how wide the filter is: (5,3,1) for types (W,M,N) respectively
monochromatic : float, optional
Setting this to a wavelength value (in meters) will compute a monochromatic PSF at that wavelength, overriding filter and nlambda settings.
fov_arcsec : float
field of view in arcsec. Default=5
fov_pixels : int
field of view in pixels. This is an alternative to fov_arcsec.
outfile : string
Filename to write. If None, then result is returned as an HDUList
oversample, detector_oversample, fft_oversample : int
How much to oversample. Default=4. By default the same factor is used for final output pixels and intermediate optical planes, but you may optionally use different factors if so desired.
overwrite : bool
overwrite output FITS file if it already exists?
display : bool
Whether to display the PSF when done or not.
save_intermediates, return_intermediates : bool
Options for saving to disk or returning to the calling function the intermediate optical planes during the propagation. This is useful if you want to e.g. examine the intensity in the Lyot plane for a coronagraphic propagation.
normalize : string
Desired normalization for output PSFs. See doc string for OpticalSystem.calc_psf. Default is to normalize the entrance pupil to have integrated total intensity = 1.
outfits : fits.HDUList
The output PSF is returned as a fits.HDUlist object. If outfile is set to a valid filename, the output is also written to that file.
display()[source]

Display the currently configured optical system on screen

get_optical_system(*args, **kwargs)[source]

Return an OpticalSystem instance corresponding to the instrument as currently configured.