Source code for tengri.observation.filters

"""Photometric filter management via the SVO Filter Profile Service.

Downloads, caches, and loads photometric filter transmission curves from
the Spanish Virtual Observatory (SVO) Filter Profile Service:
https://svo2.cab.inta-csic.es/theory/fps/

Uses astroquery.svo_fps for downloads. Filters are cached as two-column
text files (wavelength in Angstrom, transmission) under a configurable
cache directory.

Note on ALMA / interferometers
------------------------------
ALMA and similar interferometric arrays do not use photometric bandpass
filters — observations are defined by spectral windows in GHz. SVO has
no ALMA entries. For SED fitting at (sub)mm continuum frequencies, use
``load_tophat_filter()`` to create a synthetic rectangular bandpass
centered on the observed frequency.

Available submm photometric instruments (real bandpasses on SVO):
  JCMT SCUBA-2 : 450 μm, 850 μm  → scuba2_450, scuba2_850
  APEX LABOCA  : 870 μm           → laboca_870
  APEX SABOCA  : 350 μm           → saboca_350
"""

import json
from pathlib import Path

import jax.numpy as jnp
import numpy as np

# numpy 2.0 removed np.trapz; numpy >= 1.26 provides np.trapezoid.
try:
    from numpy import trapezoid as _np_trapezoid
except ImportError:  # numpy < 1.26
    from numpy import trapz as _np_trapezoid  # type: ignore[no-redef]

from tengri.observation.photometry import FilterCurve

# ── Registry: short name -> SVO Filter Profile Service ID ─────────

_REGISTRY_PATH = Path(__file__).parent.parent / "data" / "filters_registry.json"
FILTER_REGISTRY: dict[str, str] = json.loads(_REGISTRY_PATH.read_text())

_DEFAULT_CACHE_DIR = "data/filters"

# Speed of light in Å/s — used for GHz ↔ Å conversion.
_C_AA_S = 2.99792458e18

# ALMA receiver band definitions (ALMA Cycle 11 specifications).
# Each entry maps band number → (lo_ghz, hi_ghz) at the edges of the
# receiver bandwidth.  The full band width is used as the top-hat width
# so that continuum photometry integrates over the realistic frequency
# coverage rather than an arbitrarily narrow window.
_ALMA_BANDS_GHZ: dict[int, tuple[float, float]] = {
    1: (35.0, 50.0),
    2: (67.0, 90.0),
    3: (84.0, 116.0),
    4: (125.0, 163.0),
    5: (163.0, 211.0),
    6: (211.0, 275.0),
    7: (275.0, 373.0),
    8: (385.0, 500.0),
    9: (602.0, 720.0),
    10: (787.0, 950.0),
}


# ── Filter metadata: facility and description for rich listing ────

_FACILITY_FROM_PREFIX: dict[str, str] = {
    "sdss": "SDSS",
    "lsst": "LSST/Rubin",
    "ps1": "Pan-STARRS",
    "des": "DES/DECam",
    "megacam": "CFHT/MegaCam",
    "hsc": "Subaru/HSC",
    "suprime": "Subaru/SuprimeCam",
    "galex": "GALEX",
    "xmm": "XMM-Newton/OM",
    "uvot": "Swift/UVOT",
    "2mass": "2MASS",
    "vista": "VISTA/VIRCAM",
    "ukidss": "UKIRT/WFCAM",
    "hst": "HST",
    "jwst": "JWST/NIRCam",
    "nircam25": "JWST/NIRCam2025",
    "niriss": "JWST/NIRISS",
    "nirspec": "JWST/NIRSpec",
    "miri": "JWST/MIRI",
    "roman": "Roman/WFI",
    "euclid": "Euclid",
    "irac": "Spitzer/IRAC",
    "mips": "Spitzer/MIPS",
    "wise": "WISE",
    "akari": "AKARI",
    "herschel": "Herschel",
    "scuba2": "JCMT/SCUBA-2",
    "laboca": "APEX/LABOCA",
    "saboca": "APEX/SABOCA",
    "johnson": "Generic/Johnson",
    "cousins": "Generic/Cousins",
}


