astrocook.core.spectrum module#
- class astrocook.core.spectrum.TableAdapterV2(data_dict)[source]#
Bases:
objectAdapter to wrap the V2 .t dictionary and expose Astropy Table-like properties.
- class astrocook.core.spectrum.SpectrumV2(data: SpectrumDataV2, history: list = None)[source]#
Bases:
objectComposite container for spectral data and operations.
This class wraps the raw data (flux, wavelength, error) and provides immutable operations that return new
SpectrumV2instances. It serves as the primary data structure passed between Recipes and the GUI.- Parameters:
data (SpectrumDataV2) – The core data container holding the Astropy Quantities.
history (list, optional) – A list of strings describing the modification history of this spectrum.
- __init__(data: SpectrumDataV2, history: list = None)[source]#
- property x: Quantity#
The wavelength (or velocity) array.
- property xmin: Quantity#
The lower bound of the wavelength bins.
- property xmax: Quantity#
The upper bound of the wavelength bins.
- property y: Quantity#
The flux density array.
- property dy: Quantity#
The flux error (1-sigma) array.
- property cont: Quantity | None#
Returns the ‘cont’ auxiliary column as an Astropy Quantity, or None if the column does not exist.
- property norm: Quantity | None#
Returns the normalization vector (Continuum * Telluric).
Logic:
If cont exists: - Returns cont * telluric (if telluric exists). - Returns cont (otherwise).
If cont is missing: - Returns 1.0 (if flux median is ~1, assuming pre-normalized). - Returns None (if flux is physical, signaling need for continuum fitting).
- property model: Quantity | None#
Returns the ‘model’ auxiliary column as an Astropy Quantity, or None if the column does not exist.
- property meta: Dict[str, Any]#
A dictionary of metadata (header keywords, redshifts, etc.).
- property t: TableAdapterV2#
Replaces the raw table/dictionary access. Returns a TableAdapterV2 instance to simulate an Astropy Table (V1) for compatibility with GUI/V1 logic (like gui_table.py).
- Type:
Adapter
- has_aux_column(name: str) bool[source]#
Check if an auxiliary column exists.
- Parameters:
name (str) – The name of the column to check (e.g.,
'cont').- Returns:
Trueif the column is present in the auxiliary data,Falseotherwise.- Return type:
bool
- get_column(name: str) Quantity | None[source]#
Get an auxiliary column as an Astropy Quantity.
- Parameters:
name (str) – The name of the auxiliary column to retrieve (e.g.,
'cont').- Returns:
The column data with units, or
Noneif the column does not exist.- Return type:
astropy.units.Quantity or None
- to_v1_spectrum() Spectrum[source]#
Convert the immutable SpectrumV2 back into a mutable SpectrumV1 object.
This method extracts core data arrays and auxiliary columns to reconstruct the legacy object required for saving in the
.acsformat.- Returns:
A legacy Spectrum object containing the current data and metadata.
- Return type:
astrocook.legacy.spectrum.Spectrum
- with_properties(z_em: float, resol: float = 0.0, object_name: str | None = None, overwrite_resol_col: bool = True) SpectrumV2[source]#
Return a NEW SpectrumV2 instance with updated properties.
- Parameters:
z_em (float) – The emission redshift.
resol (float, optional) – Resolution (R). Defaults to
0.0.object_name (str, optional) – The object name (header keyword).
overwrite_resol_col (bool, optional) – If
True, updates/creates the ‘resol’ auxiliary column with the new resolution value. Defaults toTrue.
- Returns:
A new instance with updated metadata.
- Return type:
- sanitize_legacy_tellurics() SpectrumV2[source]#
Patches legacy spectra where ‘cont’ included tellurics and ‘cont_no_telluric’ stored the intrinsic continuum.
Logic:
If ‘cont_no_telluric’ exists: - Overwrite ‘cont’ with ‘cont_no_telluric’ (Intrinsic). - Remove ‘cont_no_telluric’.
Returns a new SpectrumV2 instance (or self if no changes needed).
- get_region_bounds(region: str) Tuple[Quantity, Quantity][source]#
Get observed wavelength bounds for a standard spectral region.
Calculates the observed wavelength range for a preset region (e.g., Ly-alpha forest) based on the object’s emission redshift.
- Parameters:
region (str) – The region identifier. Supported presets: -
'lya_forest': Ly-beta to Ly-alpha. -'all_ly_forest': Ly-limit to Ly-alpha. -'red_side': Redward of Ly-alpha.- Returns:
A tuple
(obs_min, obs_max)as Astropy Quantities.- Return type:
tuple
- Raises:
ValueError – If
z_emis 0.0 or the region preset is unknown.
- update_column(name: str, values: ndarray, unit: Unit | None = None, description: str = '') SpectrumV2[source]#
Return a NEW SpectrumV2 instance with the specified column updated or added.
Handles both core columns (y, dy) and auxiliary columns (cont, model, etc.).
- Parameters:
name (str) – Name of the column to update (e.g.,
'cont','model').values (np.ndarray) – The numerical data for the column.
unit (astropy.units.Unit, optional) – The physical unit. If None, inherits from existing column or defaults to flux unit.
description (str, optional) – Description of the column content.
- Returns:
A new instance with the updated column data.
- Return type:
- remove_meta_keys(keys: List[str]) SpectrumV2[source]#
Remove specified keys from the metadata dictionary.
Returns a new SpectrumV2 instance with the specified metadata keys removed. This is commonly used to clean up temporary large objects (like JSON mapping dictionaries) before saving the session to disk.
- Parameters:
keys (list of str) – A list of metadata keys to remove.
- Returns:
A new instance with the cleaned metadata.
- Return type:
- remove_columns(col_names: List[str]) SpectrumV2[source]#
Return a NEW SpectrumV2 with specified auxiliary columns removed.
- Parameters:
col_names (list of str) – Names of the columns to delete.
- Returns:
A new instance with the columns removed.
- Return type:
- prepare_expression_context(expression: str, other_specs: Dict[str, SpectrumV2]) Tuple[str, Dict[str, ndarray]][source]#
Prepare context for cross-spectrum arithmetic expressions.
Parses an expression for aliases (e.g.,
"s1.F"), resolves them against provided spectrum objects, and resamples their data onto the current spectrum’s grid.- Parameters:
expression (str) – The raw expression string (e.g.,
"F - s1.F").other_specs (dict) – A dictionary mapping alias names to
SpectrumV2objects.
- Returns:
A tuple containing: 1. final_expr (str): The processed expression with safe variable names. 2. extra_vars (dict): A dictionary mapping variable names to resampled numpy arrays.
- Return type:
tuple
- apply_expression(target_col: str, expression: str, extra_vars: Dict[str, ndarray] = None) SpectrumV2[source]#
Apply a numerical expression to columns and return a new spectrum.
- Parameters:
target_col (str) – The name of the column to create or overwrite.
expression (str) – A NumExpr-compatible string (e.g.
"y / cont").extra_vars (dict, optional) – Additional variables to pass to the evaluation context.
- Returns:
A new instance with the computed column.
- Return type:
- mask_expression(target_col: str, expression: str, extra_vars: Dict[str, ndarray] = None) SpectrumV2[source]#
Mask a column based on a boolean expression (sets values to NaN).
- Parameters:
target_col (str) – The name of the column to mask.
expression (str) – A boolean expression (e.g.
"x > 500"). True values will be masked.extra_vars (dict, optional) – Additional context variables.
- Returns:
A new instance with the masked column.
- Return type:
- calc_intersection_expression(others_specs: List[SpectrumV2], z_target: float, trans_self: str, trans_others: List[str], window_kms: float) str | None[source]#
Calculate the wavelength expression for the velocity intersection of spectra.
Determines the wavelength range where this spectrum and a list of others overlap in velocity space, relative to specific transitions.
- Parameters:
others_specs (list of SpectrumV2) – List of other spectra to intersect with.
z_target (float) – Target redshift for the velocity conversion.
trans_self (str) – Transition name for the current spectrum (e.g.,
'Ly_a').trans_others (list of str) – List of transition names for the other spectra, corresponding to
others_specs.window_kms (float) – Symmetric velocity window width in km/s.
- Returns:
A boolean expression string (e.g.,
"(x > 400.0) & (x < 500.0)") defining the intersection, orNoneif no overlap exists.- Return type:
str or None
- stitch(other_specs: List[SpectrumV2], sort: bool = True) SpectrumV2[source]#
Merge this spectrum with a list of other spectra.
Concatenates x, y, dy, and all auxiliary columns.
- Parameters:
other_specs (list of SpectrumV2) – List of other spectrum objects to stitch.
sort (bool, optional) – If True, sorts the final arrays by wavelength. Defaults to
True.
- Returns:
A new merged spectrum.
- Return type:
- split(expression: str, extra_vars: Dict[str, ndarray] = None) SpectrumV2[source]#
Split the spectrum based on a boolean expression.
- Parameters:
expression (str) – A boolean expression defining the rows to keep (e.g.
"x > 121.6").extra_vars (dict, optional) – Additional context variables.
- Returns:
A new instance containing only the rows where the expression evaluates to True.
- Return type:
- scale_flux(scale_factor: float | ndarray) SpectrumV2[source]#
Scale the flux (y, dy) and relevant auxiliary columns (cont, model).
- Parameters:
scale_factor (float or np.ndarray) – Multiplicative factor.
- Returns:
A new instance with scaled flux.
- Return type:
- calculate_running_std(input_col: str, output_col: str, window_pix: int) SpectrumV2[source]#
Calculate a running RMS on a column and save it as a new column.
- Parameters:
input_col (str) – Name of the column to calculate statistics on (e.g.,
'y').output_col (str) – Name of the output column (e.g.,
'running_std').window_pix (int) – Width of the window in pixels.
- Returns:
A new instance with the added RMS column.
- Return type:
- smooth(sigma_kms: float) SpectrumV2[source]#
Apply Gaussian smoothing to all flux-like columns (y, dy, cont, model).
- Parameters:
sigma_kms (float) – Standard deviation of the Gaussian kernel in km/s.
- Returns:
A new instance with smoothed columns.
- Return type:
- smooth_column(target_col: str, sigma_kms: float, renorm_model: bool = False) SpectrumV2[source]#
Apply Gaussian smoothing to a single target column.
- Parameters:
target_col (str) – Name of the column to smooth.
sigma_kms (float) – Standard deviation of the Gaussian kernel in km/s.
renorm_model (bool, optional) – If True and
target_colis'cont', the'model'column will be re-normalized.
- Returns:
A new instance with the smoothed column.
- Return type:
- rebin(xstart: Quantity | None, xend: Quantity | None, dx: Quantity, kappa: float | None, filling: float) SpectrumV2[source]#
Rebin the spectrum to a new uniform grid.
- Parameters:
xstart (Quantity, optional) – Start wavelength. If None, uses data minimum.
xend (Quantity, optional) – End wavelength. If None, uses data maximum.
dx (Quantity) – Step size (bin width). Can be wavelength or velocity.
kappa (float, optional) – Sigma clipping threshold for rejecting outliers in bins.
filling (float) – Value to fill gaps (default: NaN).
- Returns:
A new instance with the rebinned data.
- Return type:
- equalize_to_reference(reference_spec: SpectrumV2, order: int = 0) SpectrumV2[source]#
Scale this spectrum’s flux to match a reference spectrum.
Computes a scaling polynomial by comparing the flux of this spectrum to a reference spectrum, and applies this scaling to all flux-like columns (y, dy, cont, model).
- Parameters:
reference_spec (SpectrumV2) – The reference spectrum object to match.
order (int, optional) – Polynomial order for the scaling model.
0applies a constant scalar factor,-1applies a spline interpolation,>0applies a polynomial of the given order. Defaults to0.
- Returns:
A new instance with scaled flux.
- Return type:
- coadd(others: List[SpectrumV2], xstart: Quantity | None, xend: Quantity | None, dx: Quantity, kappa: float | None = 5.0, equalize_order: int = 0) SpectrumV2[source]#
Co-add this spectrum with a list of other spectra onto a new common grid.
Handles scaling, unit conversion, “drizzle” rebinning, and converts the result back to the original wavelength unit.
- Parameters:
others (list of SpectrumV2) – List of spectra to add.
xstart (Quantity, optional) – Bounds of the new grid.
xend (Quantity, optional) – Bounds of the new grid.
dx (Quantity) – Step size of the new grid.
kappa (float, optional) – Sigma clipping threshold. Defaults to
5.0.equalize_order (int, optional) – Polynomial order for flux scaling (-1=Off, 0=Scalar). Defaults to
0.
- Returns:
The resulting co-added spectrum.
- Return type:
- resample_on_grid(target_grid_spec: SpectrumV2, fill_value: float = nan) SpectrumV2[source]#
Resample this spectrum onto the exact wavelength grid of another spectrum.
- Parameters:
target_grid_spec (SpectrumV2) – The spectrum whose grid defines the target wavelengths.
fill_value (float, optional) – Value to fill outside the interpolation range (default: NaN).
- Returns:
A new instance interpolated to the target grid.
- Return type:
- get_resampled_column(col_name: str, target_x_grid: ndarray, fill_value: float = nan) ndarray[source]#
Lightweight API to interpolate a single column onto a new x-grid.
- Parameters:
col_name (str) – Name of the column to interpolate.
target_x_grid (np.ndarray) – Target grid in nanometers.
fill_value (float, optional) – Value to fill outside the interpolation range. Defaults to
NaN.
- Returns:
Interpolated values.
- Return type:
np.ndarray
- calibrate(magnitudes: str) SpectrumV2[source]#
Calibrate (and warp) the flux to match target photometric magnitudes.
- Parameters:
magnitudes (str) – Comma-separated pairs of
Filter=Mag(e.g.,"SDSS_g=17.5, SDSS_r=17.0").- Returns:
A new instance with calibrated flux.
- Return type:
- find_absorbed(smooth_len_lya: Quantity, smooth_len_out: Quantity, kappa: float, template: bool) SpectrumV2[source]#
Identify absorbed regions using V1 logic.
- Parameters:
smooth_len_lya (Quantity) – Smoothing length for the Ly-alpha forest.
smooth_len_out (Quantity) – Smoothing length outside the forest.
kappa (float) – Sigma clipping threshold.
template (bool) – Use PCA template (Not Implemented).
- Returns:
A new instance with the
abs_maskcolumn added/updated.- Return type:
- fit_continuum(fudge: float | str, smooth_std_kms: float, mask_col: str = 'abs_mask', renorm_model: bool = False, template: bool = False) SpectrumV2[source]#
Fit a continuum to the unabsorbed regions.
Orchestrates: Interpolate -> Auto-Fudge -> Smooth -> Renormalize Model.
- Parameters:
fudge (float or str) – Fudge factor. Can be a float or
'auto'.smooth_std_kms (float) – Smoothing width for the continuum spline in km/s.
mask_col (str, optional) – Name of the mask column to use. Defaults to
'abs_mask'.renorm_model (bool, optional) – If True, re-normalizes the ‘model’ column to the new continuum.
template (bool, optional) – If True, uses a QSO template to trace Lyman-alpha emission.
- Returns:
A new instance with the
contcolumn added/updated.- Return type:
- fit_powerlaw(regions_str: str, kappa: float = 3.0) SpectrumV2[source]#
Fit a power-law continuum to specified rest-frame regions.
- Parameters:
regions_str (str) – Comma-separated list of regions (e.g.
'128.0-129.0, 131.5-132.5').kappa (float, optional) – Sigma-clipping threshold. Defaults to
3.0.
- Returns:
A new instance containing the
cont_plcolumn.- Return type:
- update_continuum_from_knots(knots_x: List[float], knots_y: List[float], renorm_model: bool = True) SpectrumV2[source]#
Update the continuum using a Cubic Spline through provided knots.
- Parameters:
knots_x (list of float) – X-coordinates of the knots.
knots_y (list of float) – Y-coordinates of the knots.
renorm_model (bool, optional) – If True, re-normalizes the ‘model’ column.
- Returns:
A new instance with the updated continuum.
- Return type:
- detect_absorptions(mask_col: str = 'abs_mask', min_pix: int = 3) SpectrumV2[source]#
Detect absorption regions based on a boolean mask.
Groups adjacent pixels marked as ‘absorbed’ into discrete regions.
- Parameters:
mask_col (str, optional) – Name of the mask column. Defaults to
'abs_mask'.min_pix (int, optional) – Minimum pixels to form a region. Defaults to
3.
- Returns:
A new instance with the
region_idcolumn updated.- Return type:
- identify_lines(multiplet_list: List[str], mask_col: str = 'abs_mask', min_pix_region: int = 3, merge_dv: float = 10.0, score_threshold: float = 0.5, bypass_scoring: bool = False, debug_rating: bool = False) SpectrumV2[source]#
Identify absorption lines by correlating regions with multiplet templates.
Uses a de-blending strategy to separate merged features.
- Parameters:
multiplet_list (list of str) – List of multiplets to search for.
mask_col (str, optional) – Name of the mask column.
min_pix_region (int, optional) – Minimum pixels for a region.
merge_dv (float, optional) – Velocity threshold for merging regions (km/s).
score_threshold (float, optional) – R^2 score threshold for identification.
bypass_scoring (bool, optional) – Accept all candidates if True.
debug_rating (bool, optional) – Show debug plots if True.
- Returns:
A new instance with the
abs_idscolumn and identification metadata.- Return type:
- get_velocity_bounds(z_target: float, trans: str) tuple[source]#
Calculate the min/max velocity coverage of this spectrum relative to a target.
- Parameters:
z_target (float) – The target redshift.
trans (str) – Transition name (e.g.,
'Ly_a').
- Returns:
(min_velocity, max_velocity) in km/s. Returns
(-inf, inf)if invalid.- Return type:
tuple
- detect_doublet_z(other_spec: SpectrumV2, trans_self: str, trans_other: str) float[source]#
Find the redshift of maximum coincidence between this spectrum and another.
Calculates the cross-correlation (product of absorption depths) to align two spectra covering different lines of a doublet.
- Parameters:
other_spec (SpectrumV2) – The other spectrum object.
trans_self (str) – Transition in this spectrum (e.g.,
'OVI_1031').trans_other (str) – Transition in the other spectrum (e.g.,
'OVI_1037').
- Returns:
The detected redshift, or 0.0 on failure.
- Return type:
float
- mask_region(mask_col: str, obs_min: Quantity, obs_max: Quantity) SpectrumV2[source]#
Update a mask column, setting it to False outside the specified bounds.
- Parameters:
mask_col (str) – Name of the mask column.
obs_min (Quantity) – Minimum observed wavelength.
obs_max (Quantity) – Maximum observed wavelength.
- Returns:
A new instance with the updated mask.
- Return type:
- merge_column(col_name: str, other_values: ndarray) SpectrumV2[source]#
Merge new values into an existing column.
Strategy: Use new value if it’s non-zero, otherwise keep old value.
- Parameters:
col_name (str) – Name of the column.
other_values (np.ndarray) – Array of new values.
- Returns:
A new instance with merged data.
- Return type:
- update_model(model_flux: ndarray) SpectrumV2[source]#
Update the ‘model’ auxiliary column with new flux values.
- Parameters:
model_flux (np.ndarray) – The new model flux array.
- Returns:
A new instance with the updated model.
- Return type:
- x_convert(z_ref: float, xunit: Unit) SpectrumV2[source]#
Convert the X-axis units and return a NEW SpectrumV2 instance.
- Parameters:
z_ref (float) – Reference redshift for velocity conversion.
xunit (astropy.units.Unit) – Target unit (e.g.,
nm,km/s).
- Returns:
A new instance with the converted X-axis.
- Return type:
- y_convert(yunit: Unit) SpectrumV2[source]#
Convert the Y-axis units and return a NEW SpectrumV2 instance.
- Parameters:
yunit (astropy.units.Unit) – Target flux unit.
- Returns:
A new instance with converted Y-axis and error.
- Return type: