MIRAGE Quick Start

Table of Contents: * Getting Started * Create ``mirage` input yaml files from an APT file <#make_yaml>`__ * Single image simulation * Create simulation with one command * Running simulator steps independently * Running multiple simulations * Running in series * Running in parallel * Example observation list file * Example yaml file


Getting Started

Important: Before proceeding, ensure you have set the MIRAGE_DATA environment variable to point to the directory that contains the reference files associated with MIRAGE.

[ ]:
import os

# For examining outputs
from glob import glob
from scipy.stats import sigmaclip
import numpy as np
from astropy.io import fits
from astropy.visualization import simple_norm
import matplotlib.pyplot as plt
%matplotlib inline
[3]:
# Import top level functions
from mirage import imaging_simulator
from mirage.apt import apt_inputs
from mirage.yaml import yaml_generator

# Import the three steps of the simulator.
from mirage.seed_image import catalog_seed_image
from mirage.dark import dark_prep
from mirage.ramp_generator import obs_generator

### Create mirage input yaml files from an APT file

For convenience, observing programs with multiple pointings and detectors can be simulated starting with the program’s APT file. The xml and pointings files must be exported from APT, and are then used as input into a tool that will generate a series of yaml input files.

[ ]:
apt_xml_file = 'my_apt_file.xml'
apt_pointing_file = 'my_apt_file.pointing'
[ ]:
catalogs = {'nircam': {'sw': ['nrc_ptsrc_1.cat', 'nrc_ptsrc_2.cat'],
                       'lw': ['nrc_ptsrc_lw_1.cat', 'nrc_ptsrc_lw_2.cat']},
            'niriss': ['niriss_ptsrc_1.cat', 'niriss_ptsrc_2.cat']}
backgrounds = ['low', 'medium']
[ ]:
# Create a series of data simulator input yaml files
# from APT files
yam = yaml_generator.SimInput()

# Point to the xml and pointing files from your APT proposal
yam.input_xml = apt_xml_file
yam.pointing_file = apt_pointing_file

# Output directory for the collection of yaml files
yam.output_dir = './yaml_files'

# Output directory for the simulated observations
yam.simdata_output_dir = './imaging_data/'

# Observation table that lists the source catalogs to use with each filter/observation
yam.observation_table = observation_list_file

# Output data type. "raw", "linear", or "linear,raw" for both
yam.datatype = 'raw'

# Optional parameters
yam.use_JWST_pipeline = True
yam.use_linearized_darks = True

# Create the yaml files
yam.reffile_setup()
yam.create_inputs()
[ ]:
# Look at the list of files generated
yaml_files = glob('./yaml_files/*yaml')

## Create simulation with one command

This will take several minutes to run. For a better idea of what Mirage is doing, skip down several cells to the Running simulator steps independently section

The imaging_simulator function will run all three steps of the simulator. This convenience function is useful when creating simulated imaging mode data. WFSS data will need to be run in a slightly different way.

[ ]:
test_yaml_file = yaml_files[0]
[ ]:
# imaging_simulator is a wrapper around all 3 steps of Mirage
img_sim = imaging_simulator.ImgSim()
img_sim.paramfile = test_yaml_file
img_sim.create()

## More detail on what’s going on: Running simulation steps independently

First generate the “seed image”

This is generally a 2D noiseless countrate image that contains only simulated astronomical sources.

A seed image is generated based on a .yaml file that contains all the necessary parameters for simulating data. An example .yaml file is shown at the bottom of this notebook.

[ ]:
# Generate the seed image
cat = catalog_seed_image.Catalog_seed()
cat.paramfile = test_yaml_file
cat.make_seed()

Look at the seed image

[ ]:
# Need to flip the image vertically in order to match what ds9 would show
def show(array, title ,min=None, max=None):
    plt.figure(figsize=(12,12))
    if min is None and max is None:
            norm = simple_norm(array, 'log', percent=99)
    else:
        if min is None:
            min = np.min(array)
        if max is None:
            max = np.max(array)
        norm = simple_norm(array, 'log', min_cut=min, max_cut=max)
    plt.imshow(array[::-1,:], norm=norm)
    plt.title(title)
    plt.colorbar().set_label('DN$^{-}$/s')
