Level 1 pipeline data output

This module takes an ETC output spectrum table (created with qmostetc.QMostObservation.expose()), adds the noise, and then converts the result into a format similar to the L1 pipeline output. Both per-arm and joined output files can be generated.

The L1 pipeline output is described in VIS-TNO-4MOST-47110-0143-0001.

Note that this is not a full simulation of the L1 pipeline. Specifically, it derives the instrument parameters (sensitivity, …) from the values given in the ETC and not from calibration files.

For the examples, we assume that we created an exposure result with these lines:

>>> import astropy.units as u
>>> from qmostetc import QMostObservatory, SEDTemplate

>>> # Prepare observation
>>> observatory = QMostObservatory('lrs')
>>> obs = observatory(45*u.deg, 1.3*u.arcsec, 'gray')

>>> # Create an exposure with a sample target
>>> target = SEDTemplate('Pickles_G0V')(15*u.ABmag, 'GAIA2r.G')
>>> obs.set_target(target, 'point')
>>> texp = 1200 * u.s
>>> tbl = obs.expose(texp)
class qmostetc.L1DXU(observatory, table, texp, binwidth=None, with_noise=True)

Simulate the exposure and L1 pipeline. The L1 DXU

Parameters:
tableastropy.table.QTable

Output table of QMostObservation.expose().

observatoryQMostObservatory

Corresponding observatory (LRS or HRS).

texpastropy.units.Quantity

Exposure time [s].

binwidthastropy.units.Quantity

New wavelength array binning [Å]. Default (None) is to use 0.25 Å for LRS and 0.05 Å for HRS.

with_noisebool

Whether to add random noise. Defaults to True, but can be set to False for debugging purposes.

Examples

Create both the by-arm and joined spectra for the example observation above and save them to disk:

>>> from qmostetc import L1DXU

>>> # Create the DXU output wrapper
>>> dxu = L1DXU(observatory, tbl, texp)

>>> # Write individual L1 files
>>> for arm_name in observatory.keys():
...    dxu.per_arm(arm_name).writeto(f'l1_obs_{arm_name}.fits', overwrite=True)

>>> # Write the joined file
>>> dxu.joined().writeto('l1_obs_joined.fits', overwrite=True)
wavelength(arm_name)

Return the arm specific wavelength array

This array is regularly binned, with the bin width given in the constructor, and covers the range of the spectrograph arm.

Parameters:
arm_namestr

Spectrograph arm name (‘blue’, ‘green’, or ‘red’)

Returns:
astropy.units.Quantity

Wavelength array [Å]

Examples

Return the wavelength array of the blue arm:

>>> dxu.wavelength('blue')
<Quantity [3700.25, 3700.5 , ..., 5539.25, 5539.5 ] Angstrom>
sensitivity(arm_name)

Return the (arm specific) sensitivity function

This function is created from

  • the ADC gain,

  • the spectrograph sensitivity,

  • the spectrograph binwidth,

  • the telescope throughput,

  • the telescope area

It converts the adu counts into the flux at the entry of the telescope.

The values of the array correspond to the wavelenghths that can be retrieved with wavelength().

Parameters:
arm_namestr

Spectrograph arm name (‘blue’, ‘green’, or ‘red’)

Returns:
astropy.units.Quantity

Sensitivity factors [erg/(cm² adu)]

Examples

Give the sensitivity array for the blue arm:

>>> dxu.sensitivity('blue')
<Quantity [3.79301664e-15, ..., 3.87516497e-15] erg / (adu cm2)>
per_arm_spectrum(arm_name)

Return the uncalibrated L1 table for an individual arm

This adds some noise to the input table and rebins to the fixed-binned wavelength array. It assumes a perfect sky subtraction, i.e. the the FLUX and FLUX_NOSS columns just differ by the simulated sky.

Note that different from the DXU description, this is a conventional table with scalar values in each column.

Parameters:
arm_namestr

Spectrograph arm name (‘blue’, ‘green’, or ‘red’)

Returns:
astropy.table.QTable
DXU like table with the following columns:
  • WAVE: Wavelength array [Å],

  • FLUX: Spectral flux without sky [adu]

  • ERR_IVAR: Inverse variance without sky [1/adu²]

  • QUAL: Quality mask (always 0)

  • FLUX_NOSS: Spectral flux including sky [adu]

  • ERR_IVAR_NOSS: Inverse variance including sky [1/adu²]

  • SENSFUNC: Sensitivity function [erg/(s cm² Å)/adu]

Examples

Show the blue arm spectrum of the simulated observation:

>>> dxu.per_arm_spectrum('blue')
<QTable length=7358>
  WAVE     FLUX   ERR_IVAR  ...          SENSFUNC
