prospect API

prospect

An interactive spectrum viewer.

prospect.coaddcam

Python routines which translate the simple camera-coaddition algorithms used in prospect webpages (js code in js/interp_grid.js and js/coadd_brz_cameras.js).

prospect.coaddcam.coadd_brz_cameras(wave_in, flux_in, noise_in)[source]

Camera-coadd brz spectra.

Translated from js/coadd_brz_cameras.js.

Parameters:
  • wave_in (array-like) – Set of three wavelength arrays corresponding to brz.

  • flux_in (array-like) – Set of three flux arrays corresponding to brz.

  • noise_in (array-like) – Noise arrays for weighting.

Returns:

The coadded wavelength solution, flux and noise.

Return type:

tuple()

Notes

  • In case of a missing arm (or 2), the data are just concatenated.

  • Need to handle case of no noise.

prospect.coaddcam.coaddcam_prospect(spectra)[source]

Camera-coaddition of brz bands in a set of DESI spectra.

This is essentially a wrapper to coadd_brz_cameras().

Parameters:

spectra (Spectra or similar) –

Returns:

wave : 1D[nwave] array of wavelengths flux : 2D[nspec, nwave] array of flux densities ivar : 2D[nspec, nwave] array of inverse variances of flux

Return type:

tuple() (wave, flux, ivar), where

Raises:

RuntimeError – If spectra.bands does not contain ‘b’, ‘r’, ‘z’, and spectra.bands is not ‘brz’

prospect.coaddcam.index_dichotomy(point, grid)[source]

Find nearest index in grid, left from point; use dichotomy method.

Translated from js/interp_grid.js.

Parameters:
  • point (float) – Value to find in grid.

  • grid (array-like) – Values to search.

Returns:

Nearest index.

Return type:

int

prospect.coaddcam.interp_grid(xval, xarr, yarr)[source]

Basic linear interpolation of [xarr, yarr] on point xval.

Translated from js/interp_grid.js.

Parameters:
  • xval (xval) – Interpolate y-value at this point.

  • xarr (array-like) – X, Y data.

  • yarr (array-like) – X, Y data.

Returns:

The y-value corresponding to xval.

Return type:

float

prospect.grid_thumbs

Grid of thumbnail plots displaying a set of spectra

prospect.grid_thumbs.grid_thumbs(spectra, thumb_width, x_range=(3400, 10000), thumb_height=None, resamp_factor=15, ncols_grid=5, titles=None)[source]

Create a bokeh gridplot of thumbnail pictures from spectra - coadd arms - smooth+resample to reduce size of embedded CDS, according to resamp_factor - titles : optional list of titles for each thumb

TODO: Not tested on Spectrum1D objects.

prospect.scripts

Code to support command-line scripts.

prospect.scripts.prospect_pages

Create static html files to visually inspect DESI spectra.

prospect.scripts.prospect_pages._filter_list(args, filter_name)[source]

Make list of ‘filter_name’ from either args.filter_names (values are explicitely given to parser) or args.filter_name_list (read file)

prospect.scripts.prospect_pages.load_spectra_zcat_from_dbentry(db_entry, args, log, with_redrock_version=True)[source]

Load (spectra, zcat, redrock_cat) from a subset of DESI data (a subdirectory). Somewhat similar to ..utilities.load_spectra_zcat_from_targets(), but this function is intended to be used only by prospect_pages script.

prospect.scripts.prospect_pages.page_subset(spectra, nspec_per_page, titlepage_prefix, viewer_params, log, zcat=None, redrock_cat=None, sort_by_targetid=False)[source]

Make a set of prospect html pages from (spectra, zcat, redrock_cat) zcat and redrock_cat must be entry-matched to spectra

prospect.scripts.prospect_std_templates

Script to produce “standard templates”, used in prospect.viewer.load_std_templates

prospect.scripts.prospect_std_templates.main()[source]

Entry-point for command-line scripts.

