pypeit.core.telluric module

Module for telluric corrections.

class pypeit.core.telluric.Telluric(wave, flux, ivar, gpm, telgridfile, obj_params, init_obj_model, eval_obj_model, log10_blaze_function=None, ech_orders=None, sn_clip=30.0, teltype='pca', tell_npca=4, airmass_guess=1.5, resln_guess=None, resln_frac_bounds=(0.3, 1.5), pix_shift_bounds=(-5.0, 5.0), pix_stretch_bounds=(0.9, 1.1), maxiter=2, sticky=True, lower=3.0, upper=3.0, seed=777, ballsize=0.0005, tol=0.001, diff_evol_maxiter=1000, popsize=30, recombination=0.7, polish=True, disp=False, sensfunc=False, debug=False)[source]

Bases: DataContainer

Simultaneously fit model object and telluric spectra to an observed spectrum.

The model telluric absorption spectra for the atmosphere are generated by HITRAN atmosphere models with a baseline atmospheric profile model of the observatory in question. One can interpolate over grids of these absorption spectra directly, which have been created for each observatory, or instead (by default) use a PCA decomposition of the telluric models across all observatories, which is more flexible and much lighter in file size.

The modeling can be applied to both multi-slit or echelle data. The code performs a differential evolution optimization using scipy.optimize.differential_evolution. Several iterations are performed with rejection of outliers (governed by the maxiter parameter below). In general it produces very good results, but is rather slow. Progress can be monitored using the disp=True option, which will print out the progress to the screen. This will look like:

differential_evolution step 101: f(x)= 105392

Here f(x) represents the value of the loss function, which is essentially \(\chi^2\) (but with the non-Gaussian tails remapped and outliers rejected). Things are going well if the result of differential evolution yields a value of f(x) comparable to the number of pixels in your spectrum (or in the order in question for echelle data), which are the number of degrees of freedom. If debug=True, it will show the residuals from the fit at each iteration. Note that for spectra with S/N > 30, the residual distribution may look narrower because a floor is added to the noise (see sn_clip) to prevent excessive rejection and poorly behaved fits. Basically these telluric models are only good to about 3%, and so we cannot fit the data within the noise better than that, which is why sn_clip is implemented.

The object model is specified by the user through two user provided functions init_obj_model and eval_obj_model. The init_obj_model function is required because differential evolution can only perform bounded optimization. The bounds for the telluric are set by the grid and initialized without user involvement. But for the object model, the bounds depend on the nature of the object model, which is why init_obj_model must be provided.

The telluric model is governed by seven parameters — for the grid approach: pressure, temperature, humidity, airmass; for the pca approach 4 PCA coefficients; plus resln, shift, and stretch — where resln is the resolution of the spectrograph and shift is a shift in pixels. The stretch is a stretch of the pixel scale. The airmass of the object will be used to initalize the fit (this helps with initalizing the object model), but the models are sufficiently flexible that often the best-fit airmass actually differs from the airmass of the spectrum. In the pca approach, the number of PCA components used can be adjusted between 1 and 10, with a tradeoff between speed and accuracy.

This code can be run on stacked spectra covering a large range of airmasses and will still provide good results. The resulting airmass (in the grid approach) will be an effective value, and as per above may not have much relation to the true airmass of the observation. The exception to this rule is extremely high signal-to-noise ratio data (S/N > 30), for which it can be difficult to obtain a good fit within the noise of the data. In such cases, the user should split up the data into chunks taken around the same airmass, fit those individually with this class, and then combine the resulting telluric corrected spectra. This will, in general, result in better fits, and will also average down the residuals from the telluric model fit in the final averaged spectrum.

The datamodel attributes are:

Version: 1.0.0

Attribute

Type

Array Type

Description

airmass

float

Airmass of the observation

delta_zqso

float

Allowed range for the QSO redshift about z_qso

exptime

float

Exposure time (s)

func

str

Polynomial function used

lbound_norm

float

Flux normalization lower bound

model

astropy.table.table.Table

Table with the best-fitting model data

npca

int

Number of PCA components

pca_file

str

Name of the QSO PCA file

polish

bool

Perform a final optimization to tweak the best solution; see scipy.optimize.differential_evolution.

popsize

int

A multiplier for setting the total population size for the differential evolution optimization.

recombination

float

The recombination constant for the differential evolution optimization. Should be in the range [0, 1].

std_cal

str

File name (or shorthand) with the standard flux data

std_dec

float

DEC of the standard source

std_name

str

Type of standard source

std_ra

float

RA of the standard source

std_src

str

Name of the standard source

telgrid

str

File containing PCA components or grid of HITRAN atmosphere models

tell_norm_thresh

float

??

tell_npca

int

Number of telluric PCA components used

teltype

str

Type of telluric model, pca or grid

tol

float

Relative tolerance for converage of the differential evolution optimization.

ubound_norm

float

Flux normalization upper bound

z_qso

float

Redshift of the QSO

Todo

  • List the elements of obj_params.

Parameters:
  • wave (numpy.ndarray) – Wavelength array. Must either have shape (nspec,) or (nspec, norders), the latter being for echelle data.

  • flux (numpy.ndarray) – Flux for the object in question. Same shape as wave.

  • ivar (numpy.ndarray) – Inverse variance for the object in question. Same shape as wave.

  • gpm (numpy.ndarray) – Good pixel gpm for the object in question. Same shape as wave.

  • telgridfile (str) – File containing grid of HITRAN atmosphere models or PCA decomposition of such models, see TelluricPar.

  • teltype (str, optional) – Method for evaluating telluric models, either pca or grid. The grid method directly interpolates a pre-computed grid of HITRAN atmospheric models with nearest grid-point interpolation, while the pca method instead fits for the coefficients of a PCA decomposition of HITRAN atmospheric models, enabling a more flexible and far more file-size efficient interpolation of the telluric absorption model space.

  • tell_npca (int, optional) – Number of telluric PCA components used, must be <=10

  • obj_params (dict) – Dictionary of parameters for initializing the object model.

  • init_obj_model (callable) –

    User defined function for initializing the object model. This function must follow the calling sequence:

    obj_dict, bounds_obj = init_obj_model(obj_params, iord, wave, flux, ivar, gpm,
                                          tellmodel)
    

    See, e.g., init_star_model() for a detailed explanation of these paramaters and return values.

  • eval_obj_model (callable) –

    User defined function for evaluating the object model at a set of parameter values (theta_obj). This function must follow the calling sequence:

    obj_model, model_gpm = eval_obj_model(theta_obj, obj_dict)
    

    Where obj_dict is one of the return values from the init_obj_model above. See, e.g., eval_star_model() for a detailed explanation of these paramaters and return values.

  • log10_blaze_function (numpy.ndarray, optional) – The log10 blaze function determined from a flat field image. If this is passed in the sensitivity function model will be a (parametric) polynomial fit multiplied into the (non-parametric) log10_blaze_function. Shape = (nspec,) or (nspec, norddet), i.e. the same as wave.

  • ech_orders (numpy.ndarray, optional) – If passed the echelle orders will be included in the output data. Must be a numpy array of integers with the shape (norders,) giving the order numbers.

  • sn_clip (float, optional) – This adds an error floor to the ivar, preventing too much rejection at high-S/N (i.e. standard stars, bright objects) using the function clip_ivar(). A small error is added to the input ivar so that the output ivar_out will never give S/N greater than sn_clip. This prevents overly aggressive rejection in high S/N ratio spectra that neverthless differ at a level greater than the formal S/N due to the fact that our telluric models are only good to about 3%.

  • airmass_guess (float, optional) – A guess for the airmass of your object. The code fits the airmass as part of the telluric model, but this initial guess is useful for initializing the object model to determine the bounds for the object model parameter optimization via init_obj_model, since typically that is done by dividing out a guess for the telluric absorption and then performing some kind of initial fit.

  • resln_guess (float, optional) – A guess for the resolution of your spectrum expressed as lambda/dlambda. The resolution is fit explicitly as part of the telluric model fitting, but this guess helps determine the bounds for the optimization (see resln_frac_bounds). If not provided, the wavelength sampling of your spectrum will be used and the resolution calculated using a typical sampling of 3 spectral pixels per resolution element.

  • resln_frac_bounds (tuple, optional) – A two-tuple with the bounds for the resolution fit optimization which is part of the telluric model. This range is in units of resln_guess. For example, (0.5, 1.5) would bound the spectral resolution fit to be within the range (0.5*resln_guess, 1.5*resln_guess).

  • pix_shift_bounds (tuple, optional) – A two-tuple with the bounds for the pixel shift optimization in telluric model fit in units of pixels. The atmosphere will be allowed to shift within this range during the fit.

  • pix_stretch_bounds (tuple, optional) – A two-tuple with the bounds for the pixel scale stretch optimization in telluric model fit in units of pixels. The atmosphere spectrum will be allowed to stretch within this range during the fit.

  • maxiter (int, optional) – Maximum number of iterations for the telluric + object model fitting. The code performs multiple iterations rejecting outliers at each step. The fit is then performed anew to the remaining good pixels. For this reason if you run with the disp=True option, you will see that the f(x) loss function gets progressively better during the iterations.

  • sticky (bool, optional) – Sticky parameter for the djs_reject() algorithm for iterative model fit rejection. If True then points rejected from a previous iteration are kept rejected, in other words the bad pixel mask is the OR of all previous iterations and rejected pixels accumulate. If False, the bad pixel mask is the mask from the previous iteration, and if the model fit changes between iterations, points can alternate from being rejected to not rejected. At present this code only performs optimizations with differential evolution and experience shows that sticky needs to be True in order for these to converge. This is because the outliers can be so large that they dominate the loss function, and one never iteratively converges to a good model fit. In other words, the deformations in the model between iterations with sticky=False are too small to approach a reasonable fit.

  • lower (float, optional) – Lower rejection threshold in units of sigma_corr*sigma, where sigma is the formal noise of the spectrum and sigma_corr is an empirically determined correction to the formal error. The distribution of input chi (defined by chi = (data - model)/sigma) values is analyzed, and a correction factor to the formal error sigma_corr is returned, which is multiplied into the formal errors. In this way, a rejection threshold of i.e. 3-sigma, will always correspond to roughly the same percentile. This renormalization is performed with renormalize_errors(), and guarantees that rejection is not too agressive in cases where the empirical errors determined from the chi-distribution differ significantly from the formal noise, which is used to determine chi.

  • upper (float, optional) – Upper rejection threshold in units of sigma_corr*sigma, where sigma is the formal noise of the spectrum and sigma_corr is an empirically determined correction to the formal error. See lower for more detail.

  • seed (int, optional) – An initial seed for the differential evolution optimization, which is a random process. The default is a seed = 777 which will be used to seed a numpy.random.Generator object. A specific seed is used because otherwise the random number generator will use the time for the seed and the results will not be reproducible. TODO: (JFH) Note that differential evolution seems to have some other source of stochasticity that I have not yet figured out.

  • ballsize (float, optional) – This parameter governs how the differential evolution random population is initialized for the object model and for subsequent iterations. The idea is that after the first iteration we can generate a new random population in model space about the optimum from the previous iteration. This significantly speeds up convergence. The new population of parameters is created by perturbing the previous optimum by theta = theta_opt + ballsize*(bounds[1]-bounds[0])*standard_normal_deviate, where bounds are the bounds for the optimization and standard_normal_deviate is a random draw from a unit variance Gaussian. If the ballsize is too small, the optimization may get stuck in a local maximum. If it is too large it will not speed up convergence. This choice is highly dependent on how the bounds are chosen using in the init_obj_model function and from the bound parameters for the atmospheric, shift, and resolution parameters. If the init_obj_model functions also have a way of isolating a good guess for the object model, then this can be put assigned to obj_dict['init_obj_opt_theta'] and the object model will similarly be initialized about this optimum using the ballsize for the first iteration. See init_sensfunc_model() for an example.

  • tol (float, optional) – Relative tolerance for converage of the differential evolution optimization. See scipy.optimize.differential_evolution for details.

  • popsize (int, optional) – A multiplier for setting the total population size for the differential evolution optimization. See scipy.optimize.differential_evolution for details.

  • recombination (float, optional) – The recombination constant for the differential evolution optimization. This should be in the range [0, 1]. See scipy.optimize.differential_evolution for details.

  • polish (bool, optional) – If True then differential evolution will perform an additional optimization at the end to polish the best fit at the end, which can improve the optimization slightly. See scipy.optimize.differential_evolution for details.

  • disp (bool, optional) – Argument for scipy.optimize.differential_evolution that prints status messages to the screen indicating the status of the optimization. See description above.

  • sensfunc (bool, optional) – This option is used for usage of this class for joint telluric fitting and sensitivity function computation. If True, the input flux is in counts and then converted to counts per angstrom, since the sensfunc is obtained by fitting counts per angstrom. TODO: This explanation is unclear.

  • debug (bool, optional) – If True, QA plots will be shown to the screen indicating the quality of the fits. Specifically, the residual distributions will be shown at each iteration, and the fit will be shown at the end (for each order). This is useful if you are running the code for the first time, but since the algorithm is slow, particularly for fitting multi-order echelle data, it will require lots of clicking to close interactive matplotlib windows which halt code flow.

assign_output(iord)[source]

Routine to assign outputs to self.model for the order in question.

Parameters:

iord (int) – The order for which the output table should bbe assigned.

datamodel = {'airmass': {'descr': 'Airmass of the observation', 'otype': <class 'float'>}, 'delta_zqso': {'descr': 'Allowed range for the QSO redshift about z_qso', 'otype': <class 'float'>}, 'exptime': {'descr': 'Exposure time (s)', 'otype': <class 'float'>}, 'func': {'descr': 'Polynomial function used', 'otype': <class 'str'>}, 'lbound_norm': {'descr': 'Flux normalization lower bound', 'otype': <class 'float'>}, 'model': {'descr': 'Table with the best-fitting model data', 'otype': <class 'astropy.table.table.Table'>}, 'npca': {'descr': 'Number of PCA components', 'otype': <class 'int'>}, 'pca_file': {'descr': 'Name of the QSO PCA file', 'otype': <class 'str'>}, 'polish': {'descr': 'Perform a final optimization to tweak the best solution; see scipy.optimize.differential_evolution.', 'otype': <class 'bool'>}, 'popsize': {'descr': 'A multiplier for setting the total population size for the differential evolution optimization.', 'otype': <class 'int'>}, 'recombination': {'descr': 'The recombination constant for the differential evolution optimization. Should be in the range [0, 1].', 'otype': <class 'float'>}, 'std_cal': {'descr': 'File name (or shorthand) with the standard flux data', 'otype': <class 'str'>}, 'std_dec': {'descr': 'DEC of the standard source', 'otype': <class 'float'>}, 'std_name': {'descr': 'Type of standard source', 'otype': <class 'str'>}, 'std_ra': {'descr': 'RA of the standard source', 'otype': <class 'float'>}, 'std_src': {'descr': 'Name of the standard source', 'otype': <class 'str'>}, 'telgrid': {'descr': 'File containing PCA components or grid of HITRAN atmosphere models', 'otype': <class 'str'>}, 'tell_norm_thresh': {'descr': '??', 'otype': <class 'float'>}, 'tell_npca': {'descr': 'Number of telluric PCA components used', 'otype': <class 'int'>}, 'teltype': {'descr': 'Type of telluric model, `pca` or `grid`', 'otype': <class 'str'>}, 'tol': {'descr': 'Relative tolerance for converage of the differential evolution optimization.', 'otype': <class 'float'>}, 'ubound_norm': {'descr': 'Flux normalization upper bound', 'otype': <class 'float'>}, 'z_qso': {'descr': 'Redshift of the QSO', 'otype': <class 'float'>}}

DataContainer datamodel.

static empty_model_table(norders, nspec, tell_npca=4, n_obj_par=0)[source]

Construct an empty astropy.table.Table for the telluric model results.

Parameters:
  • norders (int) – The number of slits/orders on the detector.

  • nspec (int) – The number of spectral pixels on the detector.

  • tell_npca (int) – Number of telluric model parameters

  • n_obj_par (int, optional) – The number of parameters used to construct the object model spectrum.

Returns:

Instance of the empty model table.

Return type:

astropy.table.Table

get_bounds_tell()[source]

Return the parameter boundaries for the telluric model.

Returns:

A list of two-tuples, where each two-tuple proives the minimum and maximum allowed model values for the PCA coefficients, (or for the grid approach: pressure, temperature, humidity, airmass) as well as the resolution, shift, and stretch parameters.

Return type:

list

get_ind_lower_upper()[source]

Utiltity routine to determine the ind_lower and ind_upper for each order. This trimming makes things faster because then we only need to convolve the portion of the telluric model that is needed for the model fit to each order, rather than convolving the entire telluric model grid.

Returns:

ind_lower, ind_upper

ind_lower (int):

Lower index into the telluric model wave_grid to trim down the telluric model.

ind_upper (int):

Upper index into the telluric model wave_grid to trim down the telluric model.

get_tell_guess()[source]

Return guess parameters for the telluric model.

This is also used to determine the bounds with the init_obj_model function.

Returns:

The guess telluric model parameters, resolution, shift, and stretch parameters.

Return type:

tuple

init_output()[source]

Method to initialize the outputs

internals = ['obj_params', 'init_obj_model', 'airmass_guess', 'eval_obj_model', 'ech_orders', 'sn_clip', 'resln_frac_bounds', 'pix_shift_bounds', 'pix_stretch_bounds', 'maxiter', 'sticky', 'lower', 'upper', 'seed', 'rng', 'ballsize', 'diff_evol_maxiter', 'disp', 'sensfunc', 'debug', 'wave_in_arr', 'flux_in_arr', 'ivar_in_arr', 'mask_in_arr', 'log10_blaze_func_in_arr', 'nspec_in', 'norders', 'tell_dict', 'wave_grid', 'ngrid', 'resln_guess', 'tell_guess', 'bounds_tell', 'flux_arr', 'ivar_arr', 'mask_arr', 'log10_blaze_func_arr', 'wave_mask_arr', 'ind_lower', 'ind_upper', 'srt_order_tell', 'obj_dict_list', 'bounds_obj_list', 'bounds_list', 'arg_dict_list', 'max_ntheta_obj', 'result_list', 'outmask_list', 'obj_model_list', 'tellmodel_list', 'theta_obj_list', 'theta_tell_list']

A list of strings identifying a set of internal attributes that are not part of the datamodel.

interpolate_inmask(mask, wave_inmask, inmask)[source]

Utitlity routine to interpolate the input mask.

run(only_orders=None)[source]

Loops over orders/slits, runs the telluric correction, and evaluates the object and telluric models.

Parameters:

only_orders (list, int, optional) – A single integer or list of integers setting the orders to correct. If None, all orders are corrected.

show_fit_qa(iord)[source]

Generates QA plot for telluric fitting

Parameters:

iord – the order being currently fit

sort_telluric()[source]

Sort the model telluric spectra by the strength of their telluric features.

This is used to set the order in which the telluric model is fit for multiple spectra (either multiple orders or multiple slits). The strength of the telluric features is determined by computing the median telluric absorption at the midpoint of parameters governing the telluric grid, specifically for the wavelengths included when fitting each spectrum.

Returns:

Array of sorted indices from strongest telluric absorption to weakest.

Return type:

numpy.ndarray

version = '1.0.0'

Datamodel version.

pypeit.core.telluric.conv_telluric(tell_model, dloglam, res)[source]

Routine to convolve the telluric model to desired resolution.

Parameters:
  • tell_model (numpy.ndarray) – Input telluric model at the native resolution of the telluric model grid. The shape of this input is in general different from the size of the raw telluric model (read in by read_telluric_* above) because it is trimmed to relevant wavelengths using ind_lower, ind_upper. See eval_telluric below.

  • dloglam (float) – Wavelength spacing of the telluric grid expressed as a dlog10(lambda), i.e. stored in the tell_dict as tell_dict[‘dloglam’]

  • res (float) – Desired resolution expressed as lambda/dlambda. Note that here dlambda is linear, whereas dloglam is the delta of the log10.

Returns:

Resolution convolved telluric model. Shape = same size as input tell_model.

Return type:

numpy.ndarray

