Introduction:
The module poppy implements an object-oriented system for modeling physical optics propagation with diffraction, particularly for telescopic and coronagraphic imaging. Right now only image and pupil planes are supported; intermediate planes are a future goal. It also
This code makes use of the python standard module logging for output information. Top-level details of the calculation are output at level logging.INFO, while details of the propagation through each optical plane are printed at level logging.DEBUG. See the Python logging documentation for an explanation of how to redirect the poppy logger to the screen, a textfile, or any other log destination of your choice.
Poppy also includes a system for modeling a complete instrument (including optical propagation, synthetic photometry, and pointing jitter), and a variety of useful utility functions for analysing and plotting PSFs, documented below.
Key Concepts:
To model optical propagation, Poppy implements an object-oriented system for representing an optical train. There are a variety of OpticalElement classes representing both physical elements as apertures, mirrors, and apodizers, and also implicit operations on wavefronts, such as rotations or tilts. Each OpticalElement may be defined either via analytic functions (e.g. a simple circular aperture) or by numerical input FITS files (e.g. the complex JWST aperture with appropriate per-segment WFE). A series of such OpticalElements is chained together in order in an OpticalSystem class. That class is capable of generating Wavefronts (another class) suitable for propagation through the desired elemeents (with correct array size and sampling), and then performing the optical propagation onto the final image plane.
There is an even higher level class Instrument which adds support for selectable instrument mechanisms (such as filter wheels, pupil stops, etc). In particular it adds support for computing via synthetic photometry the appropriate weights for multiwavelength computations through a spectral bandpass filter, and for PSF blurring due to pointing jitter (neither of which effects are modeled by OpticalSystem). Given a specified instrument configuration, an appropriate OpticalSystem is generated, the appropriate wavelengths and weights are calculated based on the bandpass filter and target source spectrum, the PSF is calculated, and optionally is then convolved with a blurring kernel due to pointing jitter. All of the WebbPSF instruments are implemented by subclassing poppy.Instrument.
Poppy presently assumes that optical propagation can be modeled using Fraunhofer diffraction (far-field), such that the relationship between pupil and image plane optics is given by two-dimensional Fourier transforms. Fresnel propagation is not currently supported.
Two different algorithmic flavors of Fourier transforms are used in Poppy. The familiar FFT algorithm is used for transformations between pupil and image planes in the general case. This algorithm is relatively fast (O(N log(N))) but imposes strict constraints on the relative sizes and samplings of pupil and image plane arrays. Obtaining fine sampling in the image plane requires very large oversized pupil plane arrays and vice versa, and image plane pixel sampling becomes wavelength dependent. To avoid these constraints, for transforms onto the final Detector plane, instead a Matrix Fourier Transform (MFT) algorithm is used (See Soummer et al. 2007 Optics Express). This allows computation of the PSF directly on the desired detector pixel scale or an arbitrarily finely subsampled version therof. For equivalent array sizes N, the MFT is slower than the FFT(O(N^3)), but in practice the ability to freely choose a more appropriate N (and to avoid the need for post-FFT interpolation onto a common pixel scale) more than makes up for this and the MFT is faster.
digraph inheritancec987e26b1b { rankdir=LR; size="8.0, 12.0"; "poppy.poppy.OpticalElement" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.IdealFieldStop" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.AnalyticOpticalElement" -> "poppy.poppy.IdealFieldStop" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.instrument.Instrument" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.CircularAperture" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.AnalyticOpticalElement" -> "poppy.poppy.CircularAperture" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.poppy.FITSOpticalElement" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.OpticalElement" -> "poppy.poppy.FITSOpticalElement" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.poppy.IdealCircularOcculter" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.AnalyticOpticalElement" -> "poppy.poppy.IdealCircularOcculter" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.poppy.HexagonAperture" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.AnalyticOpticalElement" -> "poppy.poppy.HexagonAperture" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.poppy.CompoundAnalyticOptic" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.AnalyticOpticalElement" -> "poppy.poppy.CompoundAnalyticOptic" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.poppy.Rotation" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.OpticalElement" -> "poppy.poppy.Rotation" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.poppy.IdealBarOcculter" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.AnalyticOpticalElement" -> "poppy.poppy.IdealBarOcculter" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.poppy.Detector" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.OpticalElement" -> "poppy.poppy.Detector" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.poppy.AnalyticOpticalElement" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.OpticalElement" -> "poppy.poppy.AnalyticOpticalElement" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.poppy.BandLimitedCoron" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.AnalyticOpticalElement" -> "poppy.poppy.BandLimitedCoron" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.poppy.SquareAperture" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.AnalyticOpticalElement" -> "poppy.poppy.SquareAperture" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.poppy.FQPM_FFT_aligner" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.AnalyticOpticalElement" -> "poppy.poppy.FQPM_FFT_aligner" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.poppy.IdealFQPM" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.AnalyticOpticalElement" -> "poppy.poppy.IdealFQPM" [arrowsize=0.5,style="setlinewidth(0.5)"]; "poppy.poppy.Wavefront" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; "poppy.poppy.OpticalSystem" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname=Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans,height=0.25]; }
A class representing a monochromatic wavefront that can be transformed between pupil and image planes (but not to intermediate planes, yet).
Use the wf.propagateTo() method to transform a wavefront between conjugate planes. This will update those properties as appropriate.
By default, Wavefronts are created in a pupil plane. Set pixelscale=# to make an image plane instead.
| Parameters : | wavelength : float
npix : int
diam : float, optional
pixelscale : float, optional
oversample : int, optional
dtype : numpy.dtype, optional
|
|---|
Methods
Wavelength in meters
Diameter in meters. Applies to a pupil plane only.
Field of view in arcsec. Applies to an image plane only.
Pixel scale, in arcsec/pixel or meters/pixel depending on plane type
List of strings giving a descriptive history of actions performed on the wavefront. Saved to FITS headers.
A descriptive string for where a wavefront is instantaneously located (e.g. ‘before occulter’). Used mostly for titling displayed plots.
Return a copy of the wavefront as a different object.
Set this wavefront’s total intensity to 1
Return a wavefront as a pyFITS HDUList object
| Parameters : | what : string
|
|---|
Write a wavefront to a FITS file.
| Parameters : | filename : string
what : string
clobber : bool, optional
|
|---|---|
| Returns : | outfile: file on disk :
|
Display wavefront on screen
| Parameters : | what : string
nrows : int
row : int
imagecrop : float, optional
showpadding : bool, optional
colorbar : bool
ax : matplotlib Axes
|
|---|---|
| Returns : | figure : matplotlib figure
|
Propagates a wavefront object to the next optic in the list. Modifies this wavefront object itself.
Tilt a wavefront in X and Y.
Recall from Fourier optics (although this is straightforwardly rederivable by drawing triangles...) that for a wavefront tilted by some angle theta in radians, for a point r meters from the center of the pupil has
extra_pathlength = sin(theta) * r extra_waves = extra_pathlength/ wavelength = r * sin(theta) / wavelength
So we calculate the U and V arrays (corresponding to r for the pupil, in meters from the center) and then multiply by the appropriate trig factors for the angle.
The sign convention is chosen such that positive Yangle tilts move the star upwards in theg array at the focal plane. (This is sort of an inverse of what physically happens in the propagation to or through focus, but we’re ignoring that here and trying to just work in sky coords)
| Parameters : | Xangle, Yangle : float
|
|---|
Rotate a wavefront by some amount
| Parameters : | angle : float
|
|---|
Return Y, X coordinates for this wavefront, in the manner of numpy.indices()
This function knows about the offset resulting from FFTs. Use it whenever computing anything measures in wavefront coordinates.
| Returns : | Y, X : array_like
|
|---|
A class representing a series of optical elements, either Pupil, Image, or Detector planes, through which light can be propagated.
The difference between Image and Detector planes is that Detectors have fixed pixels in terms of arcsec/pixel regardless of wavelength (computed via MFT) while Image planes have variable pixels scaled in terms of lambda/D. Pupil planes are some fixed size in meters, of course.
| Parameters : | name : string
oversample : int
verbose : bool
|
|---|
Methods
Add a pupil plane optic from file(s) giving transmission or OPD
- from file(s) giving transmission and/or OPD
[set arguments transmission=filename and/or opd=filename]
- from an analytic function
[set function=’Circle’, ‘Square’ and set additional kwargs to define shape etc.
- from an already-created OpticalElement object
[set optic=that object]
| Parameters : | optic : poppy.OpticalElement, optional
function: string, optional :
opd, transmission : string, optional
Note: Now you can use the optic arguement for either an OpticalElement or a string function name, : and it will do the right thing depending on type. Both existing arguments are left for back compatibility for now. : Any provided parameters are passed to :ref:`OpticalElement`. : |
|---|
Add an image plane optic, either
- from file(s) giving transmission or OPD
[set arguments transmission=filename and/or opd=filename]
- from an analytic function
[set function=’circle, fieldstop, bandlimitedcoron, or FQPM’ and set additional kwargs to define shape etc.
- from an already-created OpticalElement object
[set optic=that object]
| Parameters : | optic : poppy.OpticalElement
function: string :
opd, transmission : string
Note: Now you can use the optic arguement for either an OpticalElement or a string function name, : and it will do the right thing depending on type. Both existing arguments are left for back compatibility for now. : |
|---|
Add a Detector object to an optical system. By default, use the same oversampling as the rest of the optical system, but the user can override to a different value if desired by setting oversample.
Other arguments are passed to the init method for Detector().
| Parameters : | pixelscale : float
oversample : int, optional
|
|---|
Create a Wavefront object suitable for sending through a given optical system, based on the size of the first optical plane, assumed to be a pupil.
If the first optical element is an Analytic pupil (i.e. has no pixel scale) then an array of 1024x1024 will be created (not including oversampling).
Uses self.source_offset to assign an off-axis tilt, if requested.
| Parameters : | wavelength : float
|
|---|---|
| Returns : | wavefront : poppy.Wavefront instance
|
Propagate a wavefront through some number of optics. Returns a pyfits.HDUList object.
| Parameters : | wavelength : float
normalize : string, {‘first’, ‘last’}
display_intermediates : bool
save_intermediates : bool
poly_weight : float
|
|---|
Calculate a multi-wavelength PSF over some weighted sum of wavelengths, across multiple processors
This version uses Python’s multiprocessing package to span tasks across available processor cores.
Any additional kwargs will be passed on to propagate_mono()
| Parameters : | source : dict
nprocesses : int
save_intermediates : bool
|
|---|---|
| Returns : | outfits : :
|
Notes
Don’t set the number of processes too high, or your machine will exhaust its memory and grind to a halt. Figure about 600 MB per process, due to the large arrays used in the FFTing, and multiple copies etc. Yes, this is surprisingly high overhead - but the raw wavefront array is typically 2048^2 * 16 (since it’s type dcomplex) = 67 MB, so this is just ~ 8-10 copies of the array floating arround. TODO: be more clever and efficient with all this somehow? TODO: write an auto tool to optimize the number of processes automatically?
Calculate a PSF, either - multi-wavelength PSF over some weighted sum of wavelengths (if you provide a source parameter) - monochromatic (if you provide just a wavelen= parameter)
Any additional kwargs will be passed on to propagate_mono()
| Parameters : | source : dict
wavelen : float, optional
save_intermediates : bool
save_intermediate_what : string
return_intermediates: bool :
display : bool
|
|---|---|
| Returns : | outfits : :
|
Display all elements in an optical system on screen.
Any extra arguments are passed to the optic.display() methods of each element.
Base class for all optical elements, whether from FITS files or analytic functions.
If instantiated on its own, this just produces a null optical element (empty space, i.e. an identity function on transmitted wavefronts.) Use one of the many subclasses to create a nontrivial optic.
The OpticalElement class follows the behavoior of the Wavefront class, using units of meters/pixel in pupil space and arcsec/pixel in image space.
| Parameters : | name : string
verbose : bool
planetype : int
oversample : int
|
|---|
Methods
string. Descriptive Name of this optic
Compute a complex phasor from an OPD, given a wavelength.
The returned value should be the complex phasor array as appropriate for multiplying by the wavefront amplitude.
| Parameters : | wave : float or obj
|
|---|
Display plots showing an optic’s transmission and OPD
| Parameters : | what : str
ax : matplotlib.Axes
|
|---|
Defines an arbitrary optic, based on amplitude transmission and/or OPD FITS files.
This optic could be a pupil or field stop, an aberrated mirror, a phase mask, etc. The FITSOpticalElement class follows the behavior of the Wavefront class, using units of meters/pixel in pupil space and arcsec/pixel in image space.
| Parameters : | name : string
transmission, opd : string or pyfits HDUList
opdunits : string
verbose : bool
planetype : int
oversample : int
shift : tuple of floats, optional
rotation : float
pixelscale : optical str or float
*NOTE:* All mask files must be *squares*. : |
|---|
Methods
Bases: poppy.poppy.OpticalElement
Defines an abstract analytic optical element, i.e. one definable by some formula rather than by an input OPD or pupil file.
This class is useless on its own; instead use its various subclasses that implement appropriate getPhasor functions. It exists mostly to provide some behaviors & initialization common to all analytic optical elements.
| Parameters : | name, verbose, oversample, planetype : various
transmission, opd : string
|
|---|
Methods
Bases: poppy.poppy.AnalyticOpticalElement
Defines an ideal circular pupil aperture
| Parameters : | name : string
radius : float
pad_factor : float, optional
|
|---|
Methods
Bases: poppy.poppy.AnalyticOpticalElement
Defines an ideal hexagonal pupil aperture
Specify either the side length (= corner radius) or the flat-to-flat distance.
| Parameters : | name : string
flattoflat : float, optional
side : float, optional
|
|---|
Methods
Bases: poppy.poppy.AnalyticOpticalElement
Defines an ideal square pupil aperture
| Parameters : | name : string
size: float :
|
|---|
Methods
Bases: poppy.poppy.AnalyticOpticalElement
Defines an ideal square field stop
| Parameters : | name : string
size : float
angle : float
|
|---|
Methods
Bases: poppy.poppy.AnalyticOpticalElement
Defines an ideal circular occulter (opaque circle)
| Parameters : | name : string
radius : float
|
|---|
Methods
Bases: poppy.poppy.AnalyticOpticalElement
Defines an ideal bar occulter (like in MIRI’s Lyot coronagraph)
| Parameters : | name : string
width : float
angle : float
|
|---|
Methods
Bases: poppy.poppy.AnalyticOpticalElement
Defines an ideal band limited coronagraph occulting mask.
| Parameters : | name : string
kind : string
sigma : float
wavelength : float
|
|---|
Methods
Bases: poppy.poppy.AnalyticOpticalElement
Defines an ideal 4-quadrant phase mask coronagraph, with its retardance set perfectly to 0.5 waves at one specific wavelength and varying linearly on either side of that.
| Parameters : | name : string
wavelength : float
|
|---|
Methods
Bases: poppy.poppy.AnalyticOpticalElement
Helper class for modeling FQPMs accurately
Adds (or removes) a slight wavelength- and pixel-scale-dependent tilt to a pupil wavefront, to ensure the correct alignment of the image plane FFT’ed PSF with the desired quad pixel alignment for the FQPM.
This is purely a computational convenience tool to work around the pixel coordinate restrictions imposed by the FFT algorithm, not a representation of any physical optic.
| Parameters : | direction : string
|
|---|
Methods
Bases: poppy.poppy.AnalyticOpticalElement
Define a compound analytic optical element made up of the combination of two or more individual optical elements.
This is just a convenience routine for semantic organization of optics. It can be useful to keep the list of optical planes cleaner, but you can certainly just add a whole bunch of planes all in a row without using this class to group them.
All optics should be of the same plane type (pupil or image); propagation between different optics contained inside one compound is not supported.
| Parameters : | opticslist : list
|
|---|
Methods
Bases: poppy.poppy.OpticalElement
Performs a rotation of the axes in the optical train.
This is not an actual optic itself, of course, but can be used to model a rotated optic by appling a Rotation before and/or after light is incident on that optic.
This is basically a placeholder to indicate the need for a rotation at a given part of the optical train. The actual rotation computation is performed in the Wavefront object’s propagation routines.
| Parameters : | angle : float
units : ‘degrees’ or ‘radians’
|
|---|
Methods
Bases: poppy.poppy.OpticalElement
A Detector is a specialized type of OpticalElement that forces a wavefront onto a specific fixed pixelization.
| Parameters : | name : string
pixelscale : float
fov_pixels, fov_arcsec : float
oversample : int
offset : tuple (X,Y)
|
|---|
Methods
A generic astronomical instrument, composed of (1) an optical system implemented using POPPY, and (2) some defined spectral bandpass(es), 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 will likely not want to 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.
You will at a minimum want to override the following class methods:
_getOpticalSystem _getFilterList _getDefaultNLambda _getDefaultFOV _getFITSHeader
Methods
Aperture for this optical system. May be a FITS filename, pyfits.HDUList, or poppy.OpticalElement
Pupil OPD for this optical system. May be a FITS filename, or pyfits.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.
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.
| Parameters : | source_offset_r : float
source_offset_theta : float
pupil_shift_x, pupil_shift_y : float
rebin : bool
jitter : string
parity : string “even” or “odd”
|
|---|
Detector pixel scale, in arcseconds/pixel
Currently selected filter name (e.g. “F200W”)
Compute a PSF. The result can either be written to disk (set outfile=”filename”) or else will be returned as a pyfits HDUlist object.
Output sampling may be specified in one of two ways:
By default, both oversampling factors are set equal to 2.
| Parameters : | source : pysynphot.SourceSpectrum or dict
nlambda : int
monochromatic : float, optional
fov_arcsec : float
fov_pixels : int
outfile : string
oversample, detector_oversample, fft_oversample : int
rebin : bool, optional
clobber : bool
display : bool
save_intermediates, return_intermediates : bool
|
|---|---|
| Returns : | outfits : pyfits.HDUList
|
Notes
More advanced PSF computation options (pupil shifts, source positions, jitter, ...) may be set by configuring the .options dictionary attribute of this class.
Display the currently configured optical system on screen
Documentation last updated on June 25, 2012