The Roman Space Telescope Module

The galsim.roman module contains information and functionality that can be used to simulate images for the Roman Space Telescope. Some of the functionality is specific to Roman per se, but some of the routines are more generically simulating aspects of the HgCdTe detectors, which will be used on Roman. These routines might therefore be useful for simulating observations from other telescopes that will use these detectors.

The demo script demo13.py illustrates the use of most of this functionality.

Note

To use this module, you must separately import galsim.roman. These functions are not automatically imported when you import galsim.

Module-level Attributes

There are a number of attributes of the galsim.roman module, which define some numerical parameters related to the Roman geometry. Some of these parameters relate to the entire wide-field imager. Others, especially the return values of the functions to get the PSF and WCS, are specific to each SCA (Sensor Chip Assembly, the equivalent of a chip for an optical CCD) and therefore are indexed based on the SCA. All SCA-related arrays are 1-indexed, i.e., the entry with index 0 is None and the entries from 1 to n_sca are the relevant ones. This is consistent with diagrams and so on provided by the Roman project, which are 1-indexed.

The NIR detectors that will be used for Roman have a different photon detection process from CCDs. In particular, the photon detection process begins with charge generation. However, instead of being read out along columns (as for CCDs), they are read directly from each pixel. Moreover, the actual quantity that is measured is technically not charge, but rather voltage. The charge is inferred based on the capacitance. To use a common language with that for CCDs, we will often refer to quantities measured in units of e-/pixel, but for some detector non-idealities, it is important to keep in mind that it is voltage that is sensed.

gain

The gain for all SCAs (sensor chip assemblies) is expected to be the roughly the same, and we currently have no information about how different they will be, so this is a single value rather than a list of values. Once the actual detectors exist and have been characterized, it might be updated to be a dict with entries for each SCA.

pixel_scale

The pixel scale in units of arcsec/pixel. This value is approximate and does not include effects like distortion, which are included in the WCS.

diameter

The telescope diameter in meters.

obscuration

The linear obscuration of the telescope, expressed as a fraction of the diameter.

collecting_area
The actual collecting area after accounting for obscuration, struts, etc. in

units of cm^2.

exptime

The typical exposure time in units of seconds. The number that is stored is for a single dither. Each location within the survey will be observed with a total of 5-7 dithers across 2 epochs.

n_dithers

The number of dithers per filter (typically 5-7, so this is currently 6 as a reasonable effective average).

dark_current

The dark current in units of e-/pix/s.

nonlinearity_beta

The coefficient of the (counts)^2 term in the detector nonlinearity function. This will not ordinarily be accessed directly by users; instead, it will be accessed by the convenience function in this module that defines the nonlinearity function as counts_out = counts_in + beta*counts_in^2. Alternatively users can use the galsim.roman.applyNonlinearity routine, which already knows about the expected form of the nonlinearity in the detectors.

reciprocity_alpha

The normalization factor that determines the effect of reciprocity failure of the detectors for a given exposure time. Alternatively, users can use the galsim.roman.addReciprocityFailure routine, which knows about this normalization factor already, and allows users to choose an exposure time or use the default Roman exposure time.

read_noise

A total of 8.5e-. This comes from 20 e- per correlated double sampling (CDS) and a 5 e- floor, so the CDS read noise dominates. The source of CDS read noise is the noise introduced when subtracting a single pair of reads; this can be reduced by averaging over multiple reads. Also, this read_noise value might be reduced based on improved behavior of newer detectors which have lower CDS noise.

thermal_backgrounds

The thermal backgrounds (in units of e-/pix/s) are based on a temperature of 282 K, but this plan might change in future. The thermal backgrounds depend on the band, so this is not a single number; instead, it’s a dictionary that is accessed by the name of the optical band, e.g., galsim.roman.thermal_backgrounds['F184'] (where the names of the bandpasses can be obtained using the getBandpasses routine described below).

pupil_plane_file

