WavelengthCalibration1D¶
- class specreduce.wavecal1d.WavelengthCalibration1D(arc_spectra: Spectrum | Sequence[Spectrum] | None = None, obs_lines: Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str] | Sequence[Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str]] | None = None, line_lists: Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str] | None = None, unit: Unit = Unit('Angstrom'), ref_pixel: float | None = None, pix_bounds: tuple[int, int] | None = None, line_list_bounds: None | tuple[float, float] = None, n_strogest_lines: None | int = None, wave_air: bool = False)[source]¶
 Bases:
objectA class for wavelength calibration of one-dimensional spectra.
This class is designed to facilitate wavelength calibration of one-dimensional spectra, with support for both direct input of line lists and observed spectra. It uses a polynomial model for fitting the wavelength solution and offers features to incorporate catalog lines and observed line positions.
- Parameters:
 - arc_spectra
 Arc spectra provided as
Spectrumobjects for wavelength fitting, by default None. This parameter andobs_linescannot be provided simultaneously.- obs_lines
 Pixel positions of observed spectral lines for wavelength fitting, by default None. This parameter and
arc_spectracannot be provided simultaneously.- line_lists
 Catalogs of spectral line wavelengths for wavelength calibration. Provide either an array of line wavelengths or a list of PypeIt catalog names. If
None, no line lists are used. You can query the list of available catalog names viaget_available_line_catalogs.- unit
 The unit of the wavelength calibration, by default
astropy.units.Angstrom.- ref_pixel
 The reference pixel in which the wavelength solution will be centered.
- pix_bounds
 Lower and upper pixel bounds for fitting, defined as a tuple (min, max). If
obs_linesis provided, this parameter is mandatory.- line_list_bounds
 Wavelength bounds as a tuple (min, max) for filtering usable spectral lines from the provided line lists.
- n_strogest_lines
 The number of strongest lines to be included from the line lists. If
None, all are included.- wave_air
 Boolean indicating whether the input wavelengths correspond to air rather than vacuum; by default
False, meaning vacuum wavelengths.
Attributes Summary
Amplitudes of the catalog lines as a list of masked arrays.
Pixel positions of the catalog lines as a list of masked arrays.
Catalog line wavelengths as a list of masked arrays.
Amplitudes of the observed lines as a list of masked arrays.
Pixel positions of the observed lines as a list of masked arrays.
Pixel positions and amplitudes of the observed lines as a list of masked arrays.
Methods Summary
find_lines(fwhm[, noise_factor])Find lines in the provided arc spectra.
fit_dispersion(wavelength_bounds, ...[, ...])Calculate a wavelength solution using all the catalog and observed lines.
fit_lines(pixels, wavelengths[, degree, ...])Fit the pixel-to-wavelength model using provided line pairs.
match_lines([max_distance])Match the observed lines to theoretical lines.
plot_catalog_lines([frames, axes, figsize, ...])Plot the catalog lines.
plot_fit([frames, figsize, plot_labels, ...])Plot the fitted catalog and observed lines for the specified arc spectra.
plot_observed_lines([frames, axes, figsize, ...])Plot observed spectral lines for the given arc spectra.
plot_residuals([ax, space, figsize])Plot the residuals of pixel-to-wavelength or wavelength-to-pixel transformation.
refine_fit([degree, max_match_distance, ...])Refine the pixel-to-wavelength transformation fit.
Remove unmatched lines from observation and catalog line data.
rms([space])Compute the RMS of the residuals between matched lines in the pixel or wavelength space.
Attributes Documentation
- catalog_lines¶
 Catalog line wavelengths as a list of masked arrays.
- degree¶
 
- observed_lines¶
 Pixel positions and amplitudes of the observed lines as a list of masked arrays.
Methods Documentation
- find_lines(fwhm: float, noise_factor: float = 1.0) None[source]¶
 Find lines in the provided arc spectra.
Determines the spectral lines within each spectrum of the arc spectra based on the provided initial guess for the line Full Width at Half Maximum (FWHM).
- Parameters:
 - fwhm
 Initial guess for the FWHM for the spectral lines, used as a parameter in the
find_arc_linesfunction to locate and identify spectral arc lines.- noise_factor
 The factor to multiply the uncertainty by to determine the noise threshold in the
