stdpipe.photometry module

Routines for object detection and photometry.

stdpipe.photometry.dilate(image, mask)
stdpipe.photometry.make_kernel(r0=1.0, ext=1.0)[source]
stdpipe.photometry.estimate_fwhm(values, *, good=None, fwhm_range=(1.0, 20.0), min_candidates=5)[source]

Robust global FWHM estimate from per-object FWHM measurements.

Uses a sliding-window mode estimator (narrowest quartile interval) that is resistant to outliers from galaxies, blends, and non-Gaussian tails.

Parameters:
valuesarray-like

1-D array of per-object FWHM values in pixels.

goodarray-like of bool, optional

Boolean mask, True = use this object. The caller is responsible for applying any S/N, flag, or ellipticity cuts before calling, or encoding them in this mask.

fwhm_rangetuple of float

Hard clip range (lo, hi) in pixels. Values outside this range are excluded before the mode calculation.

min_candidatesint

Minimum number of surviving candidates required.

Returns:
float

Scalar FWHM estimate in pixels, or np.nan on failure.

class stdpipe.photometry.FWHMMap(coeffs, order, x0, y0, sx, sy, median, n_used, fwhm_range)[source]

Bases: object

Position-dependent FWHM estimator backed by a 2-D polynomial.

Produced by estimate_fwhm_from_objects() when spatial_order >= 1. Coordinates are normalised internally to roughly [-1, 1] for numerical stability.

Instances are callable: fwhm_map(x, y) returns FWHM in pixels at the requested pixel coordinates (x, y may be scalars or broadcastable arrays; the result has the broadcast shape). Instances are also convertible via float(fwhm_map) to the scalar median summary, so they can be dropped into APIs that still expect a single number.

Attributes:
coeffsndarray

Flat array of polynomial coefficients, ordered by (i, j) with i + j <= order and j varied fastest.

orderint

Polynomial order.

x0, y0, sx, syfloat

Coordinate normalisation: xn = (x - x0) / sx.

medianfloat

Representative scalar FWHM (median of the fitted surface over the good candidates). Used by float(self).

n_usedint

Number of candidates retained after sigma-clipping.

fwhm_rangetuple of float

Output clipping range in pixels.

Methods

__call__(x, y)

Call self as a function.

coeffs
order
x0
y0
sx
sy
median
n_used
fwhm_range
stdpipe.photometry.estimate_fwhm_from_objects(obj, *, snr_min=10.0, max_ellipticity=0.3, use_flags=True, fwhm_range=(1.0, 20.0), min_candidates=5, spatial_order=0, xy_cols=('x', 'y'), image_shape=None, verbose=False)[source]

Estimate global FWHM from an object table (SEP or SExtractor).

Extracts per-object FWHM values and builds a quality mask from the available catalog columns, then delegates to estimate_fwhm().

When a flux_radius column is present (half-light radius, as produced by get_objects_sep()), it is converted to an equivalent Gaussian FWHM via FWHM = 2 * flux_radius and used as the primary estimator. Half-light radii are 2–3× more stable than SEP’s Gaussian-core FWHM (relative IQR ~0.2 vs ~0.7) and far less susceptible to bias from sub-stellar contaminants, because the measurement integrates flux over many pixels rather than fitting a slope to a handful of core pixels.

When flux_radius is not available, the function falls back to the fwhm / FWHM_IMAGE column or moment-derived FWHM, and applies an a*b-weighted median correction if the mode appears biased (weighted median exceeds mode by > 10 %).

Parameters:
obj~astropy.table.Table or structured ndarray

Detection results. Recognized columns (checked in priority order): flux_radius (half-light radius), fwhm, FWHM_IMAGE, or derived from a/b (A_IMAGE/B_IMAGE) moments.

snr_minfloat or None

Minimum S/N for candidate selection (via magerr). Set to None to disable.

max_ellipticityfloat or None

Maximum ellipticity 1 - b/a. Set to None to disable.

use_flagsbool

If True, reject objects with nonzero flags.

fwhm_rangetuple of float