def _infer_facility(name: str) -> str:
    """Infer facility from filter short name prefix."""
    for prefix, facility in _FACILITY_FROM_PREFIX.items():
        if name.startswith(prefix):
            return facility
    return "Other"


# ── Filter property computation (pure numpy, no JAX) ──────────────


def compute_effective_wavelength(wave: np.ndarray, trans: np.ndarray) -> float:
    """Photon-counting effective wavelength (pivot wavelength): λ_eff = ∫T·λ·dλ / ∫T·dλ.

    Parameters
    ----------
    wave : array, shape (n_wave,)
        Wavelength [Angstrom].
    trans : array, shape (n_wave,)
        Transmission (dimensionless [0, 1]).

    Returns
    -------
    float
        Photon-counting effective wavelength (pivot wavelength) [Angstrom].

    Notes
    -----
    This function computes the **standard astronomical photon-counting**
    effective wavelength, also known as the **pivot wavelength**. It is
    used for filter metadata, observational work, and validation against
    external codes (e.g., FSPS).

    NOT the same as the power-weighted effective wavelength used in
    approximate photometry. See :func:`~tengri.utils.optimizations.effective_wavelength`
    for the alternative definition.

    Not JAX-compatible (uses NumPy). Intended for filter metadata
    computation, not forward model evaluation.

    """
    num = _np_trapezoid(trans * wave, wave)
    den = _np_trapezoid(trans, wave)
    if den == 0:
        return 0.0
    return float(num / den)


def compute_fwhm(wave: np.ndarray, trans: np.ndarray) -> float:
    """Full width at half maximum of the transmission curve.

    Parameters
    ----------
    wave : array, shape (n_wave,)
        Wavelength [Angstrom].
    trans : array, shape (n_wave,)
        Transmission (dimensionless [0, 1]).

    Returns
    -------
    float
        FWHM [Angstrom]. Returns 0 if the curve never exceeds half-max.

    Notes
    -----
    Not JAX-compatible (uses NumPy). Intended for filter metadata
    computation, not forward model evaluation.

    """
    peak = np.max(trans)
    if peak == 0:
        return 0.0
    half_max = peak / 2.0
    above = wave[trans >= half_max]
    if len(above) < 2:
        return 0.0
    return float(above[-1] - above[0])


def _format_wavelength(wave_aa: float) -> str:
    """Format wavelength with appropriate units."""
    if wave_aa >= 1e7:
        return f"{wave_aa / 1e8:.2f} cm"
    elif wave_aa >= 1e4:
        return f"{wave_aa / 1e4:.2f} \u03bcm"
    else:
        return f"{wave_aa:.0f} \u00c5"


def filter_info(name: str, *, cache_dir: str | None = None) -> dict:
    """Return metadata for a single filter.

    Loads the transmission curve from the local cache (downloading from
    SVO if needed) and computes derived properties.

    Parameters
    ----------
    name : str
        Short filter name from ``FILTER_REGISTRY``.
    cache_dir : str, optional
        Override cache directory.

    Returns
    -------
    dict
        Keys: ``name``, ``svo_id``, ``facility``, ``lambda_eff_aa``,
        ``fwhm_aa``, ``lambda_eff_str``, ``fwhm_str``.

    Raises
    ------
    KeyError
        If *name* is not in the registry.

    Notes
    -----
    Not JAX-compatible (uses NumPy and file I/O). Intended for filter
    metadata computation and interactive exploration of the registry,
    not for forward model evaluation.

    """
    if name not in FILTER_REGISTRY:
        raise KeyError(f"Unknown filter '{name}'. Use list_available_filters() to see options.")
    kwargs = {"cache_dir": cache_dir} if cache_dir is not None else {}
    fc = load_filter(name, **kwargs)
    wave_np = np.asarray(fc.wave)
    trans_np = np.asarray(fc.trans)
    lam_eff = compute_effective_wavelength(wave_np, trans_np)
    fwhm = compute_fwhm(wave_np, trans_np)
    return {
        "name": name,
        "svo_id": FILTER_REGISTRY[name],
        "facility": _infer_facility(name),
        "lambda_eff_aa": lam_eff,
        "fwhm_aa": fwhm,
        "lambda_eff_str": _format_wavelength(lam_eff),
        "fwhm_str": _format_wavelength(fwhm),
    }