Returns:

An integer suitable for passing to sys.exit().

Return type:

int

prospect.specutils

Conversions of DESI & SDSS spectra into specutils-compatible objects.

class prospect.specutils.Spectra(*args: Any, **kwargs: Any)[source]

Represents a grouping of spectra.

This class contains an “extended” fibermap that has information about the night and exposure of each spectrum. For each band, this class has the wavelength grid, flux, ivar, mask, and resolution arrays.

Parameters:
  • bands (list) – List of strings used to identify the bands.

  • wave (dict) – Dictionary of arrays specifying the wavelength grid.

  • flux (dict) – Dictionary of arrays specifying the flux for each spectrum.

  • ivar (dict) – Dictionary of arrays specifying the inverse variance.

  • mask (dict, optional) – Dictionary of arrays specifying the bitmask.

  • resolution_data (dict, optional) – Dictionary of arrays specifying the block diagonal resolution matrix. The object for each band must be in one of the formats supported by the Resolution class constructor.

  • fibermap – Extended fibermap to use. If not specified, a fake one is created.

  • meta (dict, optional) – Dictionary of arbitrary properties.

  • extra (dict, optional) – Optional dictionary of dictionaries containing extra floating point arrays. The top-level is a dictionary over bands and each value is a dictionary containing string keys and values which are arrays of the same size as the flux array.

  • single (bool, optional) – If True, store data in memory as single precision.

  • scores – QA scores table.

property bands

the list of valid bands.

Type:

(list)

property ftype

the data type used for floating point numbers.

Type:

(numpy.dtype)

num_spectra()[source]

Get the number of spectra contained in this group.

Returns (int):

Number of spectra contained in this group.

num_targets()[source]

Get the number of distinct targets.

Returns (int):

Number of unique targets with spectra in this object.

select(nights=None, bands=None, targets=None, fibers=None, invert=False)[source]

Select a subset of the data.

This filters the data based on a logical AND of the different criteria, optionally inverting that selection.

Parameters:
  • nights (list) – optional list of nights to select.

  • bands (list) – optional list of bands to select.

  • targets (list) – optional list of target IDs to select.

  • fibers (list) – list/array of fiber indices to select.

  • invert (bool) – after combining all criteria, invert selection.

Returns (Spectra):

a new Spectra object containing the selected data.

target_ids()[source]

Return list of unique target IDs.

The target IDs are sorted by the order that they first appear.

Returns (array):

an array of integer target IDs.

update(other)[source]

Overwrite or append new data.

Given another Spectra object, compare the fibermap information with the existing one. For spectra that already exist, overwrite existing data with the new values. For spectra that do not exist, append that data to the end of the spectral data.

Parameters:

other (Spectra) – the new data to add.

Returns:

nothing (object updated in place).

wavelength_grid(band)[source]

Return the wavelength grid for a band.

Parameters:

band (str) – the name of the band.

Returns (array):

an array containing the wavelength values.

prospect.specutils.read_frame_as_spectra(filename, night=None, expid=None, band=None, single=False)[source]

Read a FITS file containing a Frame and return a Spectra.

A Frame file is very close to a Spectra object (by design), and only differs by missing the NIGHT and EXPID in the fibermap, as well as containing only one band of data.

Parameters:

infile (str) – path to read

Options:

night (int): the night value to use for all rows of the fibermap. expid (int): the expid value to use for all rows of the fibermap. band (str): the name of this band. single (bool): if True, keep spectra as single precision in memory.

Returns (Spectra):

The object containing the data read from disk.

prospect.specutils.read_spPlate(filename, limit=None)[source]

Read a SDSS spPlate file.

Parameters:
  • filename (str) – Name of the spPlate file.

  • limit (int, optional) – If set, only return the first limit spectra.

Returns:

The spectra.

Return type:

Spectrum1D

prospect.specutils.read_spZbest(filename, limit=None)[source]

Read a SDSS spZbest file.