pypeit.core.telluric.create_bal_mask(wave, bal_wv_min_max)[source]

Example of a utility function for creating a BAL mask for QSOs with BAL features. Can also be used to mask other features that the user does not want to fit.

Parameters:

wave (numpy.ndarray) – Wavelength array for the quasar in question shape=(nspec,)

Returns:

gpm – Good pixel (non-bal pixels) mask for the fits. shape=(nspec,)

Return type:

numpy.ndarray

pypeit.core.telluric.eval_poly_model(theta, obj_dict)[source]

Routine to evaluate a star spectrum model as a true model spectrum times a polynomial.

Parameters:
  • theta (numpy.ndarray) – Array containing the polynomial coefficients. shape=(ntheta,)

  • obj_dict (dict) – Dictionary containing additional arguments needed to evaluate the star model.

Returns:

  • star_model (numpy.ndarray) – array with same shape obj_dict[‘polymodel’]

  • gpm (numpy.ndarray) – Good pixel mask indicating where the model is valid. array with same shape as the star_model.

pypeit.core.telluric.eval_qso_model(theta, obj_dict)[source]

Routine to evaluate a sensitivity function model for a QSO spectrum.

Parameters:
  • theta (numpy.ndarray) – Array containing the PCA coefficients shape=(ntheta,)

  • obj_dict (dict) – Dictionary containing additional arguments needed to evaluate the PCA model

Returns:

  • qso_pca_model (array with same shape as the PCA vectors (stored in the obj_dict[‘pca_dict’])) – PCA vectors were already interpolated onto the telluric model grid by init_qso_model

  • gpm (numpy.ndarray : array with same shape as the qso_pca_model) – Good pixel mask indicating where the model is valid

pypeit.core.telluric.eval_sensfunc_model(theta, obj_dict)[source]

Utility routine to evaluate a sensitivity function model

Parameters:
  • theta (numpy.ndarray) – Array containing the parameters that describe the model. shape is (ntheta,)

  • obj_dict (dict) – Dictionary containing additional arguments needed to evaluate the model

Returns:

  • counts_per_angstrom_model (numpy.ndarray, shape is the same as obj_dict[‘wave_star’]) – Model of the quantity being fit for the sensitivity function which is exptime*f_lam_true/sensfunc. Note that this routine is used to fit the counts per angstrom. The zeropoint is defined to be zeropoint 2.5log10(S_lam) where S_lam = F_lam/(counts/s/A) where F_lam is in units of 1e-17 erg/s/cm^2/A, and so S_lam has units of (erg/s/cm^2/A)/(counts/s/A) = erg/cm^2/counts, and the zeropoint is then the magnitude of the astronomical source that produces one count per second in a one-angstrom wide filter on the detector. And zerpoint + 2.5log10(delta_lambda/1 A) is the magnitude that produces one count per spectral pixel of width delta_lambda

  • gpm (numpy.ndarray, bool, shape is the same as obj_dict[‘wave_star’]) – Good pixel mask indicating where the model is valid

pypeit.core.telluric.eval_star_model(theta, obj_dict)[source]

Routine to evaluate a star spectrum model as a true model spectrum times a polynomial.

Parameters:
  • theta (numpy.ndarray) – Array containing the polynomial coefficients. shape (ntheta,)

  • obj_dict (dict) – Dictionary containing additional arguments needed to evaluate the star model.

Returns:

  • star_model (numpy.ndarray) – PCA vectors were already interpolated onto the telluric model grid by init_qso_model. array with same shape obj_dict[‘wave’]

  • gpm (numpy.ndarray) – Good pixel mask indicating where the model is valid. array with same shape as the star_model

pypeit.core.telluric.eval_telluric(theta_tell, tell_dict, ind_lower=None, ind_upper=None)[source]

Evaluate the telluric model.

Parameters from theta_tell are: (if teltype == ‘pca’) the PCA coefficients or (if teltype == ‘grid’) the telluric grid parameters pressure, temperature, humidity, and airmass, in both cases followed by spectral resolution, shift, and stretch.

This routine performs the following steps:

  1. summation of telluric PCA components multiplied by coefficients

  2. transformation of telluric PCA model from arsinh(tau) to transmission

  3. convolution of the atmosphere model to the resolution set by the spectral resolution.

  4. shift and stretch the telluric model.

Parameters:
  • theta_tell (numpy.ndarray) – Vector with tell_npca PCA coefficients (if teltype=’pca’) or pressure, temperature, humidity, and airmass (if teltype=’grid’), followed by spectral resolution, shift, and stretch. Final length is then tell_npca+3 or 7.

  • tell_dict (dict) – Dictionary containing the telluric data. See read_telluric_pca() if teltype==’pca’. or read_telluric_grid() if teltype==’grid’.

  • ind_lower (int, optional) – The index of the first pixel to include in the model. Selecting a wavelength region for the modeling makes things faster because we only need to convolve the portion that is needed for the current model fit.

  • ind_upper (int, optional) – The index (inclusive) of the last pixel to include in the model. Selecting a wavelength region for the modeling makes things faster because we only need to convolve the portion that is needed for the current model fit.

Returns:

Telluric model evaluated at the desired location theta_tell in model atmosphere parameter space. Shape is given by the size of wave_grid plus tell_pad_pix padding from the input tell_dict.

Return type:

numpy.ndarray

pypeit.core.telluric.general_spec_reader(specfile, ret_flam=False, chk_version=True)[source]

Read a spec1d file or a coadd spectrum file.

Parameters:
  • specfile (str) – File with the data

  • ret_flam (bool, optional) – Return FLAM instead of COUNTS.

  • chk_version (bool, optional) – When reading in existing files written by PypeIt, perform strict version checking to ensure a valid file. If False, the code will try to keep going, but this may lead to faults and quiet failures. User beware!

Returns:

Seven objects are returned. The first five are numpy.ndarray objects with the wavelength, bin centers of wavelength grid, flux, inverse variance, and good-pixel mask for each spectral pixel. The 6th is a dict of metadata. And the 7th is an astropy.io.fits.Header object with the primary header from the input file.

Return type:

tuple

pypeit.core.telluric.init_poly_model(obj_params, iord, wave, flux, ivar, mask, tellmodel)[source]

Routine used to initialize a polynomial object model for telluric fits.

Parameters:
  • obj_params (dict) – Dictionary containing parameters necessary for initializing the quasar model

  • iord (int) – Order in question. This is used here because each echelle order can have a different polynomial order

  • wave (array shape (nspec,)) – Wavelength array for the object in question

  • flux (array shape (nspec,)) – Flux array for the object in question

  • ivar (array shape (nspec,)) – Inverse variance array for the oejct in question

  • mask (array shape (nspec,)) – Good pixel mask for the object in question

  • tellmodel (array shape (nspec,)) – This is a telluric model computed on the wave wavelength grid. Initialization usually requires some initial best guess for the telluric absorption, which is computed from the midpoint of the telluric model grid parameter space using the resolution of the spectrograph and the airmass of the observations.