find_lines_thresholdroutine.
- fit_dispersion(wavelength_bounds: tuple[float, float], dispersion_bounds: tuple[float, float], higher_order_limits: Sequence[float] | None = None, degree: int = 3, popsize: int = 30, max_distance: float = 100, refine_fit: bool = True, refine_max_distance: float = 5.0, refined_fit_degree: int | None = None) WavelengthSolution1D[source]¶
 Calculate a wavelength solution using all the catalog and observed lines.
This method estimates a wavelength solution without pre-matched pixel–wavelength pairs, making it suitable for automated pipelines on stable, well-characterized spectrographs. It uses differential evolution to optimize the polynomial parameters that minimize the distance between the predicted wavelengths of the observed lines and their nearest catalog lines. The resulting solution can optionally be refined with a least-squares fit to automatically matched lines.
- Parameters:
 - wavelength_bounds
 (min, max) bounds for the wavelength at
ref_pixel; used as an optimization constraint.- dispersion_bounds
 (min, max) bounds for the dispersion d(wavelength)/d(pixel) at
ref_pixel; used as an optimization constraint.- higher_order_limits
 Absolute limits for the higher-order polynomial coefficients. Each coefficient is constrained to [-limit, limit]. If provided, the number of limits must equal (polynomial degree - 1).
- degree
 The polynomial degree for the wavelength solution.
- popsize
 Population size for
scipy.optimize.differential_evolution. Larger values can improve the chance of finding the global minimum at the cost of additional time.- max_distance
 Maximum wavelength separation used when associating observed and catalog lines in the optimization. Distances larger than this threshold are clipped to this value in the cost function to limit the impact of outliers.
- refine_fit
 If True (default), call
refine_fitafter global optimization to improve the solution using a least-squares fit on matched lines.- refine_max_distance
 Maximum allowed separation between catalog and observed lines for them to be considered a match during
refine_fit. Ignored ifrefine_fitis False.- refined_fit_degree
 The polynomial degree for the refined fit. Can be higher than
degree. IfNone, equals todegree.
- fit_lines(pixels: Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str], wavelengths: Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | complex | bytes | str | _NestedSequence[complex | bytes | str], degree: int = 3, match_obs: bool = False, match_cat: bool = False, refine_fit: bool = True, refine_max_distance: float = 5.0, refined_fit_degree: int | None = None) WavelengthSolution1D[source]¶
 Fit the pixel-to-wavelength model using provided line pairs.
This method fits the pixel-to-wavelength transformation using explicitly provided pairs of pixel coordinates and their corresponding wavelengths via a linear least-squares fit
Optionally, the provided pixel and wavelength values can be “snapped” to the nearest values present in the internally stored observed line list and catalog line list, respectively. This allows the inputs to be approximate, as the snapping step selects the nearest precise centroids and catalog values when available.
- Parameters:
 - pixels
 An array of pixel positions corresponding to known spectral lines.
- wavelengths
 An array of the same size as
pixels, containing the known wavelengths corresponding to the given pixel positions.- degree
 The polynomial degree for the wavelength solution.
- match_obs
 If True, snap the input
pixelsvalues to the nearest pixel values found inself.observed_line_locations(if available). This helps ensure the fit uses the precise centroids detected byfind_linesor provided initially.- match_cat
 If True, snap the input
wavelengthsvalues to the nearest wavelength values found inself.catalog_line_locations(if available). This ensures the fit uses the precise catalog wavelengths.- refine_fit
 If True (default), automatically call the
refine_fitmethod immediately after the global optimization to improve the solution using a least-squares fit on matched lines.- refine_max_distance
 Maximum allowed separation between catalog and observed lines for them to be considered a match during
refine_fit. Ignored ifrefine_fitis False.- refined_fit_degree
 The polynomial degree for the refined fit. Can be higher than
degree. IfNone, equals todegree.
- match_lines(max_distance: float = 5) None[source]¶
 Match the observed lines to theoretical lines.
- Parameters:
 - max_distance
 The maximum allowed distance between the catalog and observed lines for them to be considered a match.
- plot_catalog_lines(frames: int | Sequence[int] | None = None, axes: Axes | Sequence[Axes] | None = None, figsize: tuple[float, float] | None = None, plot_labels: bool | Sequence[bool] = True, map_to_pix: bool = False, label_kwargs: dict | None = None) Figure[source]¶
 Plot the catalog lines.
- Parameters:
 - frames
 Specifies the frames to be plotted. If an integer, only one frame is plotted. If a sequence, the specified frames are plotted. If None, default selection or all frames are plotted.
- axes
 The matplotlib axes where catalog data will be plotted. If provided, the function will plot on these axes. If None, new axes will be created.
- figsize
 Specifies the dimensions of the figure as (width, height). If None, the default dimensions are used.
- plot_labels
 If True, the numerical values associated with the catalog data will be displayed in the plot. If False, only the graphical representation of the lines will be shown.
- map_to_pix
 Indicates whether the catalog data should be mapped to pixel coordinates before plotting. If True, the data is converted to pixel coordinates.
- label_kwargs
 Specifies the keyword arguments for the line label text objects.
- Returns:
 - Figure
 The matplotlib figure containing the plotted catalog lines.
- plot_fit(frames: int | Sequence[int] | None = None, figsize: tuple[float, float] | None = None, plot_labels: bool = True, obs_to_wav: bool = False, cat_to_pix: bool = False, label_kwargs: dict | None = None) Figure[source]¶
 Plot the fitted catalog and observed lines for the specified arc spectra.
- Parameters:
 - frames
 The indices of the frames to plot. If
None, all frames from 0 toself.nframes - 1are plotted.- figsize
 Defines the width and height of the figure in inches. If
None, the default size is used.- plot_labels
 If
True, print line locations over the plotted lines. Can also be a list with the same length asframes.- obs_to_wav
 If
True, transform the x-axis of observed lines to the wavelength domain usingself._p2w, if available.- cat_to_pix
 If
True, transforms catalog data points to pixel values before plotting.- label_kwargs
 Specifies the keyword arguments for the line label text objects.
- Returns:
 - matplotlib.figure.Figure
 The figure object containing the generated subplots.
- plot_observed_lines(frames: int | Sequence[int] | None = None, axes: Axes | Sequence[Axes] | None = None, figsize: tuple[float, float] | None = None, plot_labels: bool | Sequence[bool] = True, map_to_wav: bool = False, label_kwargs: dict | None = None) Figure[source]¶
 Plot observed spectral lines for the given arc spectra.
- Parameters:
 - frames
 Specifies the frame(s) for which the plot is to be generated. If None, all frames are plotted. When an integer is provided, a single frame is used. For a sequence of integers, multiple frames are plotted.
- axes
 Axes object(s) to plot the spectral lines on. If None, new axes are created.
- figsize
 Dimensions of the figure to be created, specified as a tuple (width, height). Ignored if
axesis provided.- plot_labels
 If True, plots the numerical values of the observed lines at their respective locations on the graph.
- map_to_wav
 Determines whether to map the x-axis values to wavelengths.
- label_kwargs
 Specifies the keyword arguments for the line label text objects.
- Returns:
 - Figure
 The matplotlib figure containing the observed lines plot.
- plot_residuals(ax: Axes | None = None, space: Literal['pixel', 'wavelength'] = 'wavelength', figsize: tuple[float, float] | None = None) Figure[source]¶
 Plot the residuals of pixel-to-wavelength or wavelength-to-pixel transformation.
- Parameters:
 - ax
 Matplotlib Axes object to plot on. If None, a new figure and axes are created.
- space
 The reference space used for plotting residuals. Options are ‘pixel’ for residuals in pixel space or ‘wavelength’ for residuals in wavelength space.
- figsize
 The size of the figure in inches, if a new figure is created.
- Returns:
 - matplotlib.figure.Figure
 
- refine_fit(degree: None | int = None, max_match_distance: float = 5.0, max_iter: int = 5) WavelengthSolution1D[source]¶
 Refine the pixel-to-wavelength transformation fit.
Fits (or re-fits) the polynomial wavelength solution using the currently matched pixel–wavelength pairs. Optionally adjusts the polynomial degree, filters matches by a maximum pixel-space separation, and iterates the fit.
- Parameters:
 - degree
 The polynomial degree for the wavelength solution. If
None, the degree previously set by thefit_linesorfit_dispersionmethod will be used.- max_match_distance
 Maximum allowable distance used to identify matched pixel and wavelength data points. Points exceeding the bound will not be considered in the fit.
- max_iter
 Maximum number of fitting iterations.
- remove_unmatched_lines() None[source]¶
 Remove unmatched lines from observation and catalog line data.
- rms(space: Literal['pixel', 'wavelength'] = 'wavelength') float[source]¶
 Compute the RMS of the residuals between matched lines in the pixel or wavelength space.
- Parameters:
 - space
 The space in which to calculate the RMS residual. If ‘wavelength’, the calculation is performed in the wavelength space. If ‘pixel’, it is performed in the pixel space. Default is ‘wavelength’.
- Returns:
 - float