# ── Internal helpers ──────────────────────────────────────────────


def _svo_id_to_filename(svo_id: str) -> str:
    """Convert SVO filter ID to a safe filename."""
    return svo_id.replace("/", "_").replace(".", "_") + ".dat"


def _save_filter(filepath: Path, wave: np.ndarray, trans: np.ndarray) -> None:
    """Write wavelength and transmission columns to a two-column text file."""
    header = "# Wavelength(Angstrom)  Transmission"
    np.savetxt(str(filepath), np.column_stack([wave, trans]), header=header, fmt="%.6e")


def _load_filter_file(filepath: Path) -> tuple[np.ndarray, np.ndarray]:
    """Read and return wavelength and transmission columns from a text file."""
    data = np.loadtxt(str(filepath))
    if data.ndim != 2 or data.shape[1] < 2:
        raise ValueError(
            f"Filter file {filepath} must have at least 2 columns "
            f"(wavelength, transmission). Got shape {data.shape}."
        )
    return data[:, 0], data[:, 1]


def _fetch_from_svo(svo_id: str) -> tuple[np.ndarray, np.ndarray]:
    """Fetch filter curve from SVO using astroquery.svo_fps."""
    try:
        from astroquery.svo_fps import SvoFps
    except ImportError as exc:
        raise ImportError(
            "astroquery is required to download filters from SVO. "
            "Install it with:  pip install astroquery"
        ) from exc

    table = SvoFps.get_transmission_data(svo_id)
    wave = np.asarray(table["Wavelength"], dtype=float)
    trans = np.asarray(table["Transmission"], dtype=float)
    if len(wave) == 0:
        raise ValueError(
            f"SVO returned zero rows for filter '{svo_id}'. Check that the filter ID is correct."
        )
    return wave, trans


# ── Public API ────────────────────────────────────────────────────


def download_filter(
    svo_id: str,
    cache_dir: str = _DEFAULT_CACHE_DIR,
) -> tuple[np.ndarray, np.ndarray]:
    """Download a single filter from the SVO Filter Profile Service.

    If the filter is already cached on disk, it is loaded from the cache
    instead of re-downloading.

    Parameters
    ----------
    svo_id : str
        SVO filter identifier (e.g. ``"JWST/NIRCam.F200W"``).
    cache_dir : str
        Directory for cached filter files.

    Returns
    -------
    wave : ndarray
        Wavelength in Angstrom.
    trans : ndarray
        Transmission (dimensionless).

    Raises
    ------
    ImportError
        If ``astroquery`` is not installed.
    ValueError
        If SVO returns no data for the given filter ID.

    Notes
    -----
    Not JAX-compatible (uses file I/O and astroquery). Caching avoids
    redundant SVO downloads. The returned transmission values are not
    normalized — the absolute scale cancels in the photometry integral
    ``∫fλTλdλ / ∫Tλdλ``.

    """
    cache_path = Path(cache_dir)
    cache_path.mkdir(parents=True, exist_ok=True)
    filepath = cache_path / _svo_id_to_filename(svo_id)

    if filepath.exists():
        return _load_filter_file(filepath)

    wave, trans = _fetch_from_svo(svo_id)

    order = np.argsort(wave)
    wave = wave[order]
    trans = trans[order]

    _save_filter(filepath, wave, trans)
    return wave, trans