There is actually a separate file for each SCA giving the pupil plane mask for the Roman telescope as seen from the center of each SCA. When building the PSF with galsim.roman.getPSF, it will use the correct one for the given SCA. However, for backwards compatibility, if anyone needs a generic image of the pupil plane, this file is for SCA 2, near the center of the WFC field.

pupil_plane_scale

The pixel scale in meters per pixel for the image in pupil_plane_file.

stray_light_fraction

The fraction of the sky background that is allowed to contribute as stray light. Currently this is required to be <10% of the background due to zodiacal light, so its value is set to 0.1 (assuming a worst-case). This could be used to get a total background including stray light.

ipc_kernel

The 3x3 kernel to be used in simulations of interpixel capacitance (IPC), using galsim.roman.applyIPC.

persistence_coefficients

The retention fraction of the previous eight exposures in a simple, linear model for persistence.

persistence_fermi_params

The parameters in the fermi persistence model.

n_sca

The number of SCAs in the focal plane.

n_pix_tot

Each SCA has n_pix_tot x n_pix_tot pixels.

n_pix

The number of pixels that are actively used. The 4 outer rows and columns will be attached internally to capacitors rather than to detector pixels, and used to monitor bias voltage drifts. Thus, images seen by users will be n_pix x n_pix.

jitter_rms
The worst-case RMS jitter per axis for Roman in the current design (reality

will likely be much better than this). Units: arcsec.

charge_diffusion

The per-axis sigma to use for a Gaussian representing charge diffusion for Roman. Units: pixels.

For example, to get the gain value, use galsim.roman.gain. Most numbers related to the nature of the detectors are subject to change as further lab tests are done.

Roman Functions

This module also contains the following routines:

galsim.roman.getBandpasses

A utility to get a dictionary containing galsim.Bandpass objects for each of the Roman imaging bandpasses, which by default have AB zeropoints given using the GalSim zeropoint convention (see getBandpasses docstring for more details).

galsim.roman.getSkyLevel

A utility to find the expected sky level due to zodiacal light at a given position, in a given band.

galsim.roman.applyNonlinearity

A routine to apply detector nonlinearity of the type expected for Roman.

galsim.roman.addReciprocityFailure

A routine to include the effects of reciprocity failure in images at the level expected for Roman.

galsim.roman.applyIPC

A routine to incorporate the effects of interpixel capacitance in Roman images.

galsim.roman.applyPersistence

A routine to incorporate the effects of persistence - the residual images from earlier exposures after resetting.

galsim.roman.allDetectorEffects

A routine to add all sources of noise and all (implemented) detector effects to an image containing astronomical objects plus background. In principle, users can simply use this routine instead of separately using the various routines like galsim.roman.applyNonlinearity.

galsim.roman.getPSF

A routine to get a chromatic representation of the PSF in a single SCA.

galsim.roman.getWCS

A routine to get the WCS for each SCA in the focal plane, for a given target RA, dec, and orientation angle.

galsim.roman.findSCA

A routine that can take the WCS from getWCS and some sky position, and indicate in which SCA that position can be found, optionally including half of the gaps between SCAs (to identify positions that are in the focal plane array but in the gap between SCAs).

galsim.roman.allowedPos

A routine to check whether Roman is allowed to look at a given position on a given date, given the constraints on orientation with respect to the sun.

galsim.roman.bestPA

A routine to calculate the best observatory orientation angle for Roman when looking at a given position on a given date.

Another routine that may be necessary is galsim.utilities.interleaveImages. The Roman PSFs at native Roman pixel scale are undersampled. A Nyquist-sampled PSF image can be obtained by a two-step process:

  1. Call the galsim.roman.getPSF routine and convolve the PSF with the Roman pixel response to get the effective PSF.

  2. Draw the effective PSF onto an Image using drawImage routine, with a pixel scale lesser than the native pixel scale (using the ‘method=no_pixel’ option).

However, if pixel-level effects such as nonlinearity and interpixel capacitance must be applied to the PSF images, then they must drawn at the native pixel scale. A Nyquist-sampled PSF image can be obtained in such cases by generating multiple images with offsets (a dither sequence) and then combining them using galsim.utilities.interleaveImages.

