pypeit.wavemodel module

Module to create models of arc lines.

pypeit.wavemodel.addlines2spec(wavelength, wl_line, fl_line, resolution, scale_spec=1.0, debug=False)[source]

Create a spectrum with a set of (gaussian) emission lines.

Parameters:
  • wavelength (numpy.ndarray) – wavelength vector of the input spectrum

  • wl_line (numpy.ndarray) – wavelength of each individual line

  • fl_line (numpy.ndarray) – flux of each individual line

  • resolution (float) – resolution of the spectrograph. In other words, the lines will have a FWHM equal to: fwhm_line = wl_line / resolution

  • scale_spec (float) – rescale all the normalization of the final spectrum. Default scale_spec=1.

  • debug (bool) – If True will show debug plots

Returns:

line_spec – Spectrum with lines

Return type:

numpy.ndarray

pypeit.wavemodel.blackbody(wavelength, T_BB=250.0, debug=False)[source]

Given wavelength [in microns] and Temperature in Kelvin it returns the black body emission.

Parameters:
  • wavelength (numpy.ndarray) – wavelength vector in microns

  • T_BB (float) – black body temperature in Kelvin. Default is set to: T_BB = 250.

Returns:

  • blackbody (numpy.ndarray) – spectral radiance of the black body in cgs units: B_lambda = 2.*h*c^2/lambda^5.*(1./(exp(h*c/(lambda*k_b*T_BB))-1.)

  • blackbody_counts (numpy.ndarray) – Same as above but in flux density

pypeit.wavemodel.conv2res(wavelength, flux, resolution, central_wl='midpt', debug=False)[source]

Convolve an imput spectrum to a specific resolution. This is only approximate. It takes a fix FWHM for the entire spectrum given by: fwhm = wl_cent / resolution

Parameters:
  • wavelength (numpy.ndarray) – wavelength

  • flux (numpy.ndarray) – flux

  • resolution (float) – resolution of the spectrograph

  • central_wl – if ‘midpt’ the central pixel of wavelength is used, otherwise the central_wl will be used.

  • debug (bool) – If True will show debug plots

Returns:

  • flux_convolved (numpy.ndarray) – Resulting flux after convolution

  • px_sigma (float) – Size of the sigma in pixels at central_wl

  • px_bin (float) – Size of one pixel at central_wl

pypeit.wavemodel.create_OHlinelist(resolution, waveminmax=(0.8, 2.6), dlam=40.0, flgd=True, nirsky_outfile=None, fwhm=None, sigdetec=3.0, line_name='OH', file_root_name=None, iraf_frmt=False, debug=False)[source]

Create a synthetic sky spectrum at a given resolution, extract significant lines, and store them in a PypeIt compatibile file. The skymodel is built from nearIR_modelsky and includes black body at 250K, OH lines, and H2O lines (but only at lambda>2.3microns).

Parameters:
  • resolution (float) – resolution of the spectrograph

  • waveminmax (tuple) – wavelength range in microns to be covered by the model. Default is: (0.8, 2.6)

  • dlam – bin to be used to create the wavelength grid of the model. If flgd=’True’ it is a bin in velocity (km/s). If flgd=’False’ it is a bin in linear space (microns). Default is: 40.0 (with flgd=’True’)

  • flgd (bool) – if flgd=’True’ (default) wavelengths are created with equal steps in log space. If ‘False’, wavelengths will be created wit equal steps in linear space.

  • nirsky_outfile (str) – name of the fits file where the model sky spectrum will be stored. default is ‘None’ (i.e., no file will be written).

  • fwhm (float) – fwhm in pixels used for filtering out arc lines that are too wide and not considered in fits. Parameter of arc.detect_lines(). If set to ‘None’ the fwhm will be derived from the resolution as: 2. * central_wavelength / resolution

  • sigdetec (float) – sigma threshold above fluctuations for line detection. Parameter of arc.detect_lines(). Default = 2.

  • line_name (str) – name of the lines to listed in the file. Default is ‘OH’.

  • file_root_name (str) – name of the file where the identified lines will be stored. The code automatically add ‘_lines.dat’ at the end of the root name.

  • iraf_frmt (bool) – if True, the file is written in the IRAF format (i.e. wavelength, ion name, amplitude).

  • debug (bool) – If True will show debug plots

pypeit.wavemodel.create_ThArlinelist(resolution, waveminmax=(3000.0, 10500.0), dlam=40.0, flgd=True, thar_outfile=None, fwhm=None, sigdetec=3.0, line_name='ThAr', file_root_name=None, iraf_frmt=False, debug=False, convert_air_to_vac=True)[source]

Create a syntetic ThAr spectrum at a given resolution, extract significant lines, and store them in a PypeIt compatibile file. This is based on the Murphy et al. ThAr spectrum. Detailed information are here: http://astronomy.swin.edu.au/~mmurphy/thar/index.html

Parameters:
  • resolution (float) – resolution of the spectrograph

  • waveminmax (tuple) – wavelength range in ang. to be covered by the model. Default is: (3000., 10500.)

  • dlam – bin to be used to create the wavelength grid of the model. If flgd=’True’ it is a bin in velocity (km/s). If flgd=’False’ it is a bin in linear space (angstrom). Default is: 40.0 (with flgd=’True’)

  • flgd (bool) – if flgd=’True’ (default) wavelengths are created with equal steps in log space. If ‘False’, wavelengths will be created wit equal steps in linear space.

  • thar_outfile (str) – name of the fits file where the model sky spectrum will be stored. default is ‘None’ (i.e., no file will be written).

  • fwhm (float) – fwhm in pixels used for filtering out arc lines that are too wide and not considered in fits. Parameter of arc.detect_lines(). If set to ‘None’ the fwhm will be derived from the resolution as: 2. * central_wavelength / resolution

  • sigdetec (float) – sigma threshold above fluctuations for line detection. Parameter of arc.detect_lines(). Default = 2.

  • line_name (str) – name of the lines to listed in the file.

  • file_root_name (str) – name of the file where the identified lines will be stored. The code automatically add ‘_lines.dat’ at the end of the root name.

  • iraf_frmt (bool) – if True, the file is written in the IRAF format (i.e. wavelength, ion name, amplitude).

  • debug (bool) – If True will show debug plots

  • convert_air_to_vac (bool) – If True, convert the wavelengths of the created linelist from air to vacuum

pypeit.wavemodel.create_linelist(wavelength, spec, fwhm, sigdetec=2.0, cont_samp=10.0, line_name=None, file_root_name=None, iraf_frmt=False, debug=False, convert_air_to_vac=True)[source]

Create list of lines detected in a spectrum in a PypeIt compatible format. The name of the output file is file_root_name+’_lines.dat’.

Parameters:
  • wavelength (numpy.ndarray) – wavelength

  • spec (numpy.ndarray) – spectrum

  • fwhm (float) – fwhm in pixels used for filtering out arc lines that are too wide and not considered in fits. Parameter of arc.detect_lines().

  • sigdetec (float) – sigma threshold above fluctuations for line detection. Parameter of arc.detect_lines(). Default = 2.

  • cont_samp (float) – the number of samples across the spectrum used for continuum subtraction. Parameter of arc.detect_lines(). Default = 10.

  • line_name (str) – name of the lines to listed in the file.

  • file_root_name (str) – name of the file where the identified lines will be stored. The code automatically add ‘_lines.dat’ at the end of the root name.

  • iraf_frmt (bool) – if True, the file is written in the IRAF format (i.e. wavelength, ion name, amplitude).

  • convert_air_to_vac (bool) – If True, convert the wavelengths of the created linelist from air to vacuum

pypeit.wavemodel.h2o_lines()[source]

Reads in the H2O atmospheric spectrum”

Returns:

  • wavelength (numpy.ndarray) – Wavelength [in microns] of the H2O atmospheric spectrum.

  • flux (numpy.ndarray) – Flux of the H2O atmospheric spectrum.

pypeit.wavemodel.iraf_datareader(database_dir, id_file)[source]

Reads in a line identification database created with IRAF identify. These are usually locate in a directory called ‘database’. This read pixel location and wavelength of the lines that have been id with IRAF. Note that the first pixel in IRAF is ‘1’, while is ‘0’ in python. The pixel location is thus shifted of one pixel while reading the database.

Parameters:
  • database_dir (str) – directory where the id files are located.

  • id_file (str) – filename that is going to be read.

Returns:

  • pixel (numpy.ndarray) – Position of the line in pixel of the line. For IRAF output, these are usually in Ang.

  • line_id (numpy.ndarray) – ID of the line.

pypeit.wavemodel.nearIR_modelsky(resolution, waveminmax=(0.8, 2.6), dlam=40.0, flgd=True, nirsky_outfile=None, T_BB=250.0, SCL_BB=1.0, SCL_OH=1.0, SCL_H2O=10.0, WAVE_WATER=2.3, debug=False)[source]

Generate a model sky in the near-IR. This includes a continuum model to match to gemini broadband level, a black body at T_BB, OH lines, and H2O lines (but only at lambda>WAVE_WATER). Everythins is smoothed at the given resolution.

Parameters:
  • resolution (float) – resolution of the spectrograph. The OH and H2O lines will have a FWHM equal to: fwhm_line = wl_line / resolution

  • waveminmax (tuple) – wavelength range in microns to be covered by the model. Default is: (0.8, 2.6)

  • dlam – bin to be used to create the wavelength grid of the model. If flgd=’True’ it is a bin in velocity (km/s). If flgd=’False’ it is a bin in linear space (microns). Default is: 40.0 (with flgd=’True’)

  • flgd (bool) – if flgd=’True’ (default) wavelengths are created with equal steps in log space. If ‘False’, wavelengths will be created wit equal steps in linear space.

  • nirsky_outfile (str) – name of the fits file where the model sky spectrum will be stored. default is ‘None’ (i.e., no file will be written).

  • T_BB (float) – black body temperature in Kelvin. Default is set to: T_BB = 250.

  • SCL_BB (float) – scale factor for modelling the sky black body emssion. Default: SCL_BB=1.

  • SCL_OH (float) – scale factor for modelling the OH emssion. Default: SCL_OH=1.

  • SCL_H2O (float) – scale factor for modelling the H2O emssion. Default: SCL_H2O=10.

  • WAVE_WATER (float) – wavelength (in microns) at which the H2O are inclued. Default: WAVE_WATER = 2.3

  • debug (bool) – If True will show debug plots

Returns:

  • wave (numpy.ndarray) – Wavelength (in Ang.) of the final model of the sky.

  • sky_model (numpy.ndarray) – Flux of the final model of the sky.

pypeit.wavemodel.oh_lines()[source]

Reads in the Rousselot (2000) OH line list”

Returns:

  • wavelength (numpy.ndarray) – Wavelength [in microns] of the OH lines.

  • amplitude (numpy.ndarray) – Amplitude of the OH lines.

pypeit.wavemodel.optical_modelThAr(resolution, waveminmax=(3000.0, 10500.0), dlam=40.0, flgd=True, thar_outfile=None, debug=False)[source]

Generate a model of a ThAr lamp in the uvb/optical. This is based on the Murphy et al. ThAr spectrum. Detailed information are here: http://astronomy.swin.edu.au/~mmurphy/thar/index.html Everythins is smoothed at the given resolution.

Parameters:
  • resolution (float) – resolution of the spectrograph. The ThAr lines will have a FWHM equal to: fwhm_line = wl_line / resolution

  • waveminmax (tuple) – wavelength range in angstrom to be covered by the model. Default is: (3000.,10500.)

  • dlam – bin to be used to create the wavelength grid of the model. If flgd=’True’ it is a bin in velocity (km/s). If flgd=’False’ it is a bin in linear space (microns). Default is: 40.0 (with flgd=’True’)

  • flgd (bool) – if flgd=’True’ (default) wavelengths are created with equal steps in log space. If ‘False’, wavelengths will be created wit equal steps in linear space.

  • thar_outfile (str) – name of the fits file where the model sky spectrum will be stored. default is ‘None’ (i.e., no file will be written).

  • debug (bool) – If True will show debug plots

Returns:

  • wave (numpy.ndarray) – Wavelength (in Ang.) of the final model of the ThAr lamp emission.

  • thar_model (numpy.ndarray) – Flux of the final model of the ThAr lamp emission.

pypeit.wavemodel.thar_lines()[source]

Reads in the H2O atmospheric spectrum” Detailed information are here: http://astronomy.swin.edu.au/~mmurphy/thar/index.html

Returns:

  • wavelength (numpy.ndarray) – Wavelength [in microns] of the ThAr lamp spectrum.

  • flux (numpy.ndarray) – Flux of the ThAr lamp spectrum.

pypeit.wavemodel.transparency(wavelength, debug=False)[source]

Interpolate the atmospheric transmission model in the IR over a given wavelength (in microns) range.

Parameters:
  • wavelength (numpy.ndarray) – wavelength vector in microns

  • debug (bool) – If True will show debug plots

Returns:

transparency – Transmission of the sky over the considered wavelength range. 1. means fully transparent and 0. fully opaque

Return type:

numpy.ndarray