Returns:

  • obj_dict (dict) – Dictionary containing the meta-information and variables that are used for the object model evaluations.

  • bounds_obj (tuple) – Tuple of bounds for each parameter that will be fit for the object model, which are here the polynomial coefficients.

pypeit.core.telluric.init_qso_model(obj_params, iord, wave, flux, ivar, mask, tellmodel)[source]

Routine used to initialize the quasar spectrum telluric fits

Parameters:
  • obj_params (dict) – Dictionary containing parameters necessary for initializing the quasar model

  • iord (int) – Order in question. This is not used for initializing the qso model, but is kept here for compatibility with the standard function argument list

  • wave (array shape (nspec,)) – Wavelength array for the object in question

  • flux (array shape (nspec,)) – Flux array for the object in question

  • ivar (array shape (nspec,)) – Inverse variance array for the oejct in question

  • mask (array shape (nspec,)) – Good pixel mask for the object in question

  • tellmodel (array shape (nspec,)) – This is a telluric model computed on the wave wavelength grid. Initialization usually requires some initial best guess for the telluric absorption, which is computed from the mean of the telluric model grid using the resolution of the spectrograph.

Returns:

  • obj_dict (dict) – Dictionary containing the meta-information and variables that are used for the object model evaluations. For example for the quasar model which based on a PCA decomposition, this dictionary holds the number of PCA basis vectors and those basis vectors.

  • bounds_obj (tuple) – Tuple of bounds for each parameter that will be fit for the object model

pypeit.core.telluric.init_sensfunc_model(obj_params, iord, wave, counts_per_ang, ivar, gpm, tellmodel)[source]

Initializes a sensitivity function model fit for joint sensitivity function and telluric fitting by setting up the obj_dict and bounds.

Parameters:
  • obj_params (dict) – Dictionary of object parameters

  • iord (int) – The slit/order in question

  • wave (numpy.ndarray, float, shape is (nspec,)) – Wavelength array for the standard star

  • counts_per_ang (numpy.ndarray, float, shape is (nspec,)) – Counts per Angstrom array for the standard star

  • ivar (numpy.ndarray, float, shape is (nspec,)) – Inverse variance of the counts per Angstrom array for the standard star

  • gpm (numpy.ndarray, bool, shape is (nspec,)) – Good pixel mask for the counts per Angstrom array for the standard star

  • tellmodel (numpy.ndarray, float, shape is (nspec,)) – Telluric absorption model guess

Returns:

  • obj_dict (dict) – Dictionary of object paramaters for joint sensitivity function telluric model fitting

  • bounds_obj (list) – List of bounds for each parameter in the joint sensitivity function and telluric model fit.

pypeit.core.telluric.init_star_model(obj_params, iord, wave, flux, ivar, mask, tellmodel)[source]

Routine used to initialize the star spectrum model for telluric fits. The star model is the true spectrum of the standard star times a polynomial to accomodate situations where the star-model is not perfect.

Parameters:
  • obj_params (dict) – Dictionary containing parameters necessary for initializing the quasar model

  • iord (int) – Order in question. This is used here because each echelle order can have a different polynomial order

  • wave (array shape (nspec,)) – Wavelength array for the object in question

  • flux (array shape (nspec,)) – Flux array for the object in question

  • ivar (array shape (nspec,)) – Inverse variance array for the oejct in question

  • mask (array shape (nspec,)) – Good pixel mask for the object in question

  • tellmodel (array shape (nspec,)) – This is a telluric model computed on the wave wavelength grid. Initialization usually requires some initial best guess for the telluric absorption, which is computed from the midpoint of the telluric model grid parameter space using the resolution of the spectrograph and the airmass of the observations.

Returns:

  • obj_dict (dict) – Dictionary containing the meta-information and variables that are used for the object model evaluations.

  • bounds_obj (tuple) – Tuple of bounds for each parameter that will be fit for the object model, which are here the polynomial coefficients.

pypeit.core.telluric.interp_telluric_grid(theta, tell_dict)[source]

Interpolate the telluric model grid to the specified location in parameter space.

The interpolation is only performed over the 4D parameter space specified by pressure, temperature, humidity, and airmass. This routine performs nearest-gridpoint interpolation to evaluate the telluric model at an arbitrary location in this 4-d space. The telluric grid is assumed to be uniformly sampled in this parameter space.

Parameters:
  • theta (numpy.ndarray) – A 4-element vector with the telluric model parameters: pressure, temperature, humidity, and airmass.

  • tell_dict (dict) – Dictionary containing the telluric grid. See read_telluric_grid().

Returns:

Telluric model evaluated at the provided 4D position in parameter space. The telluric model is provided over all available wavelengths in tell_dict.

Return type:

numpy.ndarray

pypeit.core.telluric.poly_telluric(spec1dfile, telgridfile, telloutfile, outfile, z_obj=0.0, func='legendre', model='exp', polyorder=3, fit_wv_min_max=None, mask_lyman_a=True, teltype='pca', tell_npca=4, delta_coeff_bounds=(-20.0, 20.0), minmax_coeff_bounds=(-5.0, 5.0), only_orders=None, sn_clip=30.0, maxiter=3, tol=0.001, popsize=30, recombination=0.7, polish=True, disp=False, pix_shift_bounds=(-5.0, 5.0), debug_init=False, debug=False, show=False, chk_version=True)[source]

This needs a doc string.

FILL IN

COULD USE THE qso_telluric() doc string as a starting point

pypeit.core.telluric.qso_init_pca(filename, wave_grid, redshift, npca)[source]

This routine reads in the pickle file created by coarse_pca.create_coarse_pca. The relevant pieces are the wavelengths (wave_pca_c), the PCA components (pca_comp_c), and the Gaussian mixture model prior (mix_fit).

FILL THESE IN

Parameters:
  • filename (str) – Pickle filename

  • wave_grid (numpy.ndarray) –

  • redshift – (float): Approximate redshift of the quasar.

  • npca – Numer of pca components.

Return type:

dict

pypeit.core.telluric.qso_pca_eval(theta, qso_pca_dict)[source]

Function for evaluating the quasar PCA

Parameters:
  • theta (numpy.ndarray) – Parameter vector, where theta_pca[0] is redshift, theta_pca[1] is the normalization, and theta_pca[2:npca+1] are the PCA coefficients, where npca is the PCA dimensionality

  • qso_pca_dict (dict) – Dictionary continaing the PCA information generated by qso_init_pca

Returns:

Evaluated PCA for the QSO

Return type:

numpy.ndarray