galsim.roman.getBandpasses(AB_zeropoint=True, default_thin_trunc=True, **kwargs)[source]

Utility to get a dictionary containing the Roman ST bandpasses used for imaging.

This routine reads in a file containing a list of wavelengths and throughput for all Roman bandpasses, and uses the information in the file to create a dictionary. This file is in units of effective area (m^2), which includes the nominal mirror size and obscuration in each bandpass. We divide these by the nominal roman.collecting_area, so the bandpass objects include both filter transmission losses and the obscuration differences relevant for each bandpass. I.e. you should always use roman.collecting_area for the collecting area in any flux calculation, and the bandpass will account for the differences from this.

In principle it should be possible to replace the version of the file with another one, provided that the format obeys the following rules:

  • There is a column called ‘Wave’, containing the wavelengths in microns.

  • The other columns are labeled by the name of the bandpass.

The bandpasses can be either truncated or thinned before setting the zero points, by passing in the keyword arguments that need to get propagated through to the Bandpass.thin() and/or Bandpass.truncate() routines. Or, if the user wishes to thin and truncate using the defaults for those two routines, they can use default_thin_trunc=True. This option is the default, because the stored ‘official’ versions of the bandpasses cover a wide wavelength range. So even if thinning is not desired, truncation is recommended.

By default, the routine will set an AB zeropoint (unless AB_zeropoint=False). The zeropoint in GalSim is defined such that the flux is 1 photon/cm^2/sec through the bandpass. This differs from an instrumental bandpass, which is typically defined such that the flux is 1 photon/sec for that instrument. The difference between the two can be calculated as follows:

# Shift zeropoint based on effective collecting area in cm^2.
delta_zp = 2.5 * np.log10(galsim.roman.collecting_area)

delta_zp will be a positive number that should be added to the GalSim zeropoints to compare with externally calculated instrumental zeropoints. When using the GalSim zeropoints for normalization of fluxes, the area kwarg to drawImage can be used to get the right normalization (giving it the quantity galsim.roman.collecting_area).

This routine also loads information about sky backgrounds in each filter, to be used by the galsim.roman.getSkyLevel() routine. The sky background information is saved as an attribute in each Bandpass object.

There are some subtle points related to the filter edges, which seem to depend on the field angle at some level. This is more important for the grism than for the imaging, so currently this effect is not included in the Roman bandpasses in GalSim.

The bandpass throughput file is translated from a spreadsheet Roman_effarea_20201130.xlsx at https://roman.gsfc.nasa.gov/science/WFI_technical.html.

Example:

>>> roman_bandpasses = galsim.roman.getBandpasses()
>>> f184_bp = roman_bandpasses['F184']
Parameters:
  • AB_zeropoint – Should the routine set an AB zeropoint before returning the bandpass? If False, then it is up to the user to set a zero point. [default: True]

  • default_thin_trunc – Use the default thinning and truncation options? Users who wish to use no thinning and truncation of bandpasses, or who want control over the level of thinning and truncation, should have this be False. [default: True]

  • **kwargs – Other kwargs are passed to either Bandpass.thin or Bandpass.truncate as appropriate.

@returns A dictionary containing bandpasses for all Roman imaging filters.

galsim.roman.getSkyLevel(bandpass, world_pos=None, exptime=None, epoch=2025, date=None)[source]

Get the expected sky level for a Roman ST observation due to zodiacal light for this bandpass and position.

This routine requires Bandpass objects that were loaded by galsim.roman.getBandpasses(). That routine will have stored tables containing the sky background as a function of position on the sky for that bandpass. This routine then interpolates between the values in those tables to arbitrary positions on the sky.

The numbers that are stored in the Bandpass object bandpass are background level in units of e-/s/arcsec^2. Multiplying by the exposure time gives a result in e-/arcsec^2. The result can either be multiplied by the approximate pixel area to get e-/pix, or the result can be used with wcs.makeSkyImage() to make an image of the sky that properly includes the actual pixel area as a function of position on the detector.