Angstrom   adu    1 / adu2  ... erg / (Angstrom adu s cm2)
float32  float32  float32   ...          float32
-------- ------- ---------- ... --------------------------
 3700.25   38.05 0.02028342 ...                1.26434e-17
 3700.50   33.53 0.02023901 ...                1.26069e-17
     ...     ...        ... ...                        ...
 5539.50   66.58 0.01129934 ...                1.29172e-17
per_arm(arm_name)

HDUList with a L1 DXU compatible representation for a single arm

This basically creates a table with per_arm_spectrum() and converts it into a L1 DXU conform FITS HDUList. Different from the original table, the HDUList is column oriented: There is only one table row, and all values of a column are squashed into an array that goes into this row.

Also, the table comes with VO conform metadata, f.e. with UCD and UTYP for each column.

Parameters:
arm_namestr

Spectrograph arm name (‘blue’, ‘green’ or ‘red’)

Returns:
astropy.io.fits.HDUList

(Almost) L1 DXU conform HDUList

Examples

Create the HDU list for the blue arm L1 output:

>>> hdulist = dxu.per_arm('blue')
>>> hdulist[0].header['PRODCATG']
'SCIENCE.SPECTRUM'
>>> hdulist[1].name
'PHASE3SPECTRUM'
>>> hdulist[1].header['T*2']
TTYPE2  = 'FLUX    '           / Spectral flux
TFORM2  = '7358E   '           / Data format of field
TUNIT2  = 'adu     '           / Physical unit of field
TDISP2  = 'F7.2    '           / Display format for field
TDIM2   = '(7358)  '           / Array dimension of field
TCOMM2  = 'Spectral flux'      / Verbal description of field
TUCD2   = 'phot.flux.density;meta.main;stat.uncalib' / UCD for field
TUTYP2  = 'Spectrum.Data.FluxAxis.Value' / UType for field
>>> QTable.read(hdulist[1])  
...
<QTable length=1>
   WAVE [7512]        FLUX [7512]     ...        SENSFUNC [7512]
     Angstrom             adu         ...   erg / (adu Angstrom cm2 s)
     float32            float32       ...            float32
------------------ ------------------ ... ------------------------------
3671.50 .. 5549.25   21.07 ..   58.84 ...   1.52402e-17 ..   1.07370e-17
joined_spectrum()

Return the calibrated L1 table with the joined data.

This converts the input (CCD) values (given in adu) into the energy flux at the telescope’s entrance.

The sensitivity functions from the separate arms are not included in the joined table, but can be retrieved with the retrieved with the sensitivity() method.

Note that different from the DXU description, this is a conventional table with scalar values in each column.

Returns:
astropy.table.QTable
Calibrated spectrum table with the following columns:
  • WAVE: Wavelength array [Å],

  • FLUX: Spectral flux without sky [/adu]

  • ERR_FLUX: Inverse variance without sky [erg/(s cm² Å)/adu]

  • QUAL: Quality mask (always 0)

  • FLUX_NOSS: Spectral flux including sky [erg/(s cm² Å)/adu]

  • ERR_FLUX_NOSS: Inverse variance including sky [erg/(s cm² Å)/adu]

  • SENSFUNC: Sensitivity function [erg/(s cm² Å)/adu]

Examples

Show the joined spectrum of the simulated observation:

>>> dxu.joined_spectrum()
<QTable length=23198>
  WAVE            FLUX          ...     ERR_FLUX_NOSS
Angstrom erg / (Angstrom s cm2) ... erg / (Angstrom s cm2)
float32         float32         ...        float32
-------- ---------------------- ... ----------------------
 3700.25            4.81067e-16 ...            8.87755e-17
 3700.50            4.22690e-16 ...            8.86166e-17
     ...                    ... ...                    ...
 9499.50            3.41822e-16 ...            1.74843e-17
joined()

Create a L1 DXU HDUList with the joined spectrum

This basically creates a table with joined_spectrum() and converts it into a L1 DXU conform FITS HDUList. Different from the original table, the HDUList is column oriented: There is only one table row, and all values of a column are squashed into an array that goes into this row.

Also, the table comes with VO conform metadata, f.e. with UCD and UTYP for each column.

Returns:
astropy.io.fits.HDUList

(Almost) L1 DXU conform HDUList

Examples

Create the HDU list for the joined spectrum:

