Flux Calibration API#
End-to-end flux calibration: BOSZ template matching, photometric scaling, calibration-vector computation, and QC diagnostics. See the flux calibration design document for the full specification.
Data containers for kspecdr flux calibration.
All containers are plain dataclasses with light validation in __post_init__.
No heavy dependencies — only numpy.
- class kspecdr.fluxcal.containers.CalibrationVector(wavelength: ~numpy.ndarray, cal_factor: ~numpy.ndarray, cal_variance: ~numpy.ndarray, mask: ~numpy.ndarray, meta: ~typing.Dict = <factory>)[source]#
Bases:
objectWavelength-dependent factor converting observed counts to physical flux.
flux_calibrated[i] = counts[i] * cal_factor[i]- Parameters:
wavelength (ndarray, shape (N,)) – Wavelength axis in Angstrom.
cal_factor (ndarray, shape (N,)) – Calibration factor in erg/s/cm²/Å per count (or per count/s).
cal_variance (ndarray, shape (N,)) – Variance on
cal_factor.mask (ndarray of bool, shape (N,)) – True for reliable (unmasked) pixels. Telluric bands and other problematic regions are set to False.
meta (dict) –
Per-star metadata. Recognised keys:
star_name,fiber_id,teff,logg,feh,alpha_m,scale_factor(float) — photometric normalisation,band_residuals(dict) — synth − obs per filter,chi2,ndof,rv_kms.
- class kspecdr.fluxcal.containers.FilterCurve(name: str, wavelength: ndarray, transmission: ndarray, system: str = 'AB', vega_to_ab: float = 0.0)[source]#
Bases:
objectTransmission curve for one photometric bandpass.
- Parameters:
name (str) – Identifier matching a file in
data/filters/(without.dat).wavelength (ndarray, shape (M,)) – Wavelength in Angstrom (converted from µm on load).
transmission (ndarray, shape (M,)) – Fractional transmission, 0–1.
system (str) – Photometric system:
"AB"or"Vega".vega_to_ab (float) – Offset to add to a Vega magnitude to obtain an AB magnitude:
m_AB = m_Vega + vega_to_ab. Zero for AB-system filters.
- class kspecdr.fluxcal.containers.FluxCalibrationResult(combined_vector: ~kspecdr.fluxcal.containers.CalibrationVector, per_star_vectors: ~typing.List[~kspecdr.fluxcal.containers.CalibrationVector], per_star_residuals: ~typing.List[~numpy.ndarray], summary: ~typing.Dict = <factory>)[source]#
Bases:
objectComplete output of the flux calibration procedure for one exposure.
- Parameters:
combined_vector (CalibrationVector) – Final combined calibration curve to apply to all science fibers.
per_star_vectors (list of CalibrationVector) – Individual per-star calibration vectors before combination.
per_star_residuals (list of ndarray) – Fractional residual arrays
(Cal_star − Cal_combined) / Cal_combinedfor each star. Shape (N,) each.summary (dict) –
Exposure-level summary. Keys:
n_stars_used(int),n_stars_rejected(int),rms_scatter(float, fractional),wavelength_range(tuple),per_star_metrics(list of dict).
- combined_vector: CalibrationVector#
- per_star_vectors: List[CalibrationVector]#
- class kspecdr.fluxcal.containers.Photometry(filter_names: ~typing.List[str], magnitudes: ~numpy.ndarray, mag_errors: ~numpy.ndarray, meta: ~typing.Dict = <factory>)[source]#
Bases:
objectBroadband photometric measurements for a single standard star.
All magnitudes are stored in the AB system regardless of the native photometric system of the catalogue. The conversion from Vega to AB is applied when the object is constructed (see
photometry_from_catalog_row()).- Parameters:
filter_names (list of str) – Names matching the keys in
FILTER_INFO(e.g."ps1_g","gaia_g").magnitudes (ndarray, shape (N_bands,)) – AB magnitudes.
nanmarks unavailable bands.mag_errors (ndarray, shape (N_bands,)) – Magnitude errors (1-sigma).
nanor non-positive marks unavailable.meta (dict) – Optional metadata:
ra,dec,objid,teff_catalog,a_gaia(extinction),a_g(extinction in g-band).
- class kspecdr.fluxcal.containers.Spectrum1D(wavelength: ~numpy.ndarray, flux: ~numpy.ndarray, variance: ~numpy.ndarray, mask: ~numpy.ndarray, meta: ~typing.Dict = <factory>)[source]#
Bases:
objectWavelength-calibrated 1D spectrum with uncertainty and pixel mask.
- Parameters:
wavelength (ndarray, shape (N,)) – Wavelength axis in Angstrom.
flux (ndarray, shape (N,)) – Flux values (counts, counts/s, or physical flux density depending on pipeline stage).
variance (ndarray, shape (N,)) – Per-pixel variance, same units as
flux**2.mask (ndarray of bool, shape (N,)) – True for good (usable) pixels.
meta (dict) –
Freeform metadata. Recognised keys:
fiber_id,fiber_name,exptime,bunit,continuum(ndarray) — theoretical continuum carried from BOSZ,rv_kms— applied radial-velocity correction.
- class kspecdr.fluxcal.containers.StellarTemplate(wavelength: ndarray, flux: ndarray, continuum: ndarray, teff: float, logg: float, feh: float, alpha_m: float = 0.0, carbon_m: float = 0.0, vmicro: float = 1.0, atmos_model: str = '', source: str = 'BOSZ2024')[source]#
Bases:
objectA single BOSZ 2024 stellar model spectrum.
- Parameters:
wavelength (ndarray, shape (K,)) – Wavelength in Angstrom (native BOSZ log-λ grid at R = 10,000).
flux (ndarray, shape (K,)) – Surface flux density in erg/s/cm²/Å (= 4π × H, column 3 of file).
continuum (ndarray, shape (K,)) – Theoretical continuum in the same units (= 4π × C, column 4 of file). Useful as a first-pass normalization without fitting a spline.
teff (float) – Effective temperature in K.
logg (float) – Surface gravity log(g) in cgs.
feh (float) – Overall metallicity [M/H].
alpha_m (float) – Alpha-element enhancement [α/M].
carbon_m (float) – Carbon enhancement [C/M].
vmicro (float) – Microturbulent velocity in km/s.
atmos_model (str) – Atmosphere model code:
"mp"(MARCS plane-parallel) or"ap"(ATLAS9 plane-parallel).source (str) – Provenance string — typically the original filename.
- to_spectrum1d() Spectrum1D[source]#
Wrap as
Spectrum1D(variance = 0, mask = all good).
Photometric utilities for kspecdr flux calibration.
Provides:
- AB magnitude ↔ flux density conversions
- Filter curve loading (data/filters/*.dat)
- Synthetic photometry through a bandpass
- ATLAS Refcat2 catalog loading and extraction into Photometry
All magnitudes inside this module and in Photometry
objects are in the AB system. Vega-to-AB corrections for Gaia and 2MASS
bands are applied automatically when reading from a catalogue row.
Filter file format#
Two-column ASCII, whitespace-separated, no header:
wavelength(µm) transmission(0-1)
Internal representation uses Angstrom throughout.
- kspecdr.fluxcal.photometry.ab_mag_to_flam(mag: float | ndarray, wave_eff_ang: float, mag_err: float | ndarray | None = None) tuple[ndarray, ndarray | None][source]#
Convert AB magnitude to f_λ (erg/s/cm²/Å) at effective wavelength.
f_λ = f_ν × c / λ²
- kspecdr.fluxcal.photometry.ab_mag_to_fnu(mag: float | ndarray, mag_err: float | ndarray | None = None) tuple[ndarray, ndarray | None][source]#
Convert AB magnitude(s) to f_ν in erg/s/cm²/Hz.
m_AB = −2.5 log₁₀(f_ν) − 48.60
- kspecdr.fluxcal.photometry.estimate_teff_from_color(photometry: Photometry, method: str = 'bp_rp') tuple[float, tuple[float, float]][source]#
Rough photometric Teff estimate for F/G/K stars from broadband color.
- Parameters:
photometry (Photometry)
method (str) –
Color index to use. Supported:
"bp_rp"— Gaia BP − RP (recommended; ~100 K precision for F stars)"g_r"— PS1 g − r
- Returns:
teff (float) – Estimated effective temperature in K.
teff_range ((float, float)) – Conservative ±2-sigma search window (teff_lo, teff_hi).
- kspecdr.fluxcal.photometry.fnu_to_ab_mag(fnu: float | ndarray, fnu_err: float | ndarray | None = None) tuple[ndarray, ndarray | None][source]#
Convert f_ν (erg/s/cm²/Hz) to AB magnitude(s).
- kspecdr.fluxcal.photometry.load_filter_curve(filter_name: str) FilterCurve[source]#
Load a filter transmission curve from
data/filters/{filter_name}.dat.Results are cached in memory after the first load.
- Parameters:
filter_name (str) – One of the keys in
FILTER_INFO, e.g."ps1_g".- Returns:
Wavelength axis in Angstrom, transmission in [0, 1].
- Return type:
- Raises:
KeyError – If filter_name is not in
FILTER_INFO.FileNotFoundError – If the corresponding
.datfile is missing.
- kspecdr.fluxcal.photometry.load_filter_curves(filter_names: Sequence[str]) Dict[str, FilterCurve][source]#
Load multiple filters at once.
- Return type:
dict mapping filter name →
FilterCurve
- kspecdr.fluxcal.photometry.load_standard_star_catalog(catalog_path: str | Path)[source]#
Load an ATLAS Refcat2-format standard-star CSV catalog.
- Parameters:
catalog_path (str or Path) – Path to the CSV file. Must have a header row matching the Refcat2 column names (see
resources/comm/.../standard_star_atlas_refcat2.csv).- Returns:
Full catalog table with all photometric columns. Columns with
99.9sentinel values are not filtered here; usephotometry_from_catalog_row()to obtain validatedPhotometryobjects.- Return type:
- kspecdr.fluxcal.photometry.photometry_from_catalog_row(row, bands: List[str] | None = None) Photometry[source]#
Extract a
Photometryobject from one Refcat2 catalogue row.Vega-to-AB corrections for Gaia and 2MASS bands are applied automatically. Bands with sentinel values (magnitude ≥ 99 or error ≥ 9.9, or missing) are stored as
nanand excluded byvalid_bands().- Parameters:
row (astropy.table.Row or dict-like) – One row from a table returned by
load_standard_star_catalog().bands (list of str, optional) – Which filter names to extract. Defaults to
DEFAULT_BANDS.
- Returns:
Magnitudes in the AB system.
- Return type:
- kspecdr.fluxcal.photometry.query_photometry(ra: float, dec: float, radius_arcsec: float = 2.0, catalog: str = 'refcat2') Photometry[source]#
Query an external catalogue for standard-star photometry.
Note
Not yet implemented. Use
load_standard_star_catalog()andphotometry_from_catalog_row()with a pre-downloaded CSV file.- Raises:
- kspecdr.fluxcal.photometry.synthetic_color(spectrum: Spectrum1D, filter1: FilterCurve, filter2: FilterCurve) float[source]#
Return synthetic color
m_filter1 − m_filter2(AB mags).
- kspecdr.fluxcal.photometry.synthetic_photometry(spectrum: Spectrum1D, filter_curve: FilterCurve) float[source]#
Compute synthetic AB magnitude of spectrum through filter_curve.
Uses the photon-counting (energy-per-photon weighted) convention:
\[ \begin{align}\begin{aligned}\langle f_\nu \rangle = \frac{\int f_\nu(\lambda)\, T(\lambda)\, \lambda\, d\lambda} {c \int T(\lambda) / \lambda\, d\lambda}\\m_{\rm AB} = -2.5 \log_{10}(\langle f_\nu \rangle) - 48.60\end{aligned}\end{align} \]where \(f_\nu = f_\lambda \, \lambda^2 / c\).
- Parameters:
spectrum (Spectrum1D) – Observed or model spectrum with
fluxin erg/s/cm²/Å andwavelengthin Angstrom. Masked pixels are excluded.filter_curve (FilterCurve) – Bandpass transmission.
- Returns:
Synthetic AB magnitude. Returns
np.nanif the filter is not covered by the spectrum or if the result is non-physical.- Return type:
- Raises:
ValueError – If spectrum flux and wavelength have mismatched shapes.
BOSZ 2024 stellar template library for kspecdr flux calibration.
Provides:
- parse_bosz_filename() — extract stellar parameters from filename
- TemplateLibrary — indexed access to the downloaded BOSZ subgrid
- prepare_template() — convolve to instrument resolution and resample
- resample_spectrum() — general-purpose 1-D resampling
BOSZ file format (_resam.txt.gz)#
Two-column ASCII (whitespace-separated, no header), one row per wavelength
point. The shared wavelength grid is in a separate file
(bosz2024_wave_r10000.txt).
Column 0:
H(erg/s/cm²/Å/steradian)Column 1:
C(continuum; erg/s/cm²/Å/steradian)
Surface flux: F = 4π × H. Normalised flux: F/C = H/C.
- class kspecdr.fluxcal.templates.TemplateLibrary(library_dir: str | Path | None = None, resolution: str = 'r10000', wave_range: Tuple[float, float] = (3000.0, 11000.0))[source]#
Bases:
objectIndexed library of BOSZ 2024 stellar template spectra.
On first instantiation the library directory is scanned and an index CSV is written to
<library_dir>/index_<resolution>.csvfor fast subsequent starts.- Parameters:
library_dir (str or Path, optional) – Root of the BOSZ subgrid. Defaults to
data/templates/bosz2024/.resolution (str) – Resolution subdirectory. Default
"r10000".wave_range ((float, float), optional) – If given, trim the shared wavelength grid (and all loaded spectra) to this range in Angstrom, reducing memory. Default
(3000, 11000).
- get_template(teff: float, logg: float, feh: float, alpha_m: float | None = None) StellarTemplate[source]#
Nearest-grid-point lookup and load.
- Parameters:
- Return type:
- load_template(entry: Dict) StellarTemplate[source]#
Read a template spectrum from disk.
- Parameters:
entry (dict) – One element from
_index(must have"filepath"key).- Return type:
- query(teff_range: Tuple[float, float] | None = None, logg_range: Tuple[float, float] | None = None, feh_range: Tuple[float, float] | None = None) List[Dict][source]#
Return index entries within the given parameter box.
Entries are metadata dicts (cheap); call
load_template()to read the actual spectral data.
- kspecdr.fluxcal.templates.parse_bosz_filename(filename: str) Dict | None[source]#
Parse a BOSZ filename into a parameter dictionary.
- kspecdr.fluxcal.templates.prepare_template(template: StellarTemplate, target_wavelength: ndarray, instrument_fwhm_angstrom: float, fwhm_poly_coeffs: Sequence[float] | None = None) Spectrum1D[source]#
Convolve a template to instrument resolution and resample.
The BOSZ templates are provided at R = 10,000, which is higher than the typical KSPEC instrument resolution (R ~ 3,000–5,000). This function:
Computes the additional broadening needed (quadrature subtraction of the BOSZ native width from the requested instrument width).
Applies a Gaussian convolution (constant or wavelength-dependent).
Resamples onto the observed wavelength grid.
Carries the BOSZ continuum through (same operations applied).
- Parameters:
template (StellarTemplate) – A template loaded from
TemplateLibrary.target_wavelength (ndarray, shape (N,)) – Observed wavelength grid in Angstrom.
instrument_fwhm_angstrom (float) – Instrument spectral FWHM in Angstrom (
SPECFWHMheader keyword). Used when fwhm_poly_coeffs is None.fwhm_poly_coeffs (sequence of float, optional) – Polynomial coefficients giving FWHM(λ) in Angstrom as
FWHM(λ) = c0 + c1*λ + c2*λ² + ...(from wavecalFWHMPCnheader keywords). If provided, a wavelength-dependent convolution is applied.
- Returns:
Template on the observed grid.
variance = 0,mask = Trueeverywhere.meta['continuum']contains the resampled continuum.- Return type:
- kspecdr.fluxcal.templates.resample_spectrum(wave_in: ndarray, flux_in: ndarray, wave_out: ndarray, fill_value: float = 0.0) ndarray[source]#
Resample a 1-D spectrum onto a new wavelength grid.
Uses linear interpolation. Out-of-range pixels are filled with fill_value.
- Parameters:
wave_in (ndarray, shape (M,))
flux_in (ndarray, shape (M,))
wave_out (ndarray, shape (N,))
fill_value (float)
- Return type:
ndarray, shape (N,)
Continuum normalization for kspecdr flux calibration.
Provides iterative pseudo-continuum fitting with lower-sigma clipping to avoid absorption lines biasing the fit downward. Three methods:
"bspline"— LSQ B-spline (default; robust for broad coverage)"polynomial"— Legendre polynomial (simpler, good for narrow ranges)"running_median"— Smoothed running median (fast, non-parametric)
All methods operate on Spectrum1D objects.
For BOSZ templates, the theoretical continuum column
(StellarTemplate.continuum) can be used directly, bypassing the fit
entirely — see normalize_with_model_continuum().
- kspecdr.fluxcal.continuum.fit_continuum(spectrum: Spectrum1D, method: str = 'bspline', order: int = 3, n_knots: int = 20, sigma_lo: float = 2.0, sigma_hi: float = 4.0, n_iter: int = 5, mask_regions: List[Tuple[float, float]] | None = None, median_window: int = 151) ndarray[source]#
Fit the pseudo-continuum and return it as an array.
See
normalize_continuum()for parameter descriptions.- Returns:
Fitted continuum evaluated at every pixel.
- Return type:
ndarray, shape (N,)
- kspecdr.fluxcal.continuum.normalize_continuum(spectrum: Spectrum1D, method: str = 'bspline', order: int = 3, n_knots: int = 20, sigma_lo: float = 2.0, sigma_hi: float = 4.0, n_iter: int = 5, mask_regions: List[Tuple[float, float]] | None = None, median_window: int = 151) Tuple[Spectrum1D, ndarray][source]#
Fit and divide out the pseudo-continuum.
- Parameters:
spectrum (Spectrum1D) – Input spectrum (flux in any units).
method (
"bspline"|"polynomial"|"running_median") – Fitting method.order (int) – B-spline degree (for
"bspline") or polynomial degree (for"polynomial"). Ignored for"running_median".n_knots (int) – Number of interior knots for B-spline. Ignored for other methods.
sigma_lo (float) – Lower rejection threshold in sigma units. Pixels more than sigma_lo × σ below the current continuum fit are clipped (targeting absorption lines).
sigma_hi (float) – Upper rejection threshold. Pixels more than sigma_hi × σ above the continuum are clipped (targeting emission artefacts or cosmic-ray residuals).
n_iter (int) – Number of iterative sigma-clipping passes.
mask_regions (list of (lam_lo, lam_hi), optional) – Wavelength regions to exclude from the fit entirely (e.g. tellurics). Pixels in these regions still receive interpolated continuum values.
median_window (int) – Window size in pixels for
"running_median". Must be odd.
- Returns:
normalized (Spectrum1D) – Continuum-normalised spectrum (
flux / continuum). Variance is propagated.meta["continuum_fit"]stores the fitted continuum.continuum (ndarray, shape (N,)) – The fitted continuum array.
- kspecdr.fluxcal.continuum.normalize_with_model_continuum(spectrum: Spectrum1D, model_continuum: ndarray) Tuple[Spectrum1D, ndarray][source]#
Normalise using an externally supplied continuum (e.g. BOSZ
Ccolumn).- Parameters:
spectrum (Spectrum1D)
model_continuum (ndarray, shape (N,))
- Returns:
normalized (Spectrum1D)
continuum (ndarray)
Template selection and radial-velocity measurement for kspecdr flux calibration.
Workflow for each standard star:
Narrow the template search range using a photometric Teff estimate.
For each candidate template: a. Prepare (convolve + resample) to the observed grid. b. Cross-correlate to measure/correct radial velocity. c. Continuum-normalise both observed and template. d. Score the fit on line features (χ² or Huber metric).
Return the best-matching template, RV, and fit statistics.
- kspecdr.fluxcal.matching.cross_correlate_rv(observed: Spectrum1D, template: Spectrum1D, max_shift_kms: float = 500.0, rv_guess: float = 0.0) float[source]#
Measure radial velocity by cross-correlating observed and template.
Both spectra are resampled onto a uniform log-λ grid (= uniform velocity bins) before computing the FFT-based cross-correlation. The peak is refined by quadratic interpolation for sub-pixel precision.
- Parameters:
observed (Spectrum1D) – Continuum-normalised observed spectrum.
template (Spectrum1D) – Continuum-normalised (and convolved) template on the same grid.
max_shift_kms (float) – Maximum allowed shift in km/s.
rv_guess (float) – Initial guess (shifts the search window centre).
- Returns:
Best-fit radial velocity in km/s (positive = receding).
- Return type:
- kspecdr.fluxcal.matching.score_template_fit(observed: Spectrum1D, template: Spectrum1D, metric: str = 'chi2', mask_regions: List[Tuple[float, float]] | None = None, huber_delta: float = 1.5) Tuple[float, int][source]#
Score the agreement between normalised observed and template spectra.
- Parameters:
observed (Spectrum1D) – Continuum-normalised spectra on the same wavelength grid.
template (Spectrum1D) – Continuum-normalised spectra on the same wavelength grid.
metric (
"chi2"|"huber")mask_regions (list of (lo, hi), optional) – Additional regions to exclude from the score.
huber_delta (float) – Transition point for the Huber loss function.
- Returns:
score (float) – Reduced chi² or mean Huber loss.
ndof (int) – Number of pixels used.
- kspecdr.fluxcal.matching.select_best_template(observed: Spectrum1D, photometry: Photometry, library: TemplateLibrary, instrument_fwhm_angstrom: float, fwhm_poly_coeffs: Sequence[float] | None = None, rv_guess: float = 0.0, mask_regions: List[Tuple[float, float]] | None = None, teff_range: Tuple[float, float] | None = None, logg_range: Tuple[float, float] | None = None, feh_range: Tuple[float, float] | None = None, metric: str = 'chi2', max_rv_shift: float = 500.0, continuum_method: str = 'bspline', continuum_n_knots: int = 20) Tuple[StellarTemplate, float, Dict][source]#
Find the best-matching template for an observed standard star.
- Parameters:
observed (Spectrum1D) – Extracted, wavelength-calibrated observed spectrum.
photometry (Photometry) – Broadband magnitudes for the star (AB system).
library (TemplateLibrary) – Loaded BOSZ template library.
instrument_fwhm_angstrom (float) – Instrument spectral FWHM in Å.
fwhm_poly_coeffs (sequence of float, optional) – Wavelength-dependent FWHM polynomial.
rv_guess (float) – Initial RV guess in km/s.
mask_regions (list of (lo, hi), optional) – Wavelength regions to exclude from scoring.
teff_range ((float, float), optional) – Explicit Teff search range. If None, estimated from photometry.
logg_range ((float, float), optional) – Explicit search ranges.
feh_range ((float, float), optional) – Explicit search ranges.
metric (
"chi2"|"huber") – Scoring metric.max_rv_shift (float) – Maximum allowed RV shift in km/s for cross-correlation.
continuum_method (str) – Continuum fitting method for observed spectrum.
continuum_n_knots (int) – Number of knots for B-spline continuum fit.
- Returns:
best_template (StellarTemplate) – The best-matching template (un-convolved, native BOSZ grid).
best_rv (float) – Measured radial velocity in km/s.
fit_stats (dict) – Keys:
"best_score","ndof","runner_up_score","n_candidates","teff_range","best_params","all_scores"(list of dicts).
Calibration vector computation and application for kspecdr flux calibration.
Provides:
- scale_template_to_photometry() — absolute flux anchor via synthetic photometry
- compute_calibration_vector_for_star() — full per-star calibration pipeline
- combine_calibration_vectors() — robust combination of per-star vectors
- apply_flux_calibration() — apply calibration to all fibers with variance propagation
This is the top-level orchestration module (P2) that ties together P0 and P1.
- kspecdr.fluxcal.calibration.apply_flux_calibration(spectra: ndarray, variance: ndarray, calibration: FluxCalibrationResult, update_header: Dict | None = None) Tuple[ndarray, ndarray, Dict][source]#
Apply the combined calibration vector to all fibers.
\[ \begin{align}\begin{aligned}F_{\rm cal}(\lambda) = C_{\rm obs}(\lambda) \times \mathrm{Cal}(\lambda)\\\sigma^2_{\rm cal} = \sigma^2_{\rm obs} \times \mathrm{Cal}^2 + C_{\rm obs}^2 \times \sigma^2_{\rm Cal}\end{aligned}\end{align} \]- Parameters:
spectra (ndarray, shape
(NFIB, NPIX)or(NPIX, NFIB)) – Observed spectra in counts (or counts/s). The calibration vector length must match one of the two axes.variance (ndarray, same shape as spectra) – Variance array.
calibration (FluxCalibrationResult) – Output of
combine_calibration_vectors().update_header (dict, optional) – If provided, FITS header keywords are added/updated in place.
- Returns:
cal_spectra (ndarray, same shape as *spectra*) – Flux-calibrated spectra in erg/s/cm²/Å.
cal_variance (ndarray, same shape as *spectra*) – Propagated variance.
header_updates (dict) – Suggested FITS header updates (BUNIT, HISTORY lines).
- kspecdr.fluxcal.calibration.combine_calibration_vectors(vectors: List[CalibrationVector], method: str = 'weighted_mean', sigma_clip: float = 3.0, smooth: bool = False, smooth_window: int = 51) FluxCalibrationResult[source]#
Combine per-star calibration vectors into a single curve.
Gracefully handles N = 1 (single star: pass-through, no combination).
- Parameters:
vectors (list of CalibrationVector) – Per-star calibration vectors (must share the same wavelength grid).
method (
"weighted_mean"|"median") – Combination method.sigma_clip (float) – Rejection threshold in sigma for iterative clipping.
smooth (bool) – If True, apply Savitzky–Golay smoothing to the combined curve.
smooth_window (int) – Window length for smoothing (must be odd).
- Return type:
- kspecdr.fluxcal.calibration.compute_calibration_vector_for_star(observed: Spectrum1D, photometry: Photometry, library: TemplateLibrary, filter_curves: Dict[str, FilterCurve], instrument_fwhm_angstrom: float, fwhm_poly_coeffs: Sequence[float] | None = None, mask_regions: List[Tuple[float, float]] | None = None, metric: str = 'chi2', star_name: str = '', fiber_id: int = -1) CalibrationVector[source]#
Full per-star calibration pipeline.
Orchestrates:
Select best-matching template (
select_best_template()).Prepare template at instrument resolution.
Scale template to match observed photometry.
Compute
Cal(λ) = scaled_template(λ) / observed(λ).Propagate variance and apply mask.
- Parameters:
observed (Spectrum1D) – Extracted, wavelength-calibrated, sky-subtracted observed spectrum (in counts or counts/s).
photometry (Photometry) – Broadband AB magnitudes for this standard star.
library (TemplateLibrary) – Loaded BOSZ template library.
filter_curves (dict of FilterCurve) – Loaded filter curves.
instrument_fwhm_angstrom (float) – Instrument FWHM in Å.
fwhm_poly_coeffs (sequence of float, optional) – Wavelength-dependent FWHM polynomial.
mask_regions (list of (lo, hi), optional) – Wavelength regions to exclude.
metric (str) – Template scoring metric (
"chi2"or"huber").star_name (str) – Identifier for this star (for logging / metadata).
fiber_id (int) – Fiber index.
- Returns:
The per-star calibration factor
Cal(λ)such thatflux_phys = counts × Cal(λ).- Return type:
- kspecdr.fluxcal.calibration.scale_template_to_photometry(template_spec: Spectrum1D, photometry: Photometry, filter_curves: Dict[str, FilterCurve]) Tuple[float, float, Dict[str, float]][source]#
Scale a template so its synthetic photometry matches observed magnitudes.
The scale factor is determined by minimising the weighted mean squared difference between synthetic and observed magnitudes across all valid bands. Because magnitudes are logarithmic, the scale factor is computed in flux space:
\[s = 10^{-0.4 \, \Delta m}\]where
Δmis the weighted-mean offset (observed − synthetic).- Parameters:
template_spec (Spectrum1D) – Template spectrum on the observed wavelength grid (output of
prepare_template()). Flux in erg/s/cm²/Å (surface flux).photometry (Photometry) – Observed AB magnitudes with errors.
filter_curves (dict of FilterCurve) – Loaded filter curves keyed by filter name.
- Returns:
scale_factor (float) – Multiplicative factor to apply to the template flux.
scale_error (float) – 1-σ uncertainty on scale_factor (propagated from magnitude errors).
band_residuals (dict) –
{filter_name: synth_mag − obs_mag}after scaling.
Wavelength mask utilities for kspecdr flux calibration.
Masks are lists of (lam_lo, lam_hi) pairs in Angstrom that mark regions
to exclude from calibration fits, continuum normalisation, and template
matching (e.g. telluric absorption bands, detector bad regions).
Mask file format#
Two-column ASCII, whitespace-separated. Lines starting with # are
comments and are ignored:
6860.0 6960.0 # O2 B-band
7150.0 7340.0 # H2O
The default mask is loaded from data/masks/telluric_default.dat.
- kspecdr.fluxcal.masks.apply_mask_regions(spectrum: Spectrum1D, regions: List[tuple[float, float]]) Spectrum1D[source]#
Set
mask=Falsefor all pixels inside any exclusion region.Returns a new
Spectrum1Dwith an updated mask; the original is not modified.- Parameters:
spectrum (Spectrum1D)
regions (list of (lam_lo, lam_hi)) – Wavelength intervals to exclude (Angstrom).
- Returns:
New spectrum with pixels inside regions masked out.
- Return type:
- kspecdr.fluxcal.masks.apply_named_mask(spectrum: Spectrum1D, mask_name: str = 'telluric_default') Spectrum1D[source]#
Convenience wrapper: load mask by name and apply to spectrum.
- Parameters:
spectrum (Spectrum1D)
mask_name (str) – Mask file base name (see
load_mask_regions()).
- Return type:
- kspecdr.fluxcal.masks.combine_regions(*region_lists: List[tuple[float, float]]) List[tuple[float, float]][source]#
Merge multiple region lists and sort by start wavelength.
Overlapping or adjacent regions are not merged — they are stacked and applied independently, which produces the same result.
- kspecdr.fluxcal.masks.load_mask_array(wavelength: ndarray, mask_name: str = 'telluric_default') ndarray[source]#
Return a boolean mask array for a given wavelength axis.
- Parameters:
wavelength (ndarray, shape (N,)) – Wavelength axis in Angstrom.
mask_name (str) – Mask file to load (see
load_mask_regions()).
- Returns:
True for good (unmasked) pixels; False inside any exclusion region.
- Return type:
ndarray of bool, shape (N,)
- kspecdr.fluxcal.masks.load_mask_regions(mask_name: str = 'telluric_default') List[tuple[float, float]][source]#
Load wavelength exclusion regions from
data/masks/{mask_name}.dat.Results are cached after the first load.
- Parameters:
mask_name (str) – Base name of the mask file (without
.dat). Default:"telluric_default".- Returns:
Wavelength intervals in Angstrom to be masked (excluded).
- Return type:
list of (lam_lo, lam_hi)
- Raises:
FileNotFoundError – If the mask file does not exist.
- kspecdr.fluxcal.masks.mask_coverage_fraction(wavelength: ndarray, mask_name: str = 'telluric_default') float[source]#
Return the fraction of wavelength pixels inside the named mask.
Useful for QC: large fractions (> 0.3) may indicate a wavelength unit mismatch.
Quality-control plotting utilities for kspecdr flux calibration.
All functions return matplotlib.figure.Figure objects and are designed
for interactive use in Jupyter notebooks. They are not called by the
pipeline itself.
Typical usage:
from kspecdr.fluxcal.qc import plot_calibration_summary
fig = plot_calibration_summary(result)
- kspecdr.fluxcal.qc.plot_calibrated_spectrum(wavelength: ndarray, flux: ndarray, variance: ndarray | None = None, fiber_id: int = -1, title: str = '', figsize: tuple = (12, 4))[source]#
Quick plot of a flux-calibrated spectrum with optional error band.
- kspecdr.fluxcal.qc.plot_calibration_residuals(result: FluxCalibrationResult, figsize: tuple = (12, 4))[source]#
Fractional residuals of each per-star vector relative to the combined.
- Parameters:
result (FluxCalibrationResult)
figsize (tuple)
- Return type:
matplotlib.figure.Figure
- kspecdr.fluxcal.qc.plot_calibration_summary(result: FluxCalibrationResult, title: str = 'Flux Calibration Summary', figsize: tuple = (14, 10))[source]#
Four-panel summary of the flux calibration result.
Panels: 1. Combined calibration vector Cal(λ) 2. Per-star calibration vectors (overlaid) 3. Per-star fractional residuals 4. Histogram of fractional residuals
- Parameters:
result (FluxCalibrationResult)
title (str)
figsize (tuple)
- Return type:
matplotlib.figure.Figure
- kspecdr.fluxcal.qc.plot_per_star_vectors(result: FluxCalibrationResult, log_scale: bool = False, figsize: tuple = (12, 5))[source]#
Plot each per-star calibration vector on its own axis.
- Parameters:
result (FluxCalibrationResult)
log_scale (bool)
figsize (tuple)
- Return type:
matplotlib.figure.Figure
- kspecdr.fluxcal.qc.plot_photometric_residuals(cal_vectors: List[CalibrationVector], figsize: tuple = (8, 5))[source]#
Bar chart of per-band photometric residuals (synth − obs) for each star.
- Parameters:
cal_vectors (list of CalibrationVector) – Each must have
meta["band_residuals"](fromscale_template_to_photometry()).figsize (tuple)
- Return type:
matplotlib.figure.Figure
- kspecdr.fluxcal.qc.plot_template_match(observed: Spectrum1D, template: Spectrum1D, rv_kms: float = 0.0, title: str = '', figsize: tuple = (12, 6))[source]#
Compare observed and best-matching template spectra.
Shows raw flux (top) and continuum-normalised (bottom, if continuum is available in template meta).
- Parameters:
observed (Spectrum1D) – Observed standard-star spectrum.
template (Spectrum1D) – Prepared (convolved + resampled) template spectrum.
rv_kms (float) – Measured RV in km/s (annotated on the plot).
title (str)
figsize (tuple)
- Return type:
matplotlib.figure.Figure
- kspecdr.fluxcal.qc.summarize_calibration(result: FluxCalibrationResult) str[source]#
Return a formatted text summary of the calibration result.
Suitable for printing in a notebook cell or logging.
- Parameters:
result (FluxCalibrationResult)
- Return type:
Download BOSZ 2024 stellar template subgrid for kspecdr flux calibration.
Usage#
From the repository root:
python -m kspecdr.fluxcal.download_bosz # full download
python -m kspecdr.fluxcal.download_bosz --dry-run # preview file list only
python -m kspecdr.fluxcal.download_bosz --force # re-download existing files
Subgrid definition (F-type spectrophotometric standards)#
Resolution: R = 10,000
Teff: 5000-8000 K, step 250 K (13 values)
log(g): 3.5, 4.0, 4.5, 5.0
[M/H]: -1.00 to +0.50, step 0.25 (7 values)
[alpha/M]: +0.00 (all [M/H]); additionally +0.25 for [M/H] <= -0.50
[C/M]: +0.00
vmicro: 1 km/s
atmos: mp (MARCS plane-parallel) for Teff 5000-7250 K, ap (ATLAS9 plane-parallel) for Teff 7500-8000 K
Output#
Files are saved to:
<repo>/data/templates/bosz2024/r10000/<metallicity>/<filename>.txt.gz
The shared wavelength grid is saved to:
<repo>/data/templates/bosz2024/bosz2024_wave_r10000.txt
References
Meszaros et al. 2024: The updated BOSZ synthetic stellar spectral library
MAST HLSP BOSZ: https://archive.stsci.edu/hlsp/bosz
- kspecdr.fluxcal.download_bosz.alpha_values_for(feh: float) list[float][source]#
Return [α/M] values to include for a given [M/H].
- kspecdr.fluxcal.download_bosz.atmos_for(teff: int) str[source]#
ATLAS9 plane-parallel for Teff ≥ 7500 K; MARCS plane-parallel otherwise.