pypeit.core.telluric.qso_telluric(spec1dfile, telgridfile, pca_file, z_qso, telloutfile, outfile, npca=8, pca_lower=1220.0, pca_upper=3100.0, bal_wv_min_max=None, delta_zqso=0.1, teltype='pca', tell_npca=4, bounds_norm=(0.1, 3.0), tell_norm_thresh=0.9, sn_clip=30.0, only_orders=None, maxiter=3, tol=0.001, popsize=30, recombination=0.7, polish=True, disp=False, pix_shift_bounds=(-5.0, 5.0), debug_init=False, debug=False, show=False, chk_version=True)[source]

Fit and correct a QSO spectrum for telluric absorption.

Parameters:
  • spec1dfile (str) – File name of the input 1D spectrum. Can be a PypeIt spec1d file or the output from a PypeIt 1D coadd.

  • telgridfile (str) – File name with the grid of telluric spectra.

  • pca_file (str) – File name of the QSO PCA model fits file.

  • z_qso (float) – QSO redshift.

  • telloutfile (str) – Output file name for the best-fit telluric model.

  • outfile (str) – Output file name of the telluric corrected spectrum.

  • npca (int, optional) – Numer of pca components.

  • pca_lower (float, optional) – Wavelength lower bounds of the PCA model

  • pca_upper (float, optional) – Wavelength upper bounds of the PCA model

  • bal_wv_min_max (list, optional) – Broad absorption line feature mask. If there are several BAL features, the format for this mask is [wave_min_bal1, wave_max_bal1, wave_min_bal2, wave_max_bal2,...]. These wavelength ranges will be ignored during the fitting.

  • delta_zqso (float, optional) – During the fit, the QSO redshift is allowed to vary within +/-delta_zqso.

  • teltype (str, optional, default = ‘pca’) – Method for evaluating telluric models, either pca or grid.

  • tell_npca (int, optional, default = 4) – Number of telluric PCA components used, must be <=10

  • bounds_norm (tuple, optional) – A two-tuple with the lower and upper bounds on the fractional adjustment of the flux in the QSO model spectrum. For example, a value of (0.1, 3.0) means the best-fit QSO model flux is allowed to vary between 10% and a factor of 3 times the input QSO model.

  • tell_norm_thresh (float, optional) – When initializing the model, pixels in the reference telluric spectrum below this value are not included in an initial guess of the model flux normalization.

  • sn_clip (float, optional, default=30.0) – Errors are capped during rejection so that the S/N is never greater than sn_clip. This prevents overly aggressive rejection in high S/N ratio spectrum, which neverthless differ at a level greater than the implied S/N due to systematics. These systematics are most likely due to the inadequecy of our telluric model to fit the data at the the precision required at very high S/N ratio.

  • only_orders (list of int, optional, default = None) – Only fit these specific orders

  • maxiter (int, optional, default = 3) – Maximum number of iterations for the telluric + object model fitting. The code performs multiple iterations rejecting outliers at each step. The fit is then performed anew to the remaining good pixels. For this reason if you run with the disp=True option, you will see that the f(x) loss function gets progressively better during the iterations.

  • tol (float, default = 1e-3) – Relative tolerance for converage of the differential evolution optimization. See scipy.optimize.differential_evolution for details.

  • popsize (int, default = 30) – A multiplier for setting the total population size for the differential evolution optimization. See scipy.optimize.differential_evolution for details.

  • recombination (float, optional, default = 0.7) – The recombination constant for the differential evolution optimization. This should be in the range [0, 1]. See scipy.optimize.differential_evolution for details.

  • polish (bool, optional, default=True) – If True then differential evolution will perform an additional optimizatino at the end to polish the best fit at the end, which can improve the optimization slightly. See scipy.optimize.differential_evolution for details.

  • disp (bool, optional, default=True) – Argument for scipy.optimize.differential_evolution that will display status messages to the screen indicating the status of the optimization. See above for a description of the output and how to know if things are working well. Note: this is automatically set to True if debug is True.

  • debug_init (bool, optional, default=False) – Show plots to the screen useful for debugging model initialization

  • debug (bool, optional, default=False) – Show plots to the screen useful for debugging the telluric/object model fits.

  • show (bool, optional) – Show a QA plot of the final fit.

  • chk_version (bool, optional) – When reading in existing files written by PypeIt, perform strict version checking to ensure a valid file. If False, the code will try to keep going, but this may lead to faults and quiet failures. User beware!

Returns:

TelObj – Object with the telluric modeling results

Return type:

Telluric

pypeit.core.telluric.read_telluric_grid(filename, wave_min=None, wave_max=None, pad_frac=0.1)[source]

Reads in the telluric grid from a file. This method is no longer the preferred approach; see “read_telluric_pca” for the PCA mode.

Optionally, this method also trims the grid to be in within wave_min and wave_max and pads the data (see pad_frac).

Parameters:
  • filename (str) – Telluric grid filename

  • wave_min (float, optional) – Minimum wavelength at which the grid is desired

  • wave_max (float, optional) – Maximum wavelength at which the grid is desired.

  • pad_frac (float, optional) – Percentage padding to be added to the grid boundaries if wave_min or wave_max are input; ignored otherwise. The resulting grid will extend from (1.0 - pad_frac)*wave_min to (1.0 + pad_frac)*wave_max.

Returns:

Dictionary containing the telluric grid
  • wave_grid: Telluric model wavelength grid

  • dloglam: Wavelength sampling of telluric model

  • tell_pad_pix: Number of pixels to pad at edges for convolution

  • pressure_grid: Atmospheric pressure values in telluric grid [mbar]

  • temp_grid: Temperature values in telluric grid [degrees C]

  • h2o_grid: Humidity values in telluric grid [%]

  • airmass_grid: Airmass values in telluric grid

  • tell_grid: Grid of telluric models

  • teltype: Type of telluric model, i.e. ‘grid’

Return type:

dict

pypeit.core.telluric.read_telluric_pca(filename, wave_min=None, wave_max=None, pad_frac=0.1)[source]

Reads in the telluric PCA components from a file.

Optionally, this method also trims wavelength to be in within wave_min and wave_max and pads the data (see pad_frac).

Parameters:
  • filename (str) – Telluric PCA filename

  • wave_min (float, optional) – Minimum wavelength at which the grid is desired

  • wave_max (float, optional) – Maximum wavelength at which the grid is desired.

  • pad_frac (float, optional) – Percentage padding to be added to the grid boundaries if wave_min or wave_max are input; ignored otherwise. The resulting grid will extend from (1.0 - pad_frac)*wave_min to (1.0 + pad_frac)*wave_max.

Returns:

Dictionary containing the telluric PCA components.
  • wave_grid: Telluric model wavelength grid

  • dloglam: Wavelength sampling of telluric model

  • tell_pad_pix: Number of pixels to pad at edges for convolution

  • ncomp_tell_pca: Number of PCA components

  • tell_pca: PCA component vectors

  • bounds_tell_pca: Maximum/minimum coefficient

  • coefs_tell_pca: Set of model coefficient values (for prior in future)

  • teltype: Type of telluric model, i.e. ‘pca’