def load_filter(
    name: str,
    cache_dir: str = _DEFAULT_CACHE_DIR,
) -> FilterCurve:
    """Load a filter by its short registry name.

    Downloads from SVO if not already cached.

    Parameters
    ----------
    name : str
        Short name from ``FILTER_REGISTRY`` (e.g. ``"jwst_f200w"``).
    cache_dir : str
        Directory for cached filter files.

    Returns
    -------
    FilterCurve
        Filter with wavelength (Angstrom), raw transmission as returned
        by SVO, and name.  Transmission values are not normalized — the
        absolute scale cancels in the photometry integral
        ``∫fλTλdλ / ∫Tλdλ``.

    Raises
    ------
    KeyError
        If *name* is not in ``FILTER_REGISTRY``.

    Notes
    -----
    Not JAX-compatible (uses file I/O and astroquery for downloads).
    Preferred interface over :func:`download_filter` for user code;
    provides a JAX array return type (FilterCurve) rather than raw NumPy.

    """
    if name not in FILTER_REGISTRY:
        raise KeyError(
            f"Unknown filter '{name}'. Use list_available_filters() to see "
            f"valid names, or use load_custom_filter() for arbitrary files."
        )

    svo_id = FILTER_REGISTRY[name]
    wave, trans = download_filter(svo_id, cache_dir=cache_dir)
    return FilterCurve(wave=jnp.array(wave), trans=jnp.array(trans), name=name)


