Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
9 changes: 9 additions & 0 deletions CHANGES.rst
Original file line number Diff line number Diff line change
@@ -1,3 +1,12 @@
1.7.0 (unreleased)
------------------

New Features
^^^^^^^^^^^^

- Added a ``disp_bounds`` argument to ``tracing.FitTrace``. The argument allows for adjusting the
dispersion-axis window from which the trace peaks are estimated.

1.6.0 (2025-06-18)
------------------

Expand Down
4 changes: 2 additions & 2 deletions specreduce/tests/test_tracing.py
Original file line number Diff line number Diff line change
Expand Up @@ -460,7 +460,7 @@ def test_mask_treatment_apply(self, peak_method, expected):
# run FitTrace, with the testing-only flag _save_bin_peaks_testing set
# to True to return the bin peak values before fitting the trace
trace = FitTrace(imgg, peak_method=peak_method, _save_bin_peaks_testing=True)
x_bins, y_bins = trace._bin_peaks_testing
x_bins, y_bins = trace._bin_peaks
np.testing.assert_allclose(y_bins, expected, atol=0.1)

# check that final fit to all bins, accouting for fully-masked bins,
Expand Down Expand Up @@ -548,7 +548,7 @@ def test_mask_treatment_propagate(self, peak_method, expected):
mask_treatment="propagate",
_save_bin_peaks_testing=True,
)
x_bins, y_bins = trace._bin_peaks_testing
x_bins, y_bins = trace._bin_peaks
np.testing.assert_allclose(y_bins, expected)

# check that final fit to all bins, accouting for fully-masked bins,
Expand Down
22 changes: 15 additions & 7 deletions specreduce/tracing.py
Original file line number Diff line number Diff line change
Expand Up @@ -233,17 +233,22 @@ class FitTrace(Trace, _ImageParser):
into which to divide the image. If not set, defaults to one bin
per dispersion (wavelength) pixel in the given image. If set,
requires at least 4 or N bins for a degree N ``trace_model``,
whichever is greater. [default: None]
whichever is greater. [default: ``None``]
guess : int, optional
A guess at the trace's location in the cross-dispersion
(spatial) direction. If set, overrides the normal max peak
finder. Good for tracing a fainter source if multiple traces
are present. [default: None]
are present. [default: ``None``]
window : int, optional
Fit the trace to a region with size ``window * 2`` around the
guess position. Useful for tracing faint sources if multiple
traces are present, but potentially bad if the trace is
substantially bent or warped. [default: None]
substantially bent or warped. [default: ``None``]
disp_bounds
The lower and upper bounds of the pixel range along the dispersion
axis that is used when fitting the trace model. If ``None``,
defaults to the entire range of pixels along the dispersion axis.
[default: ``None``]
trace_model : one of `~astropy.modeling.polynomial.Chebyshev1D`,\
`~astropy.modeling.polynomial.Legendre1D`,\
`~astropy.modeling.polynomial.Polynomial1D`,\
Expand Down Expand Up @@ -276,14 +281,15 @@ class FitTrace(Trace, _ImageParser):
bins: int | None = None
guess: float | None = None
window: int | None = None
disp_bounds: tuple[float, float] | None = None
trace_model: Model = field(default=models.Polynomial1D(degree=1))
peak_method: Literal["gaussian", "centroid", "max"] = "max"
_crossdisp_axis: int = 0
_disp_axis: int = 1
mask_treatment: Literal["apply", "propagate", "apply_nan_only"] = "apply"
_valid_mask_treatment_methods = ("apply", "propagate", "apply_nan_only")
# for testing purposes only, save bin peaks if requested
_save_bin_peaks_testing: bool = False
_save_bin_peaks_testing: bool = True

def __post_init__(self):

Expand Down Expand Up @@ -322,7 +328,9 @@ def __post_init__(self):
"trace_model must be one of " f"{', '.join([m.name for m in valid_models])}."
)

cols = img.shape[self._disp_axis]
self.dlim = self.disp_bounds if self.disp_bounds is not None else (0, img.shape[
self._disp_axis])
cols = self.dlim[1]
model_deg = self.trace_model.degree
if self.bins is None:
self.bins = cols
Expand Down Expand Up @@ -394,7 +402,7 @@ def _fit_trace(self, img):
"for invalid values or use a larger window value."
)

x_bins = np.linspace(0, img.shape[self._disp_axis], self.bins + 1, dtype=int)
x_bins = np.linspace(self.dlim[0], self.dlim[1], self.bins + 1, dtype=int)
y_bins = np.tile(np.nan, self.bins)

warn_bins = []
Expand Down Expand Up @@ -485,7 +493,7 @@ def _fit_trace(self, img):

# for testing purposes only, save bin peaks if requested
if self._save_bin_peaks_testing:
self._bin_peaks_testing = (x_bins, y_bins)
self._bin_peaks = (x_bins, y_bins)

# filter non-finite bin peaks before filtering to all bin peaks
y_finite = np.where(np.isfinite(y_bins))[0]
Expand Down