The source of the tables that are being interpolated is Chris Hirata’s publicly-available Roman exposure time calculator (ETC):

Using the throughput files loaded into the bandpass object, the calculation nominally returns photons/s/arcsec^2, but the input bandpasses used internally by the ETC code include the quantum efficiency, to effectively convert to e-/s/arcsec^2. Note that in general results will depend on the adopted model for zodiacal light, and these are uncertain at the ~10% level. One must also better sample the integration in the zodiacal light calculation to match the output tables used by GalSim here.

Positions should be specified with the world_pos keyword, which must be a CelestialCoord object. If no world_pos is supplied, then the routine will use a default position that looks sensibly away from the sun.

Parameters:
  • bandpass – A Bandpass object.

  • world_pos – A position, given as a CelestialCoord object. If None, then the routine will use an ecliptic longitude of 90 degrees with respect to the sun position (as a fair compromise between 0 and 180), and an ecliptic latitude of 30 degrees with respect to the sun position (decently out of the plane of the Earth-sun orbit). [default: None]

  • exptime – Exposure time in seconds. If None, use the default Roman exposure time. [default: None]

  • epoch – The epoch to be used for estimating the obliquity of the ecliptic when converting world_pos to ecliptic coordinates. This keyword is only used if date is None, otherwise date is used to determine the epoch. [default: 2025]

  • date – The date of the observation, provided as a python datetime object. If None, then the conversion to ecliptic coordinates assumes the sun is at ecliptic coordinates of (0,0), as it is at the vernal equinox. [default: None]

Returns:

the expected sky level in e-/arcsec^2.

galsim.roman.getPSF(SCA, bandpass, SCA_pos=None, pupil_bin=4, wcs=None, n_waves=None, extra_aberrations=None, wavelength=None, gsparams=None, logger=None, high_accuracy=None, approximate_struts=None)[source]

Get a single PSF for Roman ST observations.

The user must provide the SCA and bandpass; the latter is used when setting up the pupil plane configuration and when interpolating chromatic information, if requested.

This routine carries out linear interpolation of the aberrations within a given SCA, based on the Roman Cycle 9 specification of the aberrations as a function of focal plane position, more specifically from the WebbPSF data files from webbpsf-data-1.2.1.tar.gz downloaded from https://webbpsf.readthedocs.io/en/latest/installation.html#data-install. The abberation file is webbpsf-data/WFI/wim_zernikes_cycle9.csv.

The mask images for the Roman pupil plane are available in the same WebbPSF data files. There are separate files for each SCA, since the view of the spider pattern varies somewhat across the field of view of the wide field camera. Users usually don’t need to worry about any of this, as GalSim will select the correct pupil image automatically based on the SCA and bandpass provided.

The full pupil plane images are 4096 x 4096, which use a lot of memory and are somewhat slow to use, so we normally bin them by a factor of 4 (resulting in 1024 x 1024 images). This provides enough detail for most purposes and is much faster to render than using the full pupil plane images. This bin factor is a settable parameter, called pupil_bin. If you want the more accurate, slower calculation using the full images, you can set it to 1. In the other direction, using pupil_bin=8 (resulting in a 512 x 512 image) still provides fairly reasonable results and is even faster to render. It is not generally recommended to use higher binning than that, as the diffraction spikes will become noticeably degraded.

Note

This function will cache the aperture calculation, so repeated calls with the same SCA and bandpass should be much faster after the first call, as the pupil plane will already be loaded. If you need to clear the cache for memory reasons, you may call:

galsim.roman.roman_psfs._make_aperture.clear()

to recover any memory currently being used for this cache. Of course, subsequent calls to getPSF will need to rebuild the aperture at that point.

The PSF that is returned by default will be oriented with respect to the SCA coordinates, not world coordinates as is typical in GalSim. The pupil plane has a fixed orientation with respect to the focal plane, so the PSF rotates with the telescope. To obtain a PSF in world coordinates, which can be convolved with galaxies (that are normally described in world coordinates), you may pass in a wcs parameter to this function. This will project the PSF into world coordinates according to that WCS before returning it. Otherwise, the return value is equivalent to using wcs=galim.PixelScale(galsim.roman.pixel_scale).

The calculation takes advantage of the fact that the diffraction limit and aberrations have a simple, understood wavelength-dependence. (The Roman abberation data for Cycle 9 does in fact provide aberrations as a function of wavelength, but the deviation from the expected chromatic dependence is sub-percent so we neglect it here.) For reference, the script used to parse the Zernikes given on the webpage and create the files in the GalSim repository can be found in devel/external/parse_roman_zernikes_1217.py. The resulting chromatic object can be used to draw into any of the Roman bandpasses, though the pupil plane configuration will only be correct for those bands in the same range (i.e., long- or short-wavelength bands).

For applications that require very high accuracy in the modeling of the PSF, with very limited aliasing, you may want to lower the folding_threshold in the gsparams. Otherwise very bright stars will show some reflections in the spider pattern and possibly some boxiness at the outskirts of the PSF. Using gsparams = GSParams(folding_threshold=2.e-3) generally provides good results even for very bright (e.g. mag=10) stars. In these cases, you probably also want to reduce pupil_bin somewhat from the default value of 4.

By default, no additional aberrations are included above the basic design. However, users can provide an optional keyword extra_aberrations that will be included on top of those that are part of the design. This should be in the same format as for the ChromaticOpticalPSF class, with units of waves at the fiducial wavelength, 1293 nm. Currently, only aberrations up to order 22 (Noll convention) are simulated. For Roman, the tolerance for additional aberrations was a total of 90 nanometers RMS as of mid-2015, distributed largely among coma, astigmatism, trefoil, and spherical aberrations (NOT defocus). This information might serve as a guide for reasonable extra_aberrations inputs. The reference for that number is an earlier Cycle 5 document:

http://roman.gsfc.nasa.gov/science/sdt_public/wps/references/instrument/README_AFTA_C5_WFC_Zernike_and_Field_Data.pdf

However, the default (non-extra) aberrations are from Cycle 7 material linked earlier in this docstring.

Jitter and charge diffusion are, by default, not included. Users who wish to include these can find some guidelines for typical length scales of the Gaussians that can represent these effects, and convolve the ChromaticOpticalPSF with appropriate achromatic Gaussians.

The PSFs are always defined assuming the user will specify length scales in arcsec.

Users may find they do not have to call getPSF for all objects in their simulations; for a given SCA and position within the SCA, and a given pupil plane configuration and wavelength information, it should be possible to reuse the PSFs.

Parameters:
  • SCA – Single value specifying the SCA for which the PSF should be loaded.

  • bandpass – Single string specifying the bandpass to use when defining the pupil plane configuration and/or interpolation of chromatic PSFs. You may also pass a string ‘long’ or ‘short’ for this argument, in which case, the correct pupil plane configuration will be used for long- or short-wavelength bands (F184 is long, all else is short). In this case, no interpolation can be used, since it is defined using the extent of the chosen bandpass. If wavelength is given, then bandpass may be None, which will use the short-wavelength pupil plane image.

  • SCA_pos – Single galsim.PositionD indicating the position within the SCA for which the PSF should be created. If None, the exact center of the SCA is chosen. [default: None]

  • pupil_bin – The binning to apply to the pupil plane image. (See discussion above.) [default: 4]

  • wcs – The WCS to use to project the PSF into world coordinates. [default: galsim.PixelScale(galsim.roman.pixel_scale)]

  • n_waves – Number of wavelengths to use for setting up interpolation of the chromatic PSF objects, which can lead to much faster image rendering. If None, then no interpolation is used. Note that users who want to interpolate can always set up the interpolation later on even if they do not do so when calling getPSF. [default: None]

  • extra_aberrations – Array of extra aberrations to include in the PSF model, on top of those that are part of the Roman design. These should be provided in units of waves at the fiducial wavelength of 1293 nm, as an array of length 23 with entries 4 through 22 corresponding to defocus through the 22nd Zernike in the Noll convention. [default: None]

  • wavelength – An option to get an achromatic PSF for a single wavelength, for users who do not care about chromaticity of the PSF. If None, then the fully chromatic PSF is returned. Alternatively the user should supply either (a) a wavelength in nanometers, and they will get achromatic OpticalPSF objects for that wavelength, or (b) a bandpass object, in which case they will get achromatic OpticalPSF objects defined at the effective wavelength of that bandpass. [default: False]

  • gsparams – An optional GSParams argument. See the docstring for GSParams for details. [default: None]