Passed through to estimate_fwhm().

min_candidatesint

Passed through to estimate_fwhm().

spatial_orderint

Polynomial order for a position-dependent FWHM model. 0 (default) preserves the historical behaviour and returns a scalar. >=1 attempts to fit a 2-D polynomial FWHM(x, y) to the per-object values and returns a FWHMMap callable. If too few candidates survive the quality cuts for the requested order, the function falls back to the scalar path.

xy_colstuple of str

Column names used for pixel coordinates when spatial_order >= 1. Defaults to ('x', 'y') (SEP-style). Use ('X_IMAGE', 'Y_IMAGE') for SExtractor output. Falls back to the scalar path if either column is missing.

image_shapetuple of int, optional

(H, W) image shape used to normalise pixel coordinates to roughly [-1, 1] when spatial_order >= 1. If not provided, the min/max of the candidate coordinates are used.

verbosebool or callable

If True, log diagnostic information about the estimation.

Returns:
float or FWHMMap

Scalar FWHM in pixels when spatial_order == 0. A FWHMMap callable otherwise, or a scalar on fallback. Returns np.nan if no usable candidates survive.

stdpipe.photometry.get_objects_sep(image, header=None, mask=None, mask_detect=None, err=None, thresh=4.0, aper=3.0, bkgann=None, r0=0.5, gain=1, edge=0, minnthresh=2, minarea=5, relfluxradius=2.0, wcs=None, bg_size=64, fwhm=False, optimal=False, group_sources=True, group_radius_factor=1.0, group_halo_factor=1.2, centroid=False, centroid_sigma_scale=1.0, centroid_maxstep_scale=1.0, centroid_maxshift_scale=0.5, use_mask_large=False, subtract_bg=True, npix_large=100, sn=10.0, get_segmentation=False, deblend_fwhm=0, deblend_method='watershed', clip_sigma=3.0, clip_iters=5, fwhm_spatial_order=0, verbose=True, **kwargs)[source]

Object detection and aperture/optimal photometry using SEP.

Signature is as similar as possible to get_objects_sextractor().

Algorithm: background estimation → object detection with optional smoothing and deblending (watershed by default) → optional windowed centroiding → edge rejection → per-object FWHM estimation → aperture or optimal extraction photometry → S/N quality cut.

Detection flags are documented at https://sep.readthedocs.io/en/v1.1.x/reference.html — they differ from SExtractor flags.

Parameters:
imagenumpy.ndarray

Input 2-D image.

headerastropy.io.fits.Header, optional

Image header.

masknumpy.ndarray of bool, optional

Pixel mask (True = masked). Used for background estimation and photometry. Objects whose detection footprint overlaps this mask are flagged with 0x100 but not excluded from detection.

mask_detectnumpy.ndarray of bool, optional

Detection mask (True = exclude from detection). Combined (OR) with non-finite pixel mask. Use this to exclude entire regions (e.g. bad columns, satellite trails). If None, only non-finite pixels are masked during detection.

errnumpy.ndarray, optional

Noise map.

threshfloat

Detection threshold in sigmas above local background.

aperfloat

Aperture radius in pixels, or in FWHM units when fwhm is provided/estimated.

bkganntuple of float, optional

Background annulus (r_in, r_out) in pixels (or FWHM units). Uses arithmetic mean of unmasked pixels inside the annulus. If None, the global background model is used.

r0float

Smoothing kernel sigma (pixels) for detection.

gainfloat

Detector gain, e-/ADU.

edgeint

Reject objects closer to image edge than this many pixels.

minnthreshint

Minimum number of pixels above threshold per detection.

minareaint

Minimum number of pixels in an object footprint.

relfluxradiusfloat

Multiplier for rough FWHM when computing flux_radius search radius (internal use).

wcsastropy.wcs.WCS, optional

WCS for sky coordinate assignment.

bg_sizeint

Background grid cell size in pixels.

fwhmbool or float

FWHM handling. False (default) — per-object FWHM for output only. True — estimate global FWHM and scale aper/bkgann (interpreted as FWHM units). Numeric — use this value directly and scale aper/bkgann.