Return type:

dict

pypeit.core.telluric.save_coadd1d_tofits(outfile, wave, flux, ivar, gpm, wave_grid_mid=None, spectrograph=None, telluric=None, obj_model=None, header=None, ex_value='OPT', overwrite=True)[source]

Write final spectrum to disk using the OneSpec write command

Parameters:
  • outfile (str) – File to output telluric corrected spectrum to.

  • wave (numpy.ndarray) – wavelength in units of Angstrom

  • flux (numpy.ndarray) – flux array

  • ivar (numpy.ndarray) – inverse variance array.

  • gpm (numpy.ndarray) – good pixel mask for your spectrum.

  • wave_grid_mid (numpy.ndarray) – wavelength grid evaluated at the bin centers, uniformly-spaced in lambda or log10-lambda/velocity. See core.wavecal.wvutils.py for more info. Default is None because not all files that uses telluric.py has this keyword parameter (spec1dfile does not but coadd1d files do).

  • spectrograph (str, optional) – spectrograph name

  • telluric (numpy.ndarray, optional) – telluric model

  • obj_model (numpy.ndarray, optional) – object model used for telluric fitting

  • header (astropy.io.fits.Header, optional) – Primary header

  • ex_value (str, optional) – extraction mode, OPT or BOX

  • overwrite (bool, optional) – Overwrite existing file?

pypeit.core.telluric.sensfunc_telluric(wave, counts, counts_ivar, counts_mask, exptime, airmass, std_dict, telgridfile, log10_blaze_function=None, ech_orders=None, polyorder=8, tell_npca=4, teltype='pca', mask_hydrogen_lines=True, mask_helium_lines=False, hydrogen_mask_wid=10.0, resln_guess=None, resln_frac_bounds=(0.6, 1.4), pix_shift_bounds=(-5.0, 5.0), delta_coeff_bounds=(-20.0, 20.0), minmax_coeff_bounds=(-5.0, 5.0), sn_clip=30.0, ballsize=0.0005, only_orders=None, maxiter=3, lower=3.0, upper=3.0, tol=0.001, popsize=30, recombination=0.7, polish=True, disp=False, debug_init=False, debug=False)[source]

Compute a sensitivity function from a standard star spectrum by simultaneously fitting a polynomial sensitivity function and a telluric model for atmospheric absorption.

This method is primarily used with SensFunc to compute sensitivity functions.

The sensitivity function is defined to be \(S_\lambda = F_\lambda/({\rm counts/s/A})\) where \(F_\lambda\) is in units of 1e-17 \({\rm erg/s/cm^2/A}\), and so the sensitivity function has units of \({\rm (erg/s/cm^2/A)/(counts/s/A) = erg/cm^2/counts}\)

Parameters:
  • wave (numpy.ndarray, shape is (nspec,) or (nspec, norders)) – Wavelength array. If array is 2D, the data is likely from an echelle spectrograph.

  • counts (numpy.ndarray, shape must match wave) – Counts for the object in question.

  • counts_ivar (numpy.ndarray, shape must match wave) – Inverse variance for the object in question.

  • counts_mask (numpy.ndarray, shape must match wave) – Good pixel mask for the object in question.

  • exptime (float) – Exposure time

  • airmass (float) – Airmass of the observation

  • std_dict (dict) – Dictionary containing the information for the true flux of the standard star.

  • log10_blaze_function (numpy.ndarray , optional) – The log10 blaze function determined from a flat field image. If this is passed in the sensitivity function model will be a (parametric) polynomial fit multiplied into the (non-parametric) log10_blaze_function. Shape must match wave, i.e. (nspec,) or (nspec, norddet).

  • telgridfile (str) – File containing grid of HITRAN atmosphere models; see TelluricPar.

  • ech_orders (numpy.ndarray, shape is (norders,), optional) – If passed, provides the true order numbers for the spectra provided.

  • polyorder (int, optional, default = 8) – Polynomial order for the sensitivity function fit.

  • teltype (str, optional, default = ‘pca’) – Method for evaluating telluric models, either pca or grid.

  • tell_npca (int, optional, default = 4) – Number of telluric PCA components used, must be <= 10

  • mask_hydrogen_lines (bool, optional) – If True, mask stellar hydrogen absorption lines before fitting sensitivity function. Default = True

  • mask_helium_lines (bool, optional) – If True, mask stellar helium absorption lines before fitting sensitivity function. Default = False

  • hydrogen_mask_wid (float, optional) – Parameter describing the width of the mask for or stellar absorption lines (i.e., mask_hydrogen_lines=True) in Angstroms. A region equal to hydrogen_mask_wid on either side of the line center is masked. Default = 10A

  • resln_guess (float, optional) – A guess for the resolution of your spectrum expressed as lambda/dlambda. The resolution is fit explicitly as part of the telluric model fitting, but this guess helps determine the bounds for the optimization (see resln_frac_bounds). If not provided, the wavelength sampling of your spectrum will be used and the resolution calculated using a typical sampling of 3 spectral pixels per resolution element.

  • resln_frac_bounds (tuple, optional) – A two-tuple with the bounds for the resolution fit optimization which is part of the telluric model. This range is in units of resln_guess. For example, (0.5, 1.5) would bound the spectral resolution fit to be within the range (0.5*resln_guess, 1.5*resln_guess).

  • delta_coeff_bounds (tuple, optional, default = (-20.0, 20.0)) – Parameters setting the polynomial coefficient bounds for sensfunc optimization.

  • minmax_coeff_bounds (tuple, optional, default = (-5.0, 5.0)) –

    Parameters setting the polynomial coefficient bounds for sensfunc optimization. Bounds are currently determined as follows. We compute an initial fit to the sensitivity function using init_sensfunc_model(), which determines a set of coefficients. The bounds are then determined according to:

    [(np.fmin(np.abs(this_coeff)*obj_params['delta_coeff_bounds'][0],
            obj_params['minmax_coeff_bounds'][0]),
    np.fmax(np.abs(this_coeff)*obj_params['delta_coeff_bounds'][1],
            obj_params['minmax_coeff_bounds'][1]))]
    

  • sn_clip (float, optional, default=30.0) – Errors are capped during rejection so that the S/N is never greater than sn_clip. This prevents overly aggressive rejection in high S/N ratio spectrum, which neverthless differ at a level greater than the implied S/N due to systematics. These systematics are most likely due to the inadequecy of our telluric model to fit the data at the the precision required at very high S/N ratio.

  • only_orders (list of int, optional, default = None) – Only fit these specific orders

  • tol (float, optional, default = 1e-3) – Relative tolerance for convergence of the differential evolution optimization. See scipy.optimize.differential_evolution for details.

  • popsize (int, optional, default = 30) – A multiplier for setting the total population size for the differential evolution optimization. See scipy.optimize.differential_evolution for details.

  • recombination (float, optional, default = 0.7) – The recombination constant for the differential evolution optimization. This should be in the range [0, 1]. See scipy.optimize.differential_evolution for details.

  • polish (bool, optional, default=True) – If True then differential evolution will perform an additional optimization at the end to polish the best fit, which can improve the optimization slightly. See scipy.optimize.differential_evolution for details.

  • disp (bool, optional, default=True) – Argument for scipy.optimize.differential_evolution, which will display status messages to the screen indicating the status of the optimization. See above for a description of the output and how to know if things are working well.

  • debug_init (bool, optional, default=False) – Show plots to the screen useful for debugging model initialization

  • debug (bool, optional, default=False) – Show plots to the screen useful for debugging the telluric/object model fits.

Returns:

TelObj – Best-fitting telluric model

Return type:

Telluric

pypeit.core.telluric.shift_telluric(tell_model, loglam, dloglam, shift, stretch)[source]

Routine to apply a shift to the telluric model. Note that the shift can be sub-pixel, i.e this routine interpolates. Note that the shift units are pixels of the telluric model.

Parameters:
  • tell_model (numpy.ndarray) – Input telluric model. The shape of this input is in general different from the size of the telluric models (read in by read_telluric_* above) because it is trimmed to relevant wavelengths using ind_lower, ind_upper. See eval_telluric_* below.

  • loglam (numpy.ndarray) – The log10 of the wavelength grid on which the tell_model is evaluated.

  • dloglam (float) – Wavelength spacing of the telluric grid expressed as a a dlog10(lambda), i.e. stored in the tell_dict as tell_dict[‘dloglam’]

  • shift (float) – Desired shift. Note that this shift can be sub-pixel.

  • stretch (float) – Desired stretch.

Returns:

Shifted telluric model. Shape = same size as input tell_model.

Return type:

numpy.ndarray

pypeit.core.telluric.star_telluric(spec1dfile, telgridfile, telloutfile, outfile, star_type=None, star_mag=None, star_ra=None, star_dec=None, func='legendre', model='exp', polyorder=5, teltype='pca', tell_npca=4, mask_hydrogen_lines=True, mask_helium_lines=False, hydrogen_mask_wid=10.0, delta_coeff_bounds=(-20.0, 20.0), minmax_coeff_bounds=(-5.0, 5.0), only_orders=None, sn_clip=30.0, maxiter=3, tol=0.001, popsize=30, recombination=0.7, polish=True, disp=False, pix_shift_bounds=(-5.0, 5.0), debug_init=False, debug=False, show=False, chk_version=True)[source]

This needs a doc string.

USE the one from qso_telluric() as a starting point

Returns:

TelObj – Object with the telluric modeling results

Return type:

Telluric

pypeit.core.telluric.tellfit(flux, thismask, arg_dict, init_from_last=None)[source]

Routine to perform the object + telluric model fitting for telluric corrections. This is a general abstracted routine that performs the fits for any object model that the user provides.

Parameters:
  • flux (numpy.ndarray) – The flux of the object being fit

  • thismask (numpy.ndarray, boolean) – A mask indicating which values are to be fit. This is a good pixel mask, i.e. True=Good

  • arg_dict (dict) –

    A dictionary containing the parameters needed to evaluate the telluric model and the object model. The required keys are:

    • arg_dict['flux_ivar']: Inverse variance for the flux array

    • arg_dict['tell_dict']: Dictionary containing the telluric model and its parameters read in by read_telluric_pca or read_telluric_grid

    • arg_dict['ind_lower']: Lower index into the telluric model wave_grid to trim down the telluric model.

    • arg_dict['ind_upper']: Upper index into the telluric model wave_grid to trim down the telluric model.

    • arg_dict['obj_model_func']: User provided function for evaluating the object model

    • arg_dict['obj_dict']: Dictionary containing the object model arguments which is passed to the obj_model_func

  • init_from_last (object, optional) – Optional. Result object returned by the differential evolution optimizer for the last iteration. If this is passed the code will initialize from the previous best-fit for faster convergence.

Returns:

Returns three objects:

  • result (obj): Result object returned by the differential evolution optimizer

  • modelfit (numpy.ndarray): Modelfit to the input flux. This has the same size as the flux

  • ivartot (numpy.ndarray): Corrected inverse variances for the flux. This has the same size as the flux. The errors are renormalized using the renormalize_errors function by a correction factor, i.e. ivartot = flux_ivar/sigma_corr**2

Return type:

tuple

pypeit.core.telluric.tellfit_chi2(theta, flux, thismask, arg_dict)[source]

Loss function which is optimized by differential evolution to perform the object + telluric model fitting for telluric corrections. This is a general abstracted routine that provides the loss function for any object model that the user provides.

Parameters:
  • theta (numpy.ndarray) –

    Parameter vector for the object + telluric model.

    This is actually two concatenated parameter vectors, one for the object and one for the telluric, i.e.:

    (in PCA mode)

    theta_obj = theta[:-(tell_npca+3)] theta_tell = theta[-(tell_npca+3):]

    (in grid mode)

    theta_obj = theta[:-7] theta_tell = theta[-7:]

    The telluric model theta_tell includes a either user-specified number of PCA coefficients (in PCA mode) or ambient pressure, temperature, humidity, and airmass (in grid mode) followed by spectral resolution, shift, and stretch.

    That is, in PCA mode,

    pca_coeffs = theta_tell[:tell_npca]

    while in grid mode,

    pressure = theta_tell[0] temperature = theta_tell[1] humidity = theta_tell[2] airmass = theta_tell[3]

    with the last three indices of the array corresponding to

    resolution = theta_tell[-3] shift = theta_tell[-2] stretch = theta_tell[-1]

    The object model theta_obj can have an arbitrary size and is provided as an argument to obj_model_func

  • flux (numpy.ndarray) – The flux of the object being fit

  • thismask (numpy.ndarray, boolean) – A mask indicating which values are to be fit. This is a good pixel mask, i.e. True=Good

  • arg_dict (dict) – A dictionary containing the parameters needed to evaluate the telluric model and the object model. See documentation of tellfit for a detailed description.

Returns:

The value of the loss function at the location in parameter space theta. This is loss function is the thing that is minimized to perform the fit.

Return type:

float

pypeit.core.telluric.unpack_orders(sobjs, ret_flam=False)[source]

Utility function to unpack the sobjs object and return the arrays necessary for telluric fitting.

Parameters:
  • sobjs (SpecObjs) – SpecObjs object

  • ret_flam (bool) – If true return the FLAM, otherwise return COUNTS

Returns:

wave, flam, flam_ivar, flam_mask

wave (numpy.ndarray) Wavelength grids; flam (numpy.ndarray) Flambda or counts; flam_ivar (numpy.ndarray) Inverse variance (of Flambda or counts); flam_mask (numpy.ndarray) Good pixel mask. True=Good.

All return values have shape (nspec, norders)

Return type:

tuple