Parameters:
  • filename (str) – Name of the spZbest file.

  • limit (int, optional) – If set, only return the first limit spectra.

Returns:

A Table containing the redshift values and a Spectrum1D object containing the best-fit models.

Return type:

tuple

prospect.specutils.read_spectra(infile, single=False, coadd=None)[source]

Read Spectra object from FITS file.

This reads data written by the write_spectra function. A new Spectra object is instantiated and returned.

Parameters:
  • infile (str) – path to read

  • single (bool) – if True, keep spectra as single precision in memory.

  • coadd (array-like) – if set, coadd all spectra from the provided targetids.

Returns (Spectra):

The object containing the data read from disk.

prospect.specutils.write_spectra(outfile, spec, units=None)[source]

Write Spectra object to FITS file.

This places the metadata into the header of the (empty) primary HDU. The first extension contains the fibermap, and then HDUs are created for the different data arrays for each band.

Floating point data is converted to 32 bits before writing.

Parameters:
  • outfile (str) – path to write

  • spec (Spectra) – the object containing the data

  • units (str) – optional string to use for the BUNIT key of the flux HDUs for each band.

Returns:

The absolute path to the file that was written.

prospect.utilities

Utility functions for prospect.

prospect.utilities._coadd(wave, flux, ivar, rdat)[source]

Return weighted coadd of spectra

Parameters:
  • wave (array-like) – 1D[nwave] array of wavelengths.

  • flux (array-like) – 2D[nspec, nwave] array of flux densities.

  • ivar (array-like) – 2D[nspec, nwave] array of inverse variances of flux.

  • rdat (array-like) – 3D[nspec, ndiag, nwave] sparse diagonals of resolution matrix.

Returns:

The coadded spectrum (wave, outflux, outivar, outrdat).

Return type:

tuple

prospect.utilities.coadd_targets(spectra, targetids=None)[source]

Coadds individual exposures of the same targets; returns new Spectra object

Parameters:
  • spectra (desispec.spectra.Spectra) –

  • targetids (array-like, optional) – Subset of target IDs to keep.

Returns:

Where individual spectra of each target have been combined into a single spectrumper camera.

Return type:

desispec.spectra.Spectra

Notes

Coadds per camera but not across cameras.

prospect.utilities.create_subsetdb(datadir, dirtree_type=None, spectra_type='coadd', tiles=None, nights=None, expids=None, survey_program=None, petals=None, pixels=None, with_zcat=True)[source]

Create a ‘mini-db’ of DESI spectra files, in a given directory tree.

Supports tile-, exposure- and healpix-based directory trees for daily, andes, … to iron. This routine parses directory trees, but does not open any file: it just checks they exist.

Parameters:
  • datadir (string) – Base directory, eg. /global/cfs/cdirs/desi/spectro/redux/iron/healpix

  • dirtree_type (string) – The directory tree and file names must be as listed in the notes below.

  • spectra_type (string, optional) – [c/s]frames are only supported when dirtree_type=’exposures’

  • petals (list, optional) – Filter a list of petal numbers.

  • tiles (list, optional) – Filter a list of tiles.

  • nights (list, optional) – Filter a list of nights (only if dirtree_type=’pernight’ or ‘exposures’).

  • expids (list, optional) – Filter a list of exposures (only if dirtree_type=’perexp’ or ‘exposures’).

  • survey_program (list, optional) – Filter a [survey, program], only if dirtree_type=’healpix’.

  • pixels (list, optional) – Filter a list of Healpix pixels (only if dirtree_type=’healpix’).

  • with_zcat (bool, optional) – If True, filter spectra for which a ‘redrock’ (or ‘zbest’) fits file exists at the same location.

Returns:

