The METIS-WCU mode in Scopesim#
This notebook documents the implementation of the METIS WCU in Scopesim. It is largely a stand-alone effect that requires no further Source object. In addition to the normal “sky” modes there are new wcu modes that combine the WCU with the existing subsystems of METIS. Currently, these are
wcu_img_lmwcu_img_nwcu_lss_lwcu_lss_mwcu_lss_nwcu_lms
The most basic way to run this mode is as follows
import scopesim as sim
cmd = sim.UserCommands(use_instrument="METIS", set_modes=["wcu_lss_m"])
metis = sim.OpticalTrain(cmd)
metis.observe()
hdul = metis.readout()[0]
However, this sequence will saturate the detector. In the following we will go through the sequence and show how to change parameters to get a satisfactory result. We will use wcu_img_lm for that purpose.
Setting up the OpticalTrain#
import numpy as np
from astropy import units as u
from matplotlib import pyplot as plt
import scopesim as sim
# Edit this path if you have a custom install directory, otherwise comment it out. [For ReadTheDocs only]
sim.link_irdb("../../../")
#sim.download_packages(["METIS"])
cmd = sim.UserCommands(use_instrument="METIS", set_modes=["wcu_img_lm"])
metis = sim.OpticalTrain(cmd)
Let’s have look at the effects that are included in the OpticalTrain:
metis['reference_pixel_mask'].include = False # switch off to check consistency with previous versions
metis.effects.pprint_all()
There are four effects associated with the WCU; these replace the atmosphere and telescope effects that are used in the sky modes. wcu_relay_optics is a static effect that describes the mirrors within the WCU and does not need to be modified by most users. pupil_masks is used to select a pupil mask; this changes the flux throughput factor but more importantly changes the PSF that will be used. This is described in the notebook demo_metis_wcu_psfs.ipynb and we will not look at this effect here.
The effect that we are concerned with here is wcu_source, which has many user-settable parameters. For convenience we will assign a variable to reference it, but all operations below could be done directly to metis['wcu_source'].
wcu = metis['wcu_source']
An overview of the setup of the WCU can be obtained with
print(wcu)
As can be seen, the default configuration uses the black-body continuum lamp at a temperature of 1000 Kelvin, with the flux-controlling mask (“BlackBody aperture”) fully open. No focal-plane mask is inserted, which means that we are observing a flat field. The thermal backgrounds from the integrating sphere and the WCU itself are both set at a temperature of 300 Kelvin, so rather warm compared to typical atmospheric temperatures.
No Source is required (unlike for most sky observations) and we can create the (noise-less) image plane just ahead of the detector with
metis.observe()
fig, ax = plt.subplots()
img = ax.imshow(metis.image_planes[0].data, origin="lower")
fig.colorbar(img, label="ImagePlane [ph/s]");
As expected, no structure is seen in the ImagePlane. The expected photon flux per pixel is around 6e8 photons per second. This flux would saturate the Hawaii2RG detector of the LM imager within a minimum DIT (1.3 seconds in slow mode, 0.04 seconds in fast mode), so a neutral-density filter has to be inserted. Repeating the observation gives a manageable image-plane flux:
metis['nd_filter_wheel'].change_filter("ND_OD4")
metis.observe()
print(f"ImagePlane signal: {metis.image_planes[0].data.mean():.2g} ph/s")
hdul = metis.readout(exptime=1.3)[0]
mean = hdul[1].data.mean()
std = hdul[1].data.std()
print(f"Mean: {mean:7.1f} ADU")
print(f"Standard deviation: {std:7.1f} ADU")
fig, ax = plt.subplots()
img = ax.imshow(hdul[1].data, vmin=mean-2*std, vmax=mean+2*std)
fig.colorbar(img, label="ADU");
Changing the black-body temperature#
The method .set_temperature is used to set the temperature of the black-body source as well as the temperatures of the integrating sphere and the ambient temperature of the WCU. The latter two determine the thermal background, the former the actual signal that we want to record. In this illustrative example we remove the backgrounds and see how the temperature affects the signal level.
wcu.set_temperature(is_temp=0*u.K, wcu_temp=0*u.K)
temps = np.array([400, 600, 800, 1000, 1200]) * u.K
signal = np.zeros(5)
for i, temp in enumerate(temps):
wcu.set_temperature(bb_temp=temp)
metis.observe()
signal[i] = np.median(metis.image_planes[0].data)
_, ax = plt.subplots()
ax.plot(temps, signal, 'o')
ax.set_xlabel("Black-body temperature [K]")
ax.set_ylabel("Image-plane flux [ph/s/pixel]");
Using the flux-controlling masks#
The flux-controlling masks sit in front of the black-body source and control the flux that is entering the integrating sphere. These masks will be used for linearity measurements. It is currently not clear how the diameters of the masks translate to the amount of flux passed through. In the current implementation the method .set_bb_aperture accepts a number between 0 and 1 to describe directly the fraction of flux let through (values larger than 1 or less than 0 are clipped to 1 and 0, respectively).
wcu.set_temperature(bb_temp=1000*u.K, is_temp=300*u.K, wcu_temp=300*u.K)
bb_ap = np.array([0., 0.2, 0.4, 0.6, 0.8, 1.])
signal = np.zeros(6)
for i, ap in enumerate(bb_ap):
wcu.set_bb_aperture(ap)
metis.observe()
signal[i] = metis.image_planes[0].data.mean()
_, ax = plt.subplots()
ax.plot(bb_ap, signal, 'o')
ax.set_xlabel("Fraction of flux transmitted into IS")
ax.set_ylabel("Image-plane flux [ph/s/pixel]");
# check
assert np.all(np.diff(signal) > 0)
Changing the focal-plane mask#
So far, we have produced flat-field images with the ‘open’ focal-plane mask. The WCU has masks for a single pinhole and grids of pinholes. These can be applied using the .set_fpmask method. This method also allows to rotate and shift the mask (which in METIS would be done with the derotator and the internal chopper, respectively).
Information on the currently inserted mask can be obtained with
print(wcu.fpmask)
wcu.set_fpmask("grid_lm", angle=20, shift=(1, 0.5))
print(wcu.fpmask)
metis.observe()
_, ax = plt.subplots()
ax.imshow(metis.image_planes[0].data, norm="log", origin="lower");
In the example above, the name of a pre-defined mask was used (these are provided in the METIS instrument package as wcu/fp_mask_*.dat). To simplify testing of mask designs the method also accepts a valid path to a filename. To create your own mask start by copying one of the existing pinhole or grid masks and enter the hole positions in the x and y columns (in arcsec using a plate scale of 3.319 mm/arcsec) and hole diameter in the diam column (also in arcsec). Note that holes are currently assumed to be non-resolved and are modelled as a single pixel in a 2047 by 2047 image (diam scales the brightness of the flux passed through the hole).
Changing the lamp#
In addition to the black-body source (bb), the WCU includes three laser sources. As there is no observing mode that covers more than one of these lasers at a time, the implementation in Scopesim treats all three as a single laser lamp with lines in the L, M and N bands. The behaviour of the tunable laser in the M band is currently not clear and the implementation in Scopesim is not much more than a preliminary place holder. To change the lamp use the .set_lamp method:
wcu.set_lamp('laser')
To show the behaviour we shall switch to the long-slit L-band mode with the grid mask and simulate spectra with both the black-body and the laser sources. To subtract the thermal background emission, we also take an exposure with the black-body source closed.
cmd = sim.UserCommands(use_instrument='METIS', set_modes=['wcu_lss_l'])
metis = sim.OpticalTrain(cmd)
wcu = metis['wcu_source']
# Using the black-body source
wcu.set_fpmask('grid_lm')
print(wcu)
metis.observe()
implane_bb = metis.image_planes[0].data
# Background, no source
wcu.set_bb_aperture(0)
print(wcu)
metis.observe()
implane_off = metis.image_planes[0].data
wcu.set_bb_aperture(1)
# Using the laser source
wcu.set_lamp('laser')
print(wcu)
metis.observe()
implane_laser = metis.image_planes[0].data
fig, (ax0, ax1) = plt.subplots(1, 2, figsize=(10, 5))
ax0.imshow(implane_bb - implane_off, origin="lower", norm="symlog")
ax0.set_title("LSS_L, black-body (ImagePlane)")
ax1.imshow(implane_laser - implane_off, origin="lower", norm="symlog")
ax1.set_title("LSS_L, laser (ImagePlane)");
Using a configuration file#
The default configuration for the WCU is read from the file metis_wcu_config.yaml in the METIS instrument package. This is a yaml file that contains both the user-settable parameters described above and also physical quantities that describe the WCU itself:
# ------------- User-settable parameters
lamps: ["bb", "laser", "off"] # available lamps
current_lamp: "bb" # the lamp currently in use
bb_temp: 1000 # [K] temperature of BB source # Kelvin!
is_temp: 300 # [K] temperature of the integratig sphere
wcu_temp: 300 # [K] ambient temperature in the WCU
bb_aperture: 1.0 # aperture of flux-controlling mask
fpmasks: ["open", "pinhole_lm", "pinhole_n", "grid_lm"]
fpmask_filename_format: "wcu/fp_mask_{}.dat"
current_fpmask: "open"
fpmask_angle: 0
fpmask_shift: [0, 0]
# ------------- Data needed to describe the WCU
bb_to_is: "wcu/WCU_BB_to_IS_throughput.fits"
is_reflect: "wcu/WCU_IS_reflectivity.dat"
tube_reflect: "wcu/WCU_tube_reflectivity.dat"
mask_reflect: "wcu/WCU_mask_reflectivity.dat"
emiss_bb: 0.98 # E-REP-MPIA-MET-1203
diam_is_in: 25.4 # [mm] diameter of integrating sphere entrance port
diam_is: 250 # [mm] diameter of integrating sphere
diam_is_out: 100 # [mm] diameter if integrating sphere output port
rho_is: 0.95 # [] reflectivity of integrating sphere
rho_tube: 0.95 # [] reflectivity of tube between BB and IS
emiss_mask: 1.00 # [] emissivity of focal-plane mask
To use your own configuration file, copy the default file into your working directory and instantiate the OpticalTrain as follows:
cmd = sim.UserCommands(use_instrument="METIS", set_modes=["wcu_img_lm"])
cmd["!WCU.config_file") = "my_config.yaml"
metis = sim.OpticalTrain(cmd)
Some shortcomings#
The gaps in the connection between the black-body source and the entrance port of the integrating sphere are not yet included (model is available in Roy van Boekel’s document)
The PSF defaults to the standard SCAO PSF included in the METIS instrument package. At the very least a PSF without atmosphere contribution should be available. Ideally the PSF should be tied to the pupil-plane mask (this also applies to sky observations).
The focal-plane grid mask design is preliminary. The one included here has been taken from E-SPE-UZK-MET-1015 (v2), taking the Figure 5-6 as being to scale. It appears that this is not the case, though.
Pupil imaging could maybe be included. This could be included as an additional focal-plane mask with the appropriate plate scale onto the detectors. There are some differences as to which optical elements are included in the optical train as compared to the normal modes, but these are presumably of secondary importance.
Sanity checks#
In the following we add a few sanity checks, mostly for regression testing. These do not add any more information for users.
# Check that background from integrating sphere and focal-plane mask (opaque parts) are consistent.
cmd = sim.UserCommands(use_instrument="METIS", set_modes=["wcu_img_lm"])
metis = sim.OpticalTrain(cmd)
metis['reference_pixel_mask'].include = False
wcu = metis['wcu_source']
wcu.set_lamp('off')
wcu.set_fpmask('open')
metis.observe()
flux_is = np.median(metis.image_planes[0].data)
wcu.set_fpmask('pinhole_lm')
metis.observe()
flux_pin = np.median(metis.image_planes[0].data)
print(f"Open: {flux_is:.1f}")
print(f"Mask: {flux_pin:.1f}")
assert np.allclose(flux_is/flux_pin, wcu.meta['rho_is']/wcu.meta['emiss_mask'], rtol=0.01)