optimalbool

Use optimal extraction (sep.sum_circle_optimal) instead of aperture photometry. Requires FWHM. SEP 1.4+ only.

group_sourcesbool

Grouped optimal extraction for overlapping sources. Default True (recommended). SEP 1.4+ only.

group_radius_factorfloat

Connectivity scale for grouped optimal extraction. The default 1.0 groups sources whose apertures overlap.

group_halo_factorfloat or None

Local context halo scale for grouped optimal extraction. The default 1.2 keeps a modest expanded solve context without widening source connectivity. None lets SEP use group_radius_factor.

centroidbool

Refine positions via sep.winpos. Default False (windowed positions can be biased in crowded fields). With SEP 1.4+, uses maxstep = 0.2 * FWHM.

centroid_sigma_scalefloat

Multiplicative scale applied to the Gaussian sigma passed to sep.winpos. The default 1.0 keeps the historical behavior sigma = FWHM / 2.355. Values below 1 use a tighter window; values above 1 use a broader window.

centroid_maxstep_scalefloat

Multiplicative scale applied to the SEP 1.4+ maxstep limit used by sep.winpos. The default 1.0 keeps the historical cap maxstep = 0.2 * FWHM. Values below 1 restrict the allowed centroid motion more strongly; values above 1 loosen it.

centroid_maxshift_scalefloat

Multiplicative scale applied to the SEP 1.4+ maxshift cap on the total centroid drift across all winpos iterations (in units of FWHM). The default 0.5 means a centroid can move at most 0.5 * FWHM from its initial position, complementing the per-iteration maxstep cap. Set to a large value to effectively disable. Ignored when SEP is older than 1.4 or when FWHM is not available.

use_mask_largebool

Filter out objects with footprints larger than npix_large.

npix_largeint

Pixel-count threshold for large-object rejection.

subtract_bgbool

Subtract background before detection (default True).

snfloat

Minimum S/N (magerr < 1/sn) for output.

get_segmentationbool

Also return the segmentation map and add seg_id column.

deblend_fwhmfloat

Fixed Gaussian FWHM (pixels) for deterministic deblending. 0 = adaptive shapes with stochastic assignment. SEP 1.4+ only.

deblend_method{‘watershed’, ‘threshold’}

Deblending algorithm. SEP 1.4+ only.

clip_sigmafloat

Sigma-clipping threshold for annulus background. SEP 1.4+ only.

clip_itersint

Max sigma-clipping iterations for annulus background. 0 = disable. SEP 1.4+ only.

fwhm_spatial_orderint

When FWHM is auto-estimated (fwhm=True or via the optimal path), fit a 2-D polynomial FWHM(x, y) of this order and apply it per source. 0 (default) keeps the historical scalar behaviour. Per-source values are used for aperture/background-annulus scaling and for the optimal-extraction PSF width; the windowed centroider still uses the scalar median. Ignored when a numeric fwhm is supplied directly. See FWHMMap.

verbosebool or callable

Print progress messages.

Returns:
astropy.table.Table or tuple

Table of detected objects. If get_segmentation is True, returns (table, segmentation_map).

Columns: x, y, xerr, yerr, flux, fluxerr, mag, magerr, flags, ra, dec, bg, fwhm, a, b, theta, and optionally seg_id.

Metadata: aper, bkgann, optimal, group_sources, fwhm_phot.

Notes

  • FWHM estimation prefers SEP’s built-in method (ported from SExtractor). Falls back to 2nd moments for older SEP versions.

  • Optimal extraction assumes a Gaussian PSF.

  • Grouped extraction fits nearby sources simultaneously — dramatically improves accuracy for close pairs.

  • Watershed deblending often gives better completeness for pairs closer than ~2 FWHM.

stdpipe.photometry.get_objects_sextractor(image, header=None, mask=None, mask_detect=None, err=None, thresh=2.0, aper=3.0, r0=0.0, gain=1, edge=0, minarea=5, wcs=None, sn=3.0, bg_size=None, sort=True, reject_negative=True, mask_to_nans=True, checkimages=[], extra_params=[], extra={}, psf=None, catfile=None, _workdir=None, _tmpdir=None, _exe=None, verbose=False)[source]

