Skip to content

Commit

Permalink
Add helper functions to control spectral shifting
Browse files Browse the repository at this point in the history
The setters for Spectrum1d.redshift and Spectrum1d.radial_velocity shift
the spectral axis, which from astropy#820 is not always what users expect. This
adds three helpers, one which does the current shifting (and is
explicitly documented to do that), and two which simply change the
values of redshift and radial_velocity.
  • Loading branch information
aragilar committed Apr 12, 2022
1 parent 26d3444 commit 08353d6
Show file tree
Hide file tree
Showing 2 changed files with 115 additions and 17 deletions.
80 changes: 67 additions & 13 deletions specutils/spectra/spectrum1d.py
Original file line number Diff line number Diff line change
Expand Up @@ -595,12 +595,6 @@ def redshift(self):
"""
return self.spectral_axis.redshift

@redshift.setter
def redshift(self, val):
new_spec_coord = self.spectral_axis.with_radial_velocity_shift(
-self.spectral_axis.radial_velocity).with_radial_velocity_shift(val)
self._spectral_axis = new_spec_coord

@property
def radial_velocity(self):
"""
Expand All @@ -613,15 +607,75 @@ def radial_velocity(self):
"""
return self.spectral_axis.radial_velocity

def set_redshift_to(self, redshift):
"""
This sets the redshift of the spectrum to be `redshift` *without*
changing the values of the `spectral_axis`.
If you want to shift the `spectral_axis` based on this value, use
`shift_spectrum_to`.
"""
new_spec_coord = self.spectral_axis.replicate(redshift=redshift)
self._spectral_axis = new_spec_coord

def set_radial_velocity_to(self, radial_velocity):
"""
This sets the radial velocity of the spectrum to be `radial_velocity`
*without* changing the values of the `spectral_axis`.
If you want to shift the `spectral_axis` based on this value, use
`shift_spectrum_to`.
"""
new_spec_coord = self.spectral_axis.replicate(
radial_velocity=radial_velocity
)
self._spectral_axis = new_spec_coord

def shift_spectrum_to(self, *, redshift=None, radial_velocity=None):
"""
This shifts in-place the values of the `spectral_axis`, given either a
redshift or radial velocity.
If you do *not* want to change the `spectral_axis`, use
`set_redshift_to` or `set_radial_velocity_to`.
"""
if redshift is not None and radial_velocity is not None:
raise ValueError(
"Only one of redshift or radial_velocity can be used."
)
if redshift is not None:
new_spec_coord = self.spectral_axis.with_radial_velocity_shift(
-self.spectral_axis.radial_velocity
).with_radial_velocity_shift(redshift)
self._spectral_axis = new_spec_coord
elif radial_velocity is not None:
if radial_velocity is not None:
if not radial_velocity.unit.is_equivalent(u.km/u.s):
raise u.UnitsError("Radial velocity must be a velocity.")

new_spectral_axis = self.spectral_axis.with_radial_velocity_shift(
-self.spectral_axis.radial_velocity
).with_radial_velocity_shift(radial_velocity)
self._spectral_axis = new_spectral_axis
else:
raise ValueError("One of redshift or radial_velocity must be set.")

@redshift.setter
def redshift(self, val):
warnings.warn(
"Setting the redshift of a spectrum is ambiguous, use either "
"set_redshift_to or shift_spectrum_to to be explicit."
)
self.shift_spectrum_to(redshift=val)

@radial_velocity.setter
def radial_velocity(self, val):
if val is not None:
if not val.unit.is_equivalent(u.km/u.s):
raise u.UnitsError("Radial velocity must be a velocity.")

new_spectral_axis = self.spectral_axis.with_radial_velocity_shift(
-self.spectral_axis.radial_velocity).with_radial_velocity_shift(val)
self._spectral_axis = new_spectral_axis
warnings.warn(
"Setting the radial velocity of a spectrum is ambiguous, use "
"either set_radial_velocity_to or shift_spectrum_to to be "
"explicit."
)
self.shift_spectrum_to(radial_velocity=val)

def __add__(self, other):
if not isinstance(other, NDCube):
Expand Down
52 changes: 48 additions & 4 deletions specutils/tests/test_spectral_axis.py
Original file line number Diff line number Diff line change
Expand Up @@ -130,7 +130,11 @@ def test_change_radial_velocity():

assert spec.radial_velocity == 0 * u.km/u.s

spec.radial_velocity = 1 * u.km / u.s
with pytest.warns(
UserWarning,
match="Setting the radial velocity of a spectrum is ambiguous"
):
spec.radial_velocity = 1 * u.km / u.s

assert spec.radial_velocity == 1 * u.km/u.s

Expand All @@ -139,11 +143,27 @@ def test_change_radial_velocity():

assert spec.radial_velocity == 10 * u.km / u.s

spec.radial_velocity = 5 * u.km / u.s
with pytest.warns(
UserWarning,
match="Setting the radial velocity of a spectrum is ambiguous"
):
spec.radial_velocity = 5 * u.km / u.s

assert spec.radial_velocity == 5 * u.km / u.s


def test_no_change_radial_velocity():
wave = np.linspace(100, 200, 100) * u.AA
flux = np.ones(100) * u.one
spec = Spectrum1D(spectral_axis=wave, flux=flux,
radial_velocity=0 * u.km / u.s)

assert spec.radial_velocity == 0 * u.km/u.s
spec.set_radial_velocity_to(10 * u.km/u.s)
assert spec.radial_velocity == 10 * u.km/u.s
assert_quantity_allclose(spec.wavelength, wave)


def test_change_redshift():
wave = np.linspace(100, 200, 100) * u.AA
flux = np.ones(100) * u.one
Expand All @@ -153,7 +173,10 @@ def test_change_redshift():
assert_quantity_allclose(spec.redshift, u.Quantity(0))
assert type(spec.spectral_axis) == SpectralAxis

spec.redshift = 0.1
with pytest.warns(
UserWarning, match="Setting the redshift of a spectrum is ambiguous"
):
spec.redshift = 0.1

assert spec.redshift.unit.physical_type == 'dimensionless'
assert_quantity_allclose(spec.redshift, u.Quantity(0.1))
Expand All @@ -165,8 +188,29 @@ def test_change_redshift():
assert_quantity_allclose(spec.redshift, u.Quantity(0.2))
assert type(spec.spectral_axis) == SpectralAxis

spec.redshift = 0.4
with pytest.warns(
UserWarning, match="Setting the redshift of a spectrum is ambiguous"
):
spec.redshift = 0.4

assert spec.redshift.unit.physical_type == 'dimensionless'
assert_quantity_allclose(spec.redshift, u.Quantity(0.4))
assert type(spec.spectral_axis) == SpectralAxis


def test_no_change_redshift():
wave = np.linspace(100, 200, 100) * u.AA
flux = np.ones(100) * u.one
spec = Spectrum1D(spectral_axis=wave, flux=flux, redshift=0)

assert spec.redshift.unit.physical_type == 'dimensionless'
assert_quantity_allclose(spec.redshift, u.Quantity(0))
assert type(spec.spectral_axis) == SpectralAxis

spec.set_redshift_to(0.5)

assert spec.redshift.unit.physical_type == 'dimensionless'
assert_quantity_allclose(spec.redshift, u.Quantity(0.5))
assert type(spec.spectral_axis) == SpectralAxis

assert_quantity_allclose(spec.wavelength, wave)

0 comments on commit 08353d6

Please sign in to comment.