Returns:

A single PSF object (either a ChromaticOpticalPSF or an OpticalPSF depending on the inputs).

galsim.roman.getWCS(world_pos, PA=None, date=None, SCAs=None, PA_is_FPA=False)[source]

This routine returns a dict containing a WCS for each of the Roman SCAs (Sensor Chip Array, the equivalent of a chip in an optical CCD). The Roman SCAs are labeled 1-18, so these numbers are used as the keys in the dict. Alternatively the user can request a subset of the SCAs using the SCAs option. The basic instrument parameters used to create the WCS correspond to those in Cycle 6, which includes some significant updates from Cycle 5, including a 90 degree rotation of the focal plane axes relative to the payload axes, and two rows of SCAs are swapped.

The user must specify a position for observation, at which the center of the focal plane array will point. This must be supplied as a CelestialCoord world_pos. In general, only certain positions are observable on certain dates, and for a given position there is an optimal position angle for the observatory (with the solar panels pointed as directly towards the sun as possible). Users who are knowledgable about these details may choose to supply a position angle as PA, either for the observatory or for the focal plane (using PA_is_FPA to indicate this). But otherwise, the routine will simply choose the optimal position angle for a given date.

To fully understand all possible inputs and outputs to this routine, users may wish to consult the diagram on the GalSim wiki, https://github.com/GalSim-developers/GalSim/wiki/GalSim-Roman-module-diagrams

Parameters:
  • world_pos – A galsim.CelestialCoord indicating the position to observe at the center of the focal plane array (FPA). Note that if the given position is not observable on the given date, then the routine will raise an exception.

  • PA – A galsim.Angle representing the position angle of the observatory +Y axis, unless PA_is_FPA=True, in which case it’s the position angle of the FPA. For users to do not care about this, then leaving this as None will result in the routine using the supplied date and world_pos to select the optimal orientation for the observatory. Note that if a user supplies a PA value, the routine does not check whether this orientation is actually allowed. [default: None]

  • date – The date of the observation, as a python datetime object. If None, then the vernal equinox in 2025 will be used. [default: None]

  • PA_is_FPA – If True, then the position angle that was provided was the PA of the focal plane array, not the observatory. [default: False]

  • SCAs – A single number or iterable giving the SCAs for which the WCS should be obtained. If None, then the WCS is calculated for all SCAs. [default: None]

Returns:

A dict of WCS objects for each SCA.

galsim.roman.findSCA(wcs_dict, world_pos, include_border=False)[source]

This is a subroutine to take a dict of WCS (one per SCA) from galsim.roman.getWCS() and query which SCA a particular real-world coordinate would be located on. The position (world_pos) should be specified as a galsim.CelestialCoord. If the position is not located on any of the SCAs, the result will be None. Note that if wcs_dict does not include all SCAs in it, then it’s possible the position might lie on one of the SCAs that was not included.

Depending on what the user wants to do with the results, they may wish to use the include_border keyword. This keyword determines whether or not to include an additional border corresponding to half of the gaps between SCAs. For example, if a user is drawing a single image they may wish to only know whether a given position falls onto an SCA, and if so, which one (ignoring everything in the gaps). In contrast, a user who plans to make a sequence of dithered images might find it most useful to know whether the position is either on an SCA or close enough that in a small dither sequence it might appear on the SCA at some point. Use of include_border switches between these scenarios.

Parameters:
  • wcs_dict – The dict of WCS’s output from galsim.roman.getWCS().

  • world_pos – A galsim.CelestialCoord indicating the sky position of interest.

  • include_border – If True, then include the half-border around SCA to cover the gap between each sensor. [default: False]

Returns:

an integer value of the SCA on which the position falls, or None if the position is not on any SCA.

galsim.roman.allowedPos(world_pos, date)[source]

This routine can be used to check whether Roman would be allowed to look at a particular position (world_pos) on a given date. This is determined by the angle of this position relative to the Sun.

In general, Roman can point at angles relative to the Sun in the range 90+/-36 degrees. Obviously, pointing too close to the Sun would result in overly high sky backgrounds. It is less obvious why Roman cannot look at a spot directly opposite from the Sun (180 degrees on the sky). The reason is that the observatory is aligned such that if the observer is looking at some sky position, the solar panels are oriented at 90 degrees from that position. So it’s always optimal for the observatory to be pointing at an angle of 90 degrees relative to the Sun. It is also permitted to look within 36 degrees of that optimal position.

Parameters:
  • world_pos – A galsim.CelestialCoord indicating the position at which the observer wishes to look.

  • date – A python datetime object indicating the desired date of observation.

Returns:

True or False, indicating whether it is permitted to look at this position on this date.

galsim.roman.bestPA(world_pos, date)[source]

This routine determines the best position angle for the observatory for a given observation date and position on the sky.

The best/optimal position angle is determined by the fact that the solar panels are at 90 degrees to the position being observed, and it is best to have those facing the Sun as directly as possible. Note that if a given world_pos is not actually observable on the given date, then this routine will return None.

Parameters:
  • world_pos – A galsim.CelestialCoord indicating the position at which the observer wishes to look.

  • date – A python datetime object indicating the desired date of observation.

Returns:

the best position angle for the observatory, as a galsim.Angle, or None if the position is not observable.

galsim.roman.convertCenter(world_pos, SCA, PA=None, date=None, PA_is_FPA=False, tol=coord.Angle(2.42406840554768e-06, coord.radians))[source]

This is a simple helper routine that takes an input position world_pos that is meant to correspond to the position of the center of an SCA, and tells where the center of the focal plane array should be. The goal is to provide a position that can be used as an input to getWCS(), which wants the center of the focal plane array.

The results of the calculation are deterministic if given a fixed position angle (PA). If it’s not given one, it will try to determine the best one for this location and date, like getWCS() does.

Because of distortions varying across the focal plane, this routine has to iteratively correct its initial result based on empirical tests. The tol kwarg can be used to adjust how careful it will be, but it always does at least one iteration.

To fully understand all possible inputs and outputs to this routine, users may wish to consult the diagram on the GalSim wiki, https://github.com/GalSim-developers/GalSim/wiki/GalSim-Roman-module-diagrams

Parameters:
  • world_pos – A galsim.CelestialCoord indicating the position to observe at the center of the given SCA. Note that if the given position is not observable on the given date, then the routine will raise an exception.

  • SCA – A single number giving the SCA for which the center should be located at world_pos.

  • PA – galsim.Angle representing the position angle of the observatory +Y axis, unless PA_is_FPA=True, in which case it’s the position angle of the FPA. For users to do not care about this, then leaving this as None will result in the routine using the supplied date and world_pos to select the optimal orientation for the observatory. Note that if a user supplies a PA value, the routine does not check whether this orientation is actually allowed. [default: None]

  • date – The date of the observation, as a python datetime object. If None, then the vernal equinox in 2025 will be used. [default: None]

  • PA_is_FPA – If True, then the position angle that was provided was the PA of the focal plane array, not the observatory. [default: False]

  • tol – Tolerance for errors due to distortions, as a galsim.Angle. [default: 0.5*galsim.arcsec]

Returns:

A CelestialCoord object indicating the center of the focal plane array.

galsim.roman.applyNonlinearity(img)[source]

Applies the Roman nonlinearity function to the supplied image im.

For more information about nonlinearity, see the docstring for galsim.Image.applyNonlinearity. Unlike that routine, this one does not require any arguments, since it uses the nonlinearity function defined within the Roman module.

After calling this method, the Image instance img is transformed to include the nonlinearity.

Parameters:

img – The Image to be transformed.

galsim.roman.addReciprocityFailure(img, exptime=139.8)[source]

Accounts for the reciprocity failure for the Roman directors and includes it in the original Image img directly.

For more information about reciprocity failure, see the docstring for galsim.Image.addReciprocityFailure. Unlike that routine, this one does not need the parameters for reciprocity failure to be provided, though it still takes exposure time as an optional argument.

Parameters:
  • img – The Image to be transformed.

  • exptime – The exposure time (t) in seconds, which goes into the expression for reciprocity failure given in the docstring. If None, then the routine will use the default Roman exposure time in galsim.roman.exptime. [default: 139.8]

galsim.roman.applyIPC(img, edge_treatment='extend', fill_value=None)[source]

Applies the effect of interpixel capacitance (IPC) to the Image instance.

For more information about IPC, see the docstring for galsim.Image.applyIPC. Unlike that routine, this one does not need the IPC kernel to be specified, since it uses the IPC kernel defined within the Roman module.

Parameters:
  • img – The Image to be transformed.

  • edge_treatment – Specifies the method of handling edges and should be one of ‘crop’, ‘extend’ or ‘wrap’. See galsim.Image.applyIPC docstring for more information. [default: ‘extend’]

  • fill_value – Specifies the value (including nan) to fill the edges with when edge_treatment is ‘crop’. If unspecified or set to ‘None’, the original pixel values are retained at the edges. If edge_treatment is not ‘crop’, then this is ignored.

galsim.roman.applyPersistence(img, prev_exposures, method='fermi')[source]

This method applies either of the two different persistence models: ‘linear’ and ‘fermi’. Slew between pointings and consecutive resets after illumination are not considered.

‘linear’ persistence model

Applies the persistence effect to the Image instance by adding a small fraction of the previous exposures (up to 8) supplied as the ‘prev_exposures’ argument. For more information about persistence, see galsim.Image.applyPersistence. Unlike that routine, this one does not need the coefficients to be specified. However, the list of previous 8 exposures will have to be supplied. Earlier exposures, if supplied, will be ignored.

‘fermi’ persistence model

Applies the persistence effect to the Image instance by adding the accumulated persistence dark current of previous exposures supplied as the ‘prev_exposures’ argument. Unlike galsim.Image.applyPersistence, this one does not use constant coefficients but a fermi model plus a linear tail below half of saturation.

For more info about the fermi model, see:

http://www.stsci.edu/hst/wfc3/ins_performance/persistence/

Parameters:
  • img – The Image to be transformed.

  • prev_exposures – List of Image instances in the order of exposures, with the recent exposure being the first element. In the linear model, the exposures exceeding the limit (8 exposures) will be ignored.

  • method – The persistence model (‘linear’ or ‘fermi’) to be applied. [default: ‘fermi’]

galsim.roman.allDetectorEffects(img, prev_exposures=(), rng=None, exptime=139.8)[source]

This utility applies all sources of noise and detector effects for Roman that are implemented in GalSim. In terms of noise, this includes the Poisson noise due to the signal (sky + background), dark current, and read noise. The detector effects that are included are reciprocity failure, quantization, persistence, nonlinearity, and interpixel capacitance. It also includes the necessary factors of gain. In short, the user should be able to pass in an Image with all sources of signal (background plus astronomical objects), and the Image will be modified to include all subsequent steps in the image generation process for Roman that are implemented in GalSim. However, to include the effect of persistence, the user needs to provide a list of recent exposures (without the readout effects) and the routine returns an updated list of recent exposures.

Parameters:
  • img – The Image to be modified.

  • prev_exposures – List of Image instances in the order of exposures, with the recent exposure being the first element. [default: ()]

  • rng – An optional galsim.BaseDeviate to use for the addition of noise. If None, a new one will be initialized. [default: None]

  • exptime – The exposure time, in seconds. If None, then the Roman default exposure time will be used. [default: 139.8]

Returns:

Updated list of previous exposures Image instances.

Return type:

prev_exposures