>>> hdulist = dxu.joined()
>>> hdulist[1].header['T*2']
TTYPE2  = 'FLUX    '           / Spectral flux
TFORM2  = '23198E  '           / Data format of field
TUNIT2  = 'erg Angstrom-1 s-1 cm-2' / Physical unit of field
TDISP2  = 'E13.5   '           / Display format for field
TDIM2   = '(23198) '           / Array dimension of field
TCOMM2  = 'Spectral flux'      / Verbal description of field
TUCD2   = 'phot.flux.density;meta.main' / UCD for field
TUTYP2  = 'Spectrum.Data.FluxAxis.Value' / UType for field
>>> QTable.read(hdulist[1])  
<QTable length=1>
   WAVE [23411]    ...     ERR_FLUX_NOSS [23411]
     Angstrom      ...     erg / (Angstrom cm2 s)
     float32       ...            float32
------------------ ... ------------------------------
3671.50 .. 9524.00 ...   9.53169e-17 ..   1.75297e-17

Output header

The L1DXU class creates a subset of the header specified in VIS-TNO-4MOST-47110-0143-0001. Only those header that can be usefully set by the ETC are included. The user may however add or modify additional header if required:

dxu = L1DXU(observatory, tbl, texp)
hdulist = dxu.joined()
hdulist[0].header['MJD-OBS'] = 59100.125  # Additional header
hdulist[0].header['TRG_NAME'] = ...       #      -"-
...
hdulist.write('l1_joined.fits')

Here are the headers included in the output:

Primary header

  • PRODCATG: Data product category ('SCIENCE.SPECTRUM')

  • ORIGIN: Observatory or facility ('ESO-PARANAL')

  • TELESCOP: ESO telescope designation ('VISTA')

  • INSTRUME: Instrument name ('4MOST')

  • SPECSYS: Frame of reference for spectral coordinates ('HELIOCEN')

  • CONTNORM: TRUE if spectrum normalised to continuum (False)

  • TOT_FLUX: TRUE if flux data represent total source flux (False)

  • FLUXERR: Fractional uncertainty on flux scale [0-100%] (10)

  • PROCSOFT: ETC faking 4MOST QL1 pipeline ('ETC/2.6.1')

  • QMNODE: 4MOST pipeline that produced this data product ('ETC')

  • RELEASE: Data release of this data product ('')

  • DXUDOC: Data eXchange Unit document ('VIS-TNO-4MOST-47110-0143-0001')

  • DATETIME: DateTime stamp when file submitted to NXP ('')

  • EXPTIME: Total integration time per data element [s]

  • TEXPTIME: Total integration time per data element [s]

  • WAVELMIN: Min wavelength for scientific information

  • WAVELMAX: Max wavelength for scientific information

  • OBSTECH: Technique used during observation ('MOS')

  • FLUXCAL: Quality of flux calib: ABSOLUTE or UNCALIBRATED ('UNCALIBRATED')

  • LAMRMS: RMS of residuals of wavelength solution [nm] (0)

  • LAMNLIN: Number of arc lines used in wavelength solution (0)

  • SPEC_BIN: Average spectral coordinate bin size [nm] (0.025)

  • SPEC_ERR: Statistical error in spectral coordinate [nm] (0.0)

  • SPEC_SYE: Systematic error in spectral coordinate [nm] (0.0)

Table extension header

  • EXTNAME: extension name ('PHASE3SPECTRUM')

  • VOCLASS: VO data model name and version ('SPECTRUM VX.x')

  • VOPUB: Name of the publisher ('4MOST/internal')

  • TITLE: Title of the dataset ('OBJECT QL1 Fake')

  • APERTURE: Aperture angular size [deg] (0.00039)

  • SKYSUB: FLUX spectrum is sky subtracted (True)

  • TELLSUB: Telluric correction applied (False)

  • TELLCOR: Pointer to telluric correction file used ('')

  • SENSFILE: Pointer to master sensfunc file that was scaled ('')

  • DELMAGZP: Mag offset applied as sensfunc scale factor (0.0)

  • NCALIBRT: Number of calibrators used (0)

  • HELIOCR: Applied heliocentric correction [km/s] (0.0)

  • WAVCR: Applied wavelength correction [Angstrom] (0.0)

  • WAVCRRMS: RMS of wavelength offset correction [Angstrom] (0.0)

  • SKYLNOFF: Offset: sky lines - rest wavelengths [Angstrom] (0.0)

  • SKYLNRMS: RMS of SKYLINE_OFF [Angstrom] (0.0)

  • SKYSHIFT: Offset: <sky> to spec in bckgd cor [Angstrom] (0.0)

  • SKYSCALE: Scaling: to <sky> before spec sky subtract (0)

  • TELAPSE: Total elapsed time [s]

  • EXPTIME: Total integration time per data element [s]

  • SPEC_VAL: Characteristic spectral coordinate value [nm]

  • SPEC_BW: Width of the spectrum [nm]

  • TOT_FLX: Total flux of extracted object

  • FLX_PIX: Mean flux per pixel of extracted object

Additionally to these headers, the metadata of all columns are included.