Thin wrapper around the SExtractor binary.

Processes the image with optional mask and noise map and returns the list of detected objects, optionally with SExtractor checkimages.

See the SExtractor documentation for details. Detection flags are documented at https://sextractor.readthedocs.io/en/latest/Flagging.html#extraction-flags-flags. Objects with masked pixels in their footprint additionally get flag 0x100.

Parameters:
imagenumpy.ndarray

Input 2-D image.

headerastropy.io.fits.Header, optional

Image header.

masknumpy.ndarray of bool, optional

Pixel mask (True = masked). Used for SExtractor flagging (FLAG_IMAGE / IMAFLAGS_ISO). Objects overlapping this mask get flag 0x100.

mask_detectnumpy.ndarray of bool, optional

Detection mask (True = exclude). When mask_to_nans is True, these pixels (plus non-finite pixels) are set to NaN before running SExtractor.

errnumpy.ndarray, optional

Noise map.

threshfloat

Detection threshold in sigmas (DETECT_THRESH).

aperfloat or list of float

Aperture radius in pixels. If a list, flux is measured for each aperture.

r0float

Smoothing kernel sigma (pixels, i.e. FWHM / 2.355).

gainfloat

Detector gain, e-/ADU.

edgeint

Reject objects closer to edge than this many pixels.

minareaint

Minimum object area (DETECT_MINAREA).

wcsastropy.wcs.WCS, optional

WCS for sky coordinate assignment.

snfloat

Minimum S/N for output.

bg_sizeint, optional

Background grid size (BACK_SIZE).

sortbool

Sort detections by decreasing brightness.

reject_negativebool

Reject objects with negative flux.

mask_to_nansbool

Replace detection-masked pixels with NaN before running SExtractor. Default True.

checkimageslist of str

SExtractor checkimage types to return (e.g. 'BACKGROUND', 'SEGMENTATION').

extra_paramslist of str

Extra SExtractor output parameters (see sex -dp).

extradict

Extra SExtractor configuration parameters (see sex -dd).

psfstr, optional

Path to PSFEx PSF model file for PSF photometry.

catfilestr, optional

Copy the output catalogue to this path.

_workdirstr, optional

Keep temporary files in this directory (for debugging).

_tmpdirstr, optional

Create temp directory inside this path.

_exestr, optional

Full path to SExtractor executable.

verbosebool or callable

Print progress messages.

Returns:
astropy.table.Table or list

Table of detected objects. If checkimages is non-empty, returns [table, checkimage1, checkimage2, ...].

class stdpipe.photometry.ConstantBackground(value, rms, shape)[source]

Bases: object

stdpipe.photometry.get_value(val)[source]

Get array from either Quantity or ndarray

stdpipe.photometry.get_objects_photutils(image, header=None, mask=None, mask_detect=None, err=None, thresh=2.0, method='segmentation', deblend=True, minarea=5, saturation=None, npixels=5, nlevels=32, contrast=0.001, connectivity=8, fwhm=3.0, sharplo=0.2, sharphi=1.0, roundlo=-1.0, roundhi=1.0, aper=3.0, bkgann=None, bg_size=64, subtract_bg=True, edge=0, sn=3.0, wcs=None, get_segmentation=False, verbose=False, **kwargs)[source]

Detect sources in image using photutils detection algorithms.

This function provides an alternative to get_objects_sep() using photutils instead of SEP. It supports multiple detection methods including segmentation-based detection with deblending and point source detection using DAOStarFinder or IRAFStarFinder.

Parameters:
imagenumpy.ndarray

2D image array

headerastropy.io.fits.Header, optional

FITS header for WCS information

masknumpy.ndarray, optional

Boolean mask (True = masked pixel). Used for background estimation, photometry, and flagging (0x100). Not used for detection itself.

mask_detectnumpy.ndarray, optional

Boolean mask for detection (True = exclude from detection). OR-ed with non-finite pixel mask. If None, only non-finite pixels are excluded from detection.

errnumpy.ndarray, optional

Error/uncertainty map. If None, uses background RMS

threshfloat, optional

Detection threshold in sigma units (default: 2.0)

methodstr, optional

Detection method: ‘segmentation’, ‘dao’, or ‘iraf’ (default: ‘segmentation’)

deblendbool, optional

Apply deblending for segmentation method (default: True)

minareaint, optional

Minimum source area in pixels (default: 5)

saturationfloat, optional

Saturation threshold for flagging saturated sources. Sources with peak pixel values exceeding this threshold will have flag 0x004 set. If None (default), saturation detection is disabled.

npixelsint, optional

Minimum connected pixels for segmentation (default: 5)

nlevelsint, optional

Number of deblending levels (default: 32)

contrastfloat, optional

Deblending contrast threshold (default: 0.001)

connectivityint, optional

Pixel connectivity (4 or 8) for segmentation (default: 8)

fwhmfloat, optional

FWHM for StarFinder methods in pixels (default: 3.0)

sharplofloat, optional

Lower sharpness bound for StarFinder (default: 0.2)

sharphifloat, optional

Upper sharpness bound for StarFinder (default: 1.0)

roundlofloat, optional

Lower roundness bound for StarFinder (default: -1.0)

roundhifloat, optional

Upper roundness bound for StarFinder (default: 1.0)

aperfloat, optional

Aperture radius in pixels (default: 3.0)

bkganntuple, optional

Background annulus (inner, outer) radii in pixels

bg_sizeint, optional

Background estimation box size (default: 64)

subtract_bgbool, optional

Subtract background before detection (default: True)

edgeint, optional

Edge exclusion in pixels (default: 0)

snfloat, optional

Minimum S/N ratio (default: 3.0)

wcsastropy.wcs.WCS, optional

WCS for coordinate conversion. If None, extracted from header

get_segmentationbool, optional

Return segmentation map (only for method=’segmentation’) (default: False)

verbosebool or callable, optional

Verbose output (default: False)

**kwargs

Additional parameters passed to photutils functions

Returns:
tableastropy.table.Table

Table of detected sources with columns: - x, y: pixel coordinates (0-indexed) - xerr, yerr: position errors - flux, fluxerr: aperture flux and error - mag, magerr: instrumental magnitude and error - fwhm: FWHM estimate - a, b: semi-major and semi-minor axes - theta: position angle - bg: local background - flags: detection flags - ra, dec: WCS coordinates (if WCS available)

segmphotutils.segmentation.SegmentationImage, optional

Segmentation map (only if get_segmentation=True and method=’segmentation’)

Notes

The function follows the same workflow as get_objects_sep(): 1. Validate inputs and create mask from NaNs 2. Estimate background using Background2D 3. Optionally subtract background 4. Detect sources using selected method 5. Measure aperture photometry 6. Apply quality filters 7. Add WCS coordinates if available 8. Return table with standard columns

Detection methods: - ‘segmentation’: Uses detect_sources() with optional deblending.

Best for extended sources and crowded fields.

  • ‘dao’: Uses DAOStarFinder for point sources. Best for stellar fields with well-defined PSF.

  • ‘iraf’: Uses IRAFStarFinder for IRAF-compatible detection. Similar to ‘dao’ but with IRAF-style parameters.

Examples

Basic segmentation detection:

>>> obj = get_objects_photutils(image, thresh=3.0, aper=5.0)

Segmentation with deblending:

>>> obj = get_objects_photutils(
...     image, method='segmentation', deblend=True,
...     nlevels=32, contrast=0.001
... )

Point source detection:

>>> obj = get_objects_photutils(
...     image, method='dao', fwhm=3.0, thresh=5.0
... )

With background annulus:

>>> obj = get_objects_photutils(
...     image, aper=5.0, bkgann=(10.0, 15.0)
... )

Get segmentation map:

>>> obj, segm = get_objects_photutils(
...     image, method='segmentation', get_segmentation=True
... )