Content of the ‘mini-db’:

  • if dirtree_type=’healpix’: [ {‘dataset’:(survey, program), ‘subset’:’pixel’, ‘petals’:[None]}]

  • if dirtree_type=’exposures’: [ {‘dataset’:night, ‘subset’:expid, ‘petals’:[list of petals]}]

  • if dirtree_type=’perexp’: [ {‘dataset’:tile, ‘subset’:expid, ‘petals’:[list of petals]}]

  • else: [ {‘dataset’:tile, ‘subset’:night, ‘petals’:[list of petals]}]

Return type:

dict

Notes

  • dirtree_type must be one of the following:

    • dirtree_type='healpix': {datadir}/{survey}/{program}/{pixel//100}/{pixel}/{spectra_type}-{survey}-{program}-{pixel}.fits

    • dirtree_type='pernight': {datadir}/{tileid}/{night}/{spectra_type}-{petal}-{tile}-{night}.fits

    • dirtree_type='perexp': {datadir}/{tileid}/{expid}/{spectra_type}-{petal}-{tile}-exp{expid}.fits

    • dirtree_type='cumulative': {datadir}/{tileid}/{night}/{spectra_type}-{petal}-{tile}-thru{night}.fits

    • dirtree_type='exposures': {datadir}/{night}/{expid}/{spectra_type}-{band}{petal}-{expid}.fits

    • Note that ‘perexp’ and ‘exposures’ are different.

    • To use blanc/cascades ‘all’ (resp ‘deep’) coadds, use dirtree_type=’pernight’ and nights=[‘all’] (resp [‘deep’]).

    • Extensions can be fits or fits.gz

prospect.utilities.create_targetdb(datadir, subsetdb, dirtree_type=None)[source]

Create a “mini-db” of DESI targetids.

To do so, redrock (or zbest) fits files are read (faster than reading spectra).

Parameters:
  • datadir (string) – No description provided.

  • subsetdb (list) – List of spectra subsets, as produced by create_subsetdb. Format: [ {‘dataset’:dataset, ‘subset’:subset, ‘petal’:petal} ]

  • dirtree_type (string) – See documentation in create_subsetdb. dirtree_type=’exposures’ is not supported here (no redrock file available in that case). Tile-based directory trees for daily, andes, … to everest are supported. Healpix-based directory tree supported for everest.

Returns:

Content of the “mini-db”: { (dataset, subset, petal): [list of TARGETIDs] } where dataset is a tile, night, or a (survey, program) tuple; subset is a night, expid or pixel; and petal is None when dirtree_type=healpix.

Return type:

dict

prospect.utilities.create_zcat_from_redrock_cat(redrock_cat, fit_num=0)[source]

Extract a catalog with unique redshift fits from a redrock catalog containing several fit results per TARGETID

Parameters:
  • redrock_cat (Table) – Catalog with rows as defined in match_rrdetails_to_spectra()

  • fit_num (int, optional) – The (fit_num)th fit in redrock_cat is extracted (default: 0 ie. redrock’s best fit)

Returns:

Table with the following columns: TARGETID, CHI2, COEFF, Z, ZERR, ZWARN, SPECTYPE, SUBTYPE, DELTACHI2.

Return type:

Table

prospect.utilities.file_or_gz(fname)[source]

Returns fname.gz if fname does not exist but fname.gz exists. Returns fname in all other cases.

Note

This function is similar to desispec.io.util.checkgzip, but does not raise an error if files are not found.

prospect.utilities.file_or_gz_exists(fname)[source]

Returns True if fname or fname.gz exists.

prospect.utilities.frames2spectra(frames, nspec=None, startspec=None, with_scores=False, with_resolution_data=False)[source]

Convert list of frames into DESI Spectra object

Parameters:
  • frames (list) – A list of Frame.

  • nspec (int, optional) – No description provided.

  • startspec (int, optional) – If nspec is set, only spectra in range [startspec:nspec+startspec] are kept

  • with_scores (bool, optional) – If True, include merged scores from input frames

  • with_resolution_data (bool, optional) – If True, include frames.resolution_data

Returns:

No description provided.

Return type:

Spectra

prospect.utilities.get_resources(filetype)[source]

Find all HTML template or JavaScript files in the package.

Caches the results for quick access.

Parameters:

filetype ({'templates', 'js'}) – The type of file resource needed.

Returns:

A dictionary mapping filename to the contents of the file.

Return type:

dict

Raises:

ValueError – If filetype is unknown.

prospect.utilities.load_redrock_templates(template_dir=None)[source]

Load redrock templates; redirect stdout because redrock is chatty

prospect.utilities.load_spectra_zcat_from_targets(targetids, datadir, targetdb, dirtree_type=None, with_redrock_details=False, with_redrock_version=True)[source]

Get spectra, redshift catalog and optional detailed Redrock catalog matched to a set of DESI TARGETIDs.

This works using a “mini-db” of targetids, as returned by create_targetdb(). The outputs of this utility can be used directly by viewer.plotspectra(), to inspect a given list of targetids. Output spectra/catalog(s) are sorted according to the input target list. When several spectra are available for a given TARGETID, they are all included in the output, in random order.

Parameters:
  • targetids (list or numpy.ndarray) – List of TARGETIDs, must be int64.

  • datadir (string) – No description provided.

  • dirtree_type (string) – The directory tree and file names must match the types listed in the notes below.

  • targetdb (dict) – Content of the “mini-db”: { (dataset, subset, petal): [list of TARGETIDs] }, see create_targetdb().

  • with_redrock_details (bool, optional) – If True, detailed Redrock output files (.h5 files) are also read

  • with_redrock_version (bool, optional) – If True, a column ‘RRVER’ is appended to the output redshift catalog, as given by HDU0 in redrock/zbest files. This is used by viewer.plotspectra() to track Redrock version in visual inspection files.

Returns:

If with_redrock_details is False (default), returns (spectra, zcat), where spectra is ~desispec.spectra.Spectra and zcat is ~astropy.table.Table. If with_redrock_details is True, returns (spectra, zcat, redrockcat) where redrockcat is ~astropy.table.Table.

Return type:

tuple()

Notes

  • dirtree_type must be one of the following, for “coadd”, “redrock”/”zbest” (.fits/.fits.gz), and “rrdetails”/”redrock” (.h5) files: - dirtree_type='healpix': {datadir}/{survey}/{program}/{pixel//100}/{pixel}/redrock-{survey}-{program}-{pixel}.fits - dirtree_type='pernight': {datadir}/{tileid}/{night}/redrock-{petal}-{tile}-{night}.fits - dirtree_type='perexp': {datadir}/{tileid}/{expid}/redrock-{petal}-{tile}-exp{expid}.fits - dirtree_type='cumulative': {datadir}/{tileid}/{night}/redrock-{petal}-{tile}-thru{night}.fits - To use blanc/cascades ‘all’ (resp ‘deep’) coadds, use dirtree_type=’pernight’ and nights=[‘all’] (resp ‘deep’)

prospect.utilities.match_catalog_to_spectra(zcat_in, spectra, return_index=False)[source]

Creates a subcatalog, matching a set of DESI spectra

Parameters:
  • zcat_in (Table, with TARGETID keys) –

  • spectra (Spectra) –

  • return_index (bool, optional) – If True, returns the list of indices in zcat_in which match spectra

Returns:

A subtable of zcat_in, with rows matching input spectra’s TARGETIDs If return_index is True, returns (subtable, list of indices)

Return type:

Table

Raises:

RuntimeError – If a unique row in zcat_in is not found matching each of spectra’s TARGETIDs

prospect.utilities.match_rrdetails_to_spectra(redrockfile, spectra, Nfit=None)[source]

Creates a Table from a detailed Redrock output fit, matching a list of DESI spectra.

Parameters:
  • redrockfile (str, filename for the detailed Redrock output file (.h5 file)) –

  • spectra (Spectra) –

  • Nfit (int, optional) – Number of best-fits to store in output Table. By default, store all fits available in the detailed Redrock file

Returns:

Table with the following columns: TARGETID, CHI2, DELTACHI2, COEFF, Z, ZERR, ZWARN, SPECTYPE, SUBTYPE. The rows are matched to spectra’s TARGETIDs

Return type:

Table

Raises:

RuntimeError – If a set of Nfit rows in redrockfile is not found matching each of spectra’s TARGETIDs

prospect.utilities.metadata_selection(spectra, mask=None, mask_type=None, gmag_range=None, rmag_range=None, chi2_range=None, snr_range=None, clean_fiberstatus=False, select_bad_fiberstatus=False, with_dirty_mask_merge=False, zcat=None, log=None, return_index=False)[source]

Simple selection of DESI spectra based on various metadata.

Filtering based on the logical AND of requested selection criteria. Note: use X_range=[min, None] to filter X > min, X_range=[None, max] to filter X < max

Parameters:
  • spectra (Spectra) – No description provided.

  • mask (string, optional) – DESI targeting mask to select, eg ‘ELG’. Requires to set mask_type.

  • mask_type (string, optional) – DESI targeting mask category, currently supported: ‘DESI_TARGET’, ‘BGS_TARGET’, ‘MWS_TARGET’, ‘SECONDARY_TARGET’, ‘CMX_TARGET’, ‘SV[1/2/3]_DESI_TARGET’, ‘SV[1/2/3]_BGS_TARGET’, ‘SV[1/2/3]_MWS_TARGET’, ‘SV[1/2/3]_SCND_TARGET’.

  • with_dirty_mask_merge (bool, optional) – Option for specific targeting mask selection in early CMX data, see code…

  • gmag_range (list) – g magnitude range to select, gmag_range = [gmag_min, gmag_max]

  • rmag_range (list) – r magnitude range to select, rmag_range = [rmag_min, rmag_max]

  • snr_range (list) – SNR range to select, snr_range = [snr_min, snr_max]. This filter applies on all B, R and Z bands, from scores.MEDIAN_COADD_SNR_band, or scores.MEDIAN_CALIB_SNR_band if the former is not found.

  • chi2_range (list) – chi2 range to select, chi2_range = [chi2_min, chi2_max]. Requires to set zcat.

  • clean_fiberstatus (bool) – if True, remove spectra with FIBERSTATUS!=0 or COADD_FIBERSTATUS!=0

  • select_bad_fiberstatus (bool) – if True, select only spectra with FIBERSTATUS!=0 or COADD_FIBERSTATUS!=0

  • zcat (Table) – catalog with chi2 information, must be matched to spectra (needed for chi2_range filter).

  • log (optional log.) –

  • (bool) (return_index) –

Returns:

Selected spectra. If no spectra are selected, returns None If return_index is True, returns (selected spectra, array of selected indices)

Return type:

Spectra

prospect.viewer

Run a spectral viewer (plot spectra and show widgets).

Spectra can be:

prospect.viewer.create_model(spectra, zcat, archetype_fit=False, archetypes_dir=None, template_dir=None)[source]

Returns model_wave[nwave], model_flux[nspec, nwave], row matched to zcat, which can be in a different order than spectra. - zcat must be entry-matched to spectra.

prospect.viewer.plotspectra(spectra, zcatalog=None, redrock_cat=None, notebook=False, html_dir=None, title=None, colors=None, with_imaging=True, with_noise=True, with_thumb_tab=True, with_vi_widgets=True, top_metadata=None, vi_countdown=-1, with_thumb_only_page=False, with_coaddcam=True, mask_type='DESI_TARGET', model_from_zcat=True, model=None, num_approx_fits=None, with_full_2ndfit=True, template_dir=None, archetype_fit=False, archetypes_dir=None, std_template_file=None)[source]

Main prospect routine. From a set of spectra, creates a bokeh document used for VI, to be displayed as an HTML page or within a Jupyter notebook.

Parameters:
  • spectra (Spectra or Spectrum1D or SpectrumList or list of Frame) – Input spectra. Spectrum1D are assumed to be SDSS/BOSS/eBOSS. Otherwise DESI spectra or frames is assumed.

  • zcatalog (Table, optional) – Redshift values, matched one-to-one with the input spectra.

  • redrock_cat (Table, optional) – Redrock output (as defined in match_rrdetails_to_spectra()). Entries must be matched one-by-one (in order) to spectra.

  • notebook (bool, optional) – If True, bokeh outputs the viewer to a Jupyter notebook.

  • html_dir (str, optional) – Directory to store the HTML page if notebook is False.

  • title (str, optional) – Title used to name the HTML page / the bokeh figure / the VI file.

  • colors (list, optional) – Customize the curve’s colors: 3 colors should be given, associated respectively to the coadded data, the model and the noise.

  • with_imaging (bool, optional) – If False, don’t include thumb image from https://www.legacysurvey.org/viewer.

  • with_noise (bool, optional) – If False, don’t include uncertainty for each spectrum.

  • with_thumb_tab (bool, optional) – If False, don’t include a tab with spectra thumbnails.

  • with_vi_widgets (bool, optional) – Include widgets used to enter VI information. Set it to False if you do not intend to record VI files.

  • top_metadata (list, optional) – List of metadata to be highlighted in the top (most visible) table. Default values [‘TARGETID’, ‘EXPID’, ‘COADD_NUMEXP’, ‘COADD_EXPTIME’]

  • vi_countdown (int, optional) – If >0, add a countdown widget in the VI panel, with a value in minutes given by vi_countdown`.

  • with_thumb_only_page (bool, optional) – When creating a static HTML (notebook is False), a light HTML page including only the thumb gallery will also be produced.

  • with_coaddcam (bool, optional) – Include camera-coaddition, only relevant for DESI.

  • mask_type (str, optional (default: DESI_TARGET)) – Bitmask type to identify target categories in the spectra. Supported types are in ..utilities.supported_desitarget_masks

  • model_from_zcat (bool, optional) – If True, model spectra will be computed from the input zcatalog.

  • model (tuple(), optional) – If set, use this input set of model spectra instead of computing it from zcatalog. model consists of (mwave, mflux); model must be entry-matched to zcatalog.

  • num_approx_fits (int, optional) – Number of best-fit models to display, if redrock_cat is provided. By default, all best-fit models available in redrock_cat are diplayed.

  • with_full_2ndfit (bool, optional) – If True, the second best-fit model from redrock_cat will be displayed without approximation (no undersampling, full resolution).

  • template_dir (str, optional) – Redrock template directory.

  • archetype_fit (bool, optional) – If True, assume zcatalog derived from redrock --archetypes and plot model accordingly.

  • archetypes_dir (str, optional) – Directory path for archetypes if not RR_ARCHETYPE_DIR.

  • std_template_file (str, optional) – File containing standard templates to display in viewer.

prospect.viewer.cds

Class containing all bokeh’s ColumnDataSource objects needed in viewer.py

class prospect.viewer.cds.ViewerCDS[source]

Encapsulates Bokeh ColumnDataSource objects to be passed to js callback functions.

compute_median_spectra(spectra)[source]

Stores the median value for each spectrum into CDS. Simple concatenation of all values from different bands.

init_coaddcam_spec(spectra, with_noise=True)[source]

Creates column data source for camera-coadded observed spectra Do NOT store all coadded spectra in CDS obj, to reduce size of html files Except for the first spectrum, coaddition is done later in javascript

init_model(model, second_fit=False)[source]

Creates a CDS for model spectrum

init_othermodel(zcatalog)[source]

Initialize CDS for the ‘other model’ curve, from the best fit

load_fit_templates(template_dir=None, nbpts_templates=4000)[source]

Create dict for spectral templates used in Redrock fits. These are used to recompute Redrock’s Nth best-fit spectra on-the-fly in javascript. Templates are resampled in order to limit the size of html pages (and the browser’s CPU usage). This resampling is dictated by parameter nbpts_templates.

load_metadata(spectra, mask_type=None, zcatalog=None, survey='DESI')[source]

Creates column data source for target-related metadata, from fibermap, zcatalog and VI files

load_rrdetails(redrock_cat)[source]

Create dict for detailled redrock outputs. Used to recompute redrock’s Nth best fit spectra on-the-fly in javascript, and display them in a table.

load_spectra(spectra, with_noise=True)[source]

Creates column data source for observed spectra

load_std_templates(std_template_file=None)[source]

Load a dict of “standard” templates. The std template file is data/std_templates.fits. It was created from ../scripts/prospect_std_templates.py.

prospect.viewer.cds._airtovac(w)[source]

Convert air wavelengths to vacuum wavelengths. Don’t convert less than 2000 Å.

Parameters:

w (float) – Wavelength [Å] of the line in air.

Returns:

Wavelength [Å] of the line in vacuum.

Return type:

float

prospect.viewer.layouts

Full bokeh layouts for prospect

class prospect.viewer.layouts.ViewerLayout(plots, widgets, vi_widgets, with_vi_widgets=True)[source]

Layout of plots and widgets for prospect’s GUI

prospect.viewer.plots

Class containing bokeh plots needed for the viewer

class prospect.viewer.plots.ViewerPlots(colors=None)[source]

Encapsulates Bokeh plot-like objects that are part of prospect’s GUI.

add_spectral_lines(viewer_cds, figure='main', label_offset_top=100, label_offset_bottom=5)[source]

Add spectral line markers to plot.

Parameters:
  • viewer_cds (array-like) – Viewer data.

  • figure ({'main', 'zoom'}, optional) – Figure to add spectral lines to.

  • label_offset_top (float, optional) – Offset in y-position for line labels with respect to top of the figure for emission lines.

  • label_offset_bottom (float, optional) – Offset in y-position for line labels with respect to bottom of the figure for absorption lines.

prospect.viewer.plots._cross_hair_points(x_cross, y_cross)[source]

Points to make a cross-hair centered on x_cross,y_cross Return (xs,ys) to be given to bokeh.plotting.figure.multi_line()

prospect.viewer.plots._viewer_urls(spectra, pixscale=0.262, zoom=13, layer='ls-dr9')[source]

Return legacysurvey.org viewer URLs + other infos for the imaging cutout, for all spectra.

Notes

layer does not apply to the JPEG cutout service. pixscale is in arcsec/pixel (default 0.262 is the native DECam scale)

prospect.viewer.vi_widgets

Class containing bokeh widgets related to visual inspection

class prospect.viewer.vi_widgets.ViewerVIWidgets(title, viewer_cds)[source]

Encapsulates Bokeh widgets, and related callbacks, used for VI

prospect.viewer.widgets

Class containing bokeh widgets needed for the viewer (except for VI widgets)

class prospect.viewer.widgets.ViewerWidgets(plots, nspec)[source]
Encapsulates Bokeh widgets, and related callbacks, that are part of prospect’s GUI.

Except for VI widgets

add_metadata_tables(viewer_cds, show_zcat=True, top_metadata=['TARGETID', 'EXPID', 'COADD_NUMEXP', 'COADD_EXPTIME'])[source]
Display object-related informations

top_metadata: metadata to be highlighted in table_a

Note: “short” CDS, with a single row, are used to fill these bokeh tables. When changing object, js code modifies these short CDS so that tables are updated.

prospect.viewer.widgets._metadata_table(table_keys, viewer_cds, table_width=500, shortcds_name='shortcds', selectable=False)[source]

Returns bokeh’s (ColumnDataSource, DataTable) needed to display a set of metadata given by table_keys.