[docs] def load_filter_set( names: list[str], cache_dir: str = _DEFAULT_CACHE_DIR, ) -> tuple[list[jnp.ndarray], list[jnp.ndarray], list[FilterCurve]]: """Load multiple filters by short name. Parameters ---------- names : list of str Short names from ``FILTER_REGISTRY``. cache_dir : str Directory for cached filter files. Default: ``"data/filters"``. Returns ------- filter_waves : list of jnp.ndarray Wavelength arrays per filter, each shape ``(n_wave,)`` [Angstrom]. filter_trans : list of jnp.ndarray Transmission arrays per filter, each shape ``(n_wave,)`` (dimensionless [0, 1]). filter_curves : list of FilterCurve Full FilterCurve objects with wavelength, transmission, and name. Raises ------ KeyError If any name is not in ``FILTER_REGISTRY``. Notes ----- Filters are downloaded from SVO on first use and cached locally. See ``load_filter()`` for single-filter loading. Examples -------- >>> from tengri import load_filter_set >>> waves, trans, curves = load_filter_set(["sdss_r", "sdss_i", "sdss_z"]) >>> len(curves) 3 >>> curves[0].name 'sdss_r' """ filter_waves: list[jnp.ndarray] = [] filter_trans: list[jnp.ndarray] = [] filter_curves: list[FilterCurve] = [] for name in names: fc = load_filter(name, cache_dir=cache_dir) filter_waves.append(fc.wave) filter_trans.append(fc.trans) filter_curves.append(fc) return filter_waves, filter_trans, filter_curves
def load_custom_filter(filepath: str) -> FilterCurve: """Load a custom filter from a two-column text file. Parameters ---------- filepath : str Path to a text file with columns: wavelength (Angstrom), transmission. Returns ------- FilterCurve Filter with raw transmission values (not normalized). The absolute scale cancels in the photometry integral ``∫fλTλdλ / ∫Tλdλ``. Raises ------ FileNotFoundError If *filepath* does not exist. ValueError If the file format is invalid. Notes ----- Not JAX-compatible (uses file I/O). Useful when filter profiles are not available in SVO or when using synthetic, custom-defined bandpasses. File format is flexible: trailing whitespace and comment lines (starting with #) are ignored. """ path = Path(filepath) if not path.exists(): raise FileNotFoundError(f"Filter file not found: {filepath}") wave, trans = _load_filter_file(path) order = np.argsort(wave) wave, trans = wave[order], trans[order] return FilterCurve(wave=jnp.array(wave), trans=jnp.array(trans), name=path.stem) def load_tophat_filter( wave_center_aa: float, width_aa: float, name: str = "", n_points: int = 50, ) -> FilterCurve: """Create a synthetic top-hat filter (e.g. for ALMA continuum bands). Use this when the photometric measurement does not correspond to a standard bandpass on SVO — for example, an ALMA continuum flux at a given observed frequency. Parameters ---------- wave_center_aa : float Central wavelength [Angstrom]. width_aa : float Full width of the top-hat [Angstrom]. name : str Label for this filter (e.g. ``"alma_band6"``). Default: empty string. n_points : int Number of wavelength samples. Default: 50. Returns ------- FilterCurve Rectangular transmission curve with uniform transmission = 1.0. Notes ----- Useful for continuum photometry measurements (e.g., ALMA, SCUBA-2) that are defined in frequency space rather than bandpass shape. """ wave = jnp.linspace(wave_center_aa - width_aa / 2, wave_center_aa + width_aa / 2, n_points) trans = jnp.ones(n_points) return FilterCurve(wave=wave, trans=trans, name=name) def load_alma_band(band: int, name: str | None = None) -> FilterCurve: """Create a synthetic top-hat filter for an ALMA continuum band. ALMA is an interferometric array with no entries on the SVO Filter Profile Service. This function constructs a rectangular bandpass spanning the full receiver bandwidth of the requested band, which is appropriate for fitting SED continuum photometry. Parameters ---------- band : int ALMA band number (1–10). name : str, optional Label for the filter. Defaults to ``"alma_band{N}"``. Returns ------- FilterCurve Top-hat bandpass in observed-frame wavelengths (Angstrom). Examples -------- >>> fc = load_alma_band(6) # 1.23 mm continuum (211–275 GHz) >>> fc = load_alma_band(7) # 870 μm continuum (275–373 GHz) Notes ----- Band definitions follow the ALMA Cycle 11 receiver specifications. Wavelengths are in the *observed* frame — the filter should be applied at the observed frequency. For a source at redshift *z*, Band N probes rest-frame wavelength λ_rest = λ_obs / (1 + z). """ if band not in _ALMA_BANDS_GHZ: valid = sorted(_ALMA_BANDS_GHZ) raise ValueError(f"ALMA band must be one of {valid}, got {band}.") lo_ghz, hi_ghz = _ALMA_BANDS_GHZ[band] # High frequency = short wavelength and vice versa. lo_aa = _C_AA_S / (hi_ghz * 1e9) hi_aa = _C_AA_S / (lo_ghz * 1e9) center_aa = (lo_aa + hi_aa) / 2.0 width_aa = hi_aa - lo_aa label = name if name is not None else f"alma_band{band}" return load_tophat_filter(center_aa, width_aa, name=label) def list_available_filters( *, group_by: str = "facility", compute_properties: bool = False, cache_dir: str | None = None, ) -> dict[str, str]: """Print and return the filter registry, optionally grouped by facility. Parameters ---------- group_by : str Grouping key. ``"facility"`` (default) groups by telescope/instrument. ``"none"`` lists filters alphabetically without grouping. Default: ``"facility"``. compute_properties : bool If ``True``, load each filter's transmission curve and display effective wavelength and FWHM columns. This triggers SVO downloads for any filters not yet cached. Default: ``False``. cache_dir : str, optional Override cache directory for filter downloads. Default: ``None``. Returns ------- dict Copy of ``FILTER_REGISTRY`` (short name -> SVO ID). Also prints the registry to stdout in human-readable format. Notes ----- This function is primarily for interactive exploration of available filters. For programmatic use, access ``FILTER_REGISTRY`` directly. """ if group_by == "none": _print_flat_listing(compute_properties, cache_dir) else: _print_grouped_listing(compute_properties, cache_dir) return dict(FILTER_REGISTRY) def _print_flat_listing(compute_properties: bool, cache_dir: str | None) -> None: """Print filter registry as a flat (non-grouped) alphabetical list.""" if compute_properties: hdr = f"{'Name':<22s} {'SVO ID':<35s} {'lambda_eff':>12s} {'FWHM':>10s}" print(hdr) print("-" * len(hdr)) for name in sorted(FILTER_REGISTRY): kwargs = {"cache_dir": cache_dir} if cache_dir is not None else {} info = filter_info(name, **kwargs) print( f"{name:<22s} {info['svo_id']:<35s} " f"{info['lambda_eff_str']:>12s} {info['fwhm_str']:>10s}" ) else: hdr = f"{'Name':<22s} {'SVO ID':<35s}" print(hdr) print("-" * len(hdr)) for name, svo_id in sorted(FILTER_REGISTRY.items()): print(f"{name:<22s} {svo_id:<35s}") print(f"\nTotal: {len(FILTER_REGISTRY)} filters") def _print_grouped_listing(compute_properties: bool, cache_dir: str | None) -> None: """Print filter registry grouped by facility/telescope.""" groups: dict[str, list[str]] = {} for name in sorted(FILTER_REGISTRY): fac = _infer_facility(name) groups.setdefault(fac, []).append(name) for fac in sorted(groups): names = groups[fac] print(f"\n{'=' * 60}") print(f" {fac} ({len(names)} filters)") print(f"{'=' * 60}") if compute_properties: hdr = f" {'Name':<22s} {'lambda_eff':>12s} {'FWHM':>10s}" print(hdr) print(f" {'-' * (len(hdr) - 2)}") kwargs = {"cache_dir": cache_dir} if cache_dir is not None else {} for name in names: info = filter_info(name, **kwargs) print(f" {name:<22s} {info['lambda_eff_str']:>12s} {info['fwhm_str']:>10s}") else: for name in names: print(f" {name:<22s} {FILTER_REGISTRY[name]}") print(f"\nTotal: {len(FILTER_REGISTRY)} filters across {len(groups)} facilities") # ── User-facing filter discovery helpers ────────────────────────── # Thin convenience functions for discoverability on top of load_filter_set def list_filters(instrument: str | None = None) -> list[str]: """Return sorted list of filter names shipped with tengri. Parameters ---------- instrument : str, optional Filter by instrument prefix (case-insensitive, substring match). E.g., "sdss", "jwst", "hst". Default: None (return all). Returns ------- list of str Sorted filter names. Empty list if no matches. Notes ----- Filtering is permissive: matches if the instrument string appears anywhere in the filter name (case-insensitive). Examples -------- >>> all_filters = list_filters() >>> sdss_filters = list_filters(instrument="sdss") >>> len(sdss_filters) 5 """ names = sorted(FILTER_REGISTRY.keys()) if instrument is None: return names instrument_lower = instrument.lower() return [name for name in names if instrument_lower in name.lower()] def load(names: list[str]): """Load multiple filters by short name. Thin alias for ``load_filter_set`` for discoverability and consistency with the filters namespace. Parameters ---------- names : list of str Short filter names from the registry (e.g., ["sdss_r", "jwst_f200w"]). Returns ------- filter_waves : list of jnp.ndarray Wavelength arrays per filter, shape (n_wave,) [Angstrom]. filter_trans : list of jnp.ndarray Transmission arrays per filter, shape (n_wave,) (dimensionless [0, 1]). filter_curves : list of FilterCurve Full FilterCurve objects with wavelength, transmission, and name. Raises ------ KeyError If any name is not in the filter registry. Examples -------- >>> waves, trans, curves = load(["sdss_r", "sdss_i"]) >>> len(curves) 2 >>> curves[0].name 'sdss_r' Notes ----- Filters are downloaded from the SVO Filter Profile Service on first use and cached locally under data/filters/. """ return load_filter_set(names) def describe(name: str) -> str: """Return a one-line description of a filter. Computes the transmission-weighted effective wavelength and wavelength range of the transmission curve. Parameters ---------- name : str Filter short name from the registry. Returns ------- str Human-readable description. Format: "<name>: lambda_eff ~ X.XXX μm (range A–B μm)" or similar on failure. Notes ----- If the filter fails to load, returns a fallback description. Effective wavelength is computed as the transmission-weighted mean. """ try: fc = load_filter_set([name])[2][0] wave_np = np.asarray(fc.wave) trans_np = np.asarray(fc.trans) # Compute transmission-weighted effective wavelength lam_eff = compute_effective_wavelength(wave_np, trans_np) # Compute range (wavelengths with nonzero transmission) nonzero = trans_np > 0 if np.any(nonzero): wave_min = wave_np[nonzero].min() wave_max = wave_np[nonzero].max() else: wave_min, wave_max = wave_np.min(), wave_np.max() # Format wavelengths def format_wave(w): if w >= 1e4: return f"{w / 1e4:.2f}" else: return f"{w:.0f}" if lam_eff >= 1e4: unit = "μm" lam_eff_fmt = f"{lam_eff / 1e4:.3f}" else: unit = "Å" lam_eff_fmt = f"{lam_eff:.0f}" if wave_min >= 1e4: min_fmt = f"{wave_min / 1e4:.2f}" max_fmt = f"{wave_max / 1e4:.2f}" range_unit = "μm" else: min_fmt = f"{wave_min:.0f}" max_fmt = f"{wave_max:.0f}" range_unit = "Å" return f"{name}: λ_eff ~ {lam_eff_fmt} {unit} (range {min_fmt}{max_fmt} {range_unit})" except Exception: return f"{name}: (filter found; no summary available)" def suggest( redshift: float, coverage: str = "visible_to_nir", ) -> list[str]: """Suggest filters covering a rest-frame wavelength range at a redshift. Parameters ---------- redshift : float Redshift of the source (z >= 0). coverage : str Rest-frame wavelength coverage preset. Options: - "visible": 3500–9000 Å (optical) - "visible_to_nir": 3500–25000 Å (optical + near-IR) [default] - "uv_to_ir": 1200–50000 Å (UV + optical + IR) - "jwst_cover": 6000–50000 Å (rest-frame for JWST epochs) Returns ------- list of str Filter names with effective wavelength falling within the observed-frame span corresponding to the rest-frame coverage. Sorted by effective wavelength. Raises ------ ValueError If coverage is not recognized. Notes ----- Observed-frame wavelength is computed as: λ_obs = λ_rest * (1 + z). Examples -------- >>> suggest(z=3.0, coverage="visible_to_nir") # z=3 galaxies, optical→NIR ['jwst_f115w', 'jwst_f150w', ...] """ # Coverage presets (rest-frame, Angstrom) coverage_map = { "visible": (3500, 9000), "visible_to_nir": (3500, 25000), "uv_to_ir": (1200, 50000), "jwst_cover": (6000, 50000), } if coverage not in coverage_map: raise ValueError( f"Unknown coverage '{coverage}'. Must be one of {list(coverage_map.keys())}." ) lam_rest_min, lam_rest_max = coverage_map[coverage] # Convert to observed frame lam_obs_min = lam_rest_min * (1 + redshift) lam_obs_max = lam_rest_max * (1 + redshift) # Load all filters and compute effective wavelengths all_names = list_filters() if not all_names: return [] # Load all filters; skip any that fail wavelengths_by_name = {} try: for name in all_names: try: fc = load_filter_set([name])[2][0] wave_np = np.asarray(fc.wave) trans_np = np.asarray(fc.trans) lam_eff = compute_effective_wavelength(wave_np, trans_np) wavelengths_by_name[name] = lam_eff except Exception: # Skip filters that fail to load continue except Exception: return [] # Find filters within observed-frame span matches = [ name for name, lam_eff in wavelengths_by_name.items() if lam_obs_min <= lam_eff <= lam_obs_max ] # Sort by effective wavelength matches.sort(key=lambda name: wavelengths_by_name[name]) return matches