[ ]:
show(cat.seedimage,'Seed Image', max=100, min=0)

Prepare the dark current exposure

This will serve as the base of the simulated data. This step will linearize the dark current (if it is not already), and reorganize it into the requested readout pattern and number of groups.

[ ]:
d = dark_prep.DarkPrep()
d.paramfile = test_yaml_file
d.prepare()

Look at the dark current

[ ]:
# For this, we will look at an image of the final group minus the first group
exptime = d.linDark.header['NGROUPS'] * cat.frametime
diff = (d.linDark.data[0,-1,:,:] - d.linDark.data[0,0,:,:]) / exptime
show(diff,'Dark Current Countrate', max=0.1, min=0)

Create the final exposure

Turn the seed image into a exposure of the proper readout pattern, and combine it with the dark current exposure. Cosmic rays and other detector effects are added.

The output can be either this linearized exposure, or a ‘raw’ exposure where the linearized exposure is “unlinearized” and the superbias and reference pixel signals are added, or the user can request both outputs. This is controlled from within the yaml parameter file.

[ ]:
obs = obs_generator.Observation()
obs.linDark = d.prepDark
obs.seed = cat.seedimage
obs.segmap = cat.seed_segmap
obs.seedheader = cat.seedinfo
obs.paramfile = test_yaml_file
obs.create()

Examine the final output image

Again, we will look at the last group minus the first group

[ ]:
with fits.open(obs.raw_output) as h:
    lindata = h[1].data
    header = h[0].header
[ ]:
exptime = header['EFFINTTM']
diffdata = (lindata[0,-1,:,:] - lindata[0,0,:,:]) / exptime
show(diffdata,'Simulated Data',min=0,max=20)

If you have multiple exposures that will use the same dark current image (with the same readout pattern, subarray size, and number of groups), you can feed the output from the initial run of dark_prep into future runs of the obs_generator, to save time. This can be accomplished with the imaging_simulator.py code, as shown below. (Note that time savings are minimal in this case, where the readout pattern is RAPID and there are only a handful of groups. This means that no averaging/skipping of frames has to be done within dark_prep.py)

[ ]:
# Now that the linearized dark product has been created, if you want to use it
# when running the simulator with a different yaml file (or repeating the run
# with the same yaml file) you can provide the filename of the dark product, and the
# dark_prep step will be skipped.
# NOTE: if you use the same dark product for multiple exposures, those exposures
# will contain exactly the same dark signal. This may or may not be advisable, depending
# on your goals for the simulated data.
img_sim_same_dark = imaging_simulator.ImgSim()
img_sim_same_dark.paramfile = second_yaml_file
img_sim_same_dark.override_dark = 'jw44444001001_01101_00001_nrcb1_uncal_linear_dark_prep_object.fits'
img_sim_same_dark.create()

## Running Multiple Simulations

Each yaml file will simulate an exposure for a single pointing using a single detector.

### Function to simulate multiple detectors/pointings in series
[ ]:
def make_sim(paramlist):
    '''Function to run many simulations in series
    '''
    for file in paramlist:
        m = imaging_simulator.ImgSim()
        m.paramfile = file
        m.create()
[ ]:
from multiprocessing import Pool

n_procs = 3 # number of cores available

with Pool(n_procs) as pool:
    pool.map(make_sim, yaml_files)

## Example yaml input file

Entries listed as ‘config’ have default files that are present in the config directory of the repository. The scripts are set up to automatically find and use these files. The user can replace ‘config’ with a filename if they wish to override the default.

In general, if ‘None’ is placed in a field, then the step that uses that particular file will be skipped.

Note that the linearized_darkfile entry overrides the dark entry, unless linearized_darkfile is set to None, in which case the dark entry will be used.

```yaml Inst: instrument: NIRCam #Instrument name mode: imaging #Observation mode (e.g. imaging, WFSS, moving_target) use_JWST_pipeline: False #Use pipeline in data transformations

Readout: readpatt: RAPID #Readout pattern (RAPID, BRIGHT2, etc) overrides nframe,nskip unless it is not recognized ngroup: 3 #Number of groups in integration nint: 1 #Number of integrations per exposure array_name: NRCB5_FULL #Name of array (FULL, SUB160, SUB64P, etc) filter: F250M #Filter of simulated data (F090W, F322W2, etc) pupil: CLEAR #Pupil element for simulated data (CLEAR, GRISMC, etc)

Reffiles: #Set to None or leave blank if you wish to skip that step dark: None #Dark current integration used as the base linearized_darkfile: $MIRAGE_DATA/nircam/darks/linearized/B5/Linearized_Dark_and_SBRefpix_NRCNRCBLONG-DARK-60090141241_1_490_SE_2016-01-09T02h46m50_uncal.fits # Linearized dark ramp to use as input. Supercedes dark above badpixmask: $MIRAGE_DATA/nircam/reference_files/badpix/NRCB5_17161_BPM_ISIMCV3_2016-01-21_ssbspmask_DMSorient.fits # If linearized dark is used, populate output DQ extensions using this file superbias: $MIRAGE_DATA/nircam/reference_files/superbias/NRCB5_superbias_from_list_of_biasfiles.list.fits #Superbias file. Set to None or leave blank if not using linearity: $MIRAGE_DATA/nircam/reference_files/linearity/NRCBLONG_17161_LinearityCoeff_ADU0_2016-05-22_ssblinearity_v2_DMSorient.fits #linearity correction coefficients saturation: $MIRAGE_DATA/nircam/reference_files/saturation/NRCB5_17161_WellDepthADU_2016-03-10_ssbsaturation_wfact_DMSorient.fits #well depth reference files gain: $MIRAGE_DATA/nircam/reference_files/gain/NRCB5_17161_Gain_ISIMCV3_2016-02-25_ssbgain_DMSorient.fits #Gain map pixelflat: None illumflat: None #Illumination flat field file astrometric: $MIRAGE_DATA/nircam/reference_files/distortion/NRCB5_FULL_distortion.asdf #Astrometric distortion file (asdf) ipc: $MIRAGE_DATA/nircam/reference_files/ipc/NRCB5_17161_IPCDeconvolutionKernel_2016-03-18_ssbipc_DMSorient.fits #File containing IPC kernel to apply invertIPC: True #Invert the IPC kernel before the convolution. True or False. Use True if the kernel is designed for the removal of IPC effects, like the JWST reference files are. occult: None #Occulting spots correction image pixelAreaMap: $MIRAGE_DATA/nircam/reference_files/pam/NIRCam_B5_PAM_imaging.fits #Pixel area map for the detector. Used to introduce distortion into the output ramp. subarray_defs: config #File that contains a list of all possible subarray names and coordinates readpattdefs: config #File that contains a list of all possible readout pattern names and associated NFRAME/NSKIP values crosstalk: config #File containing crosstalk coefficients filtpupilcombo: config #File that lists the filter wheel element / pupil wheel element combinations. Used only in writing output file flux_cal: config #File that lists flux conversion factor and pivot wavelength for each filter. Only used when making direct image outputs to be fed into the grism disperser code.

nonlin: limit: 60000.0 #Upper singal limit to which nonlinearity is applied (ADU) accuracy: 0.000001 #Non-linearity accuracy threshold maxiter: 10 #Maximum number of iterations to use when applying non-linearity robberto: False #Use Massimo Robberto type non-linearity coefficients

cosmicRay: path: $MIRAGE_DATA/nircam/cosmic_ray_library/ #Path to CR library library: SUNMIN #Type of cosmic rayenvironment (SUNMAX, SUNMIN, FLARE) scale: 1.5 #Cosmic ray scaling factor suffix: IPC_NIRCam_B5 #Suffix of library file names seed: 2956411739 #Seed for random number generator

simSignals: pointsource: my_point_sources.cat #File containing a list of point sources to add (x,y locations and magnitudes) psfpath: $MIRAGE_DATA/nircam/webbpsf_library/ #Path to PSF library psfbasename: nircam #Basename of the files in the psf library psfpixfrac: 0.25 #Fraction of a pixel between entries in PSF library (e.g. 0.25 = files for PSF centered at 0.25 pixel intervals within pixel) psfwfe: predicted #PSF WFE value (“predicted” or “requirements”) psfwfegroup: 0 #WFE realization group (0 to 4) galaxyListFile: my_galaxies_catalog.list extended: None #Extended emission count rate image file name extendedscale: 1.0 #Scaling factor for extended emission image extendedCenter: 1024,1024 #x,y pixel location at which to place the extended image if it is smaller than the output array size PSFConvolveExtended: True #Convolve the extended image with the PSF before adding to the output image (True or False) movingTargetList: None #Name of file containing a list of point source moving targets (e.g. KBOs, asteroids) to add. movingTargetSersic: None #ascii file containing a list of 2D sersic profiles to have moving through the field movingTargetExtended: None #ascii file containing a list of stamp images to add as moving targets (planets, moons, etc) movingTargetConvolveExtended: True #convolve the extended moving targets with PSF before adding. movingTargetToTrack: None #File containing a single moving target which JWST will track during observation (e.g. a planet, moon, KBO, asteroid) This file will only be used if mode is set to “moving_target” zodiacal: None #Zodiacal light count rate image file zodiscale: 1.0 #Zodi scaling factor scattered: None #Scattered light count rate image file scatteredscale: 1.0 #Scattered light scaling factor bkgdrate: 0.0 #Constant background count rate (electrons/sec/pixel) poissonseed: 2012872553 #Random number generator seed for Poisson simulation) photonyield: True #Apply photon yield in simulation pymethod: True #Use double Poisson simulation for photon yield

Telescope: ra: 53.1 #RA of simulated pointing dec: -27.8 #Dec of simulated pointing rotation: 0.0 #y axis rotation (degrees E of N)

newRamp: dq_configfile: config #config file used by JWST pipeline sat_configfile: config #config file used by JWST pipeline superbias_configfile: config #config file used by JWST pipeline refpix_configfile: config #config file used by JWST pipeline linear_configfile: config #config file used by JWST pipeline

Output: file: jw44444024002_01101_00001_nrcb1_uncal.fits #Output filename directory: ./ # Directory in which to place output files datatype: linear,raw # Type of data to save. ‘linear’ for linearized ramp. ‘raw’ for raw ramp. ‘linear,raw’ for both format: DMS #Output file format Options: DMS, SSR(not yet implemented) save_intermediates: False #Save intermediate products separately (point source image, etc) grism_source_image: False # Create an image to be dispersed? unsigned: True #Output unsigned integers? (0-65535 if true. -32768 to 32768 if false) dmsOrient: True #Output in DMS orientation (vs. fitswriter orientation). program_number: 44444 #Program Number title: Supernovae and Black Holes Near Hyperspatial Bypasses #Program title PI_Name: Doug Adams #Proposal PI Name Proposal_category: GO #Proposal category Science_category: Cosmology #Science category observation_number: ‘024’ #Observation Number observation_label: Obs2 #User-generated observation Label visit_number: ‘002’ #Visit Number visit_group: ‘01’ #Visit Group visit_id: ‘42424024002’ #Visit ID sequence_id: ‘1’ #Sequence ID activity_id: ‘01’ #Activity ID. Increment with each exposure. exposure_number: ‘00001’ #Exposure Number obs_id: ‘V44444024002P0000000001101’ #Observation ID number date_obs: ‘2019-10-15’ #Date of observation time_obs: ‘06:29:11.852’ #Time of observation obs_template: ‘NIRCam Imaging’ #Observation template primary_dither_type: NONE #Primary dither pattern name total_primary_dither_positions: 1 #Total number of primary dither positions primary_dither_position: 1 #Primary dither position number subpix_dither_type: 2-POINT-MEDIUM-WITH-NIRISS #Subpixel dither pattern name total_subpix_dither_positions: 2 #Total number of subpixel dither positions subpix_dither_position: 2 #Subpixel dither position number xoffset: 344.284 #Dither pointing offset in x (arcsec) yoffset: 466.768 #Dither pointing offset in y (arcsec)

```

[ ]: