Skip to content

Commit

Permalink
Provide correct slope for color ratios with `SpectralGradient.from_co…
Browse files Browse the repository at this point in the history
…lor` (#392)

* Provide correct slope for color ratios with SpectralGradient.from_color

* Adapt to astropy 6 behavior; light reformatting.

u.Magnitudes (e.g., from the difference of two VEGAmags) do not
support the unary -.  Instead mulitply by -1.
  • Loading branch information
mkelley authored Aug 28, 2024
1 parent 035041e commit e566f6d
Show file tree
Hide file tree
Showing 2 changed files with 21 additions and 22 deletions.
42 changes: 20 additions & 22 deletions sbpy/spectroscopy/core.py
Original file line number Diff line number Diff line change
Expand Up @@ -9,7 +9,7 @@
import astropy.units as u


__all__ = ['Spectrum', 'SpectralModel', 'SpectralGradient']
__all__ = ["Spectrum", "SpectralModel", "SpectralGradient"]

__doctest_requires__ = {
"SpectralGradient": ["synphot"],
Expand Down Expand Up @@ -91,16 +91,14 @@ class SpectralGradient(u.SpecificTypeQuantity):
_equivalent_unit = u.meter**-1
_include_easy_conversion_members = False

def __new__(cls, value, unit=None, wave=None, wave0=None,
dtype=None, copy=True):
S = super().__new__(cls, value, unit=unit, dtype=dtype,
copy=copy)
def __new__(cls, value, unit=None, wave=None, wave0=None, dtype=None, copy=True):
S = super().__new__(cls, value, unit=unit, dtype=dtype, copy=copy)

if wave is not None:
if np.size(wave) != 2:
raise ValueError(
'Two wavelengths must be provided, got {}'
.format(np.size(wave)))
"Two wavelengths must be provided, got {}".format(np.size(wave))
)
S.wave = S._lambda_eff(wave)

if wave0 is None and wave is not None:
Expand Down Expand Up @@ -179,16 +177,15 @@ def from_color(cls, wfb, color):

lambda_eff = SpectralGradient._lambda_eff(wfb)

try:
# works for u.Magnitudes and dimensionless u.Quantity
if isinstance(color, u.Magnitude):
alpha = u.Quantity(-1 * color, u.dimensionless_unscaled)
elif u.Quantity(color).unit.is_equivalent(u.mag):
alpha = (-color).to(u.dimensionless_unscaled, u.logarithmic())
else:
alpha = u.Quantity(color, u.dimensionless_unscaled)
except u.UnitConversionError:
# works for u.mag
alpha = color.to(u.dimensionless_unscaled, u.logarithmic())

dw = lambda_eff[0] - lambda_eff[1]
S = ((2 / dw * (alpha - 1) / (alpha + 1))
.to(u.percent / hundred_nm))
dw = lambda_eff[1] - lambda_eff[0]
S = (2 / dw * (alpha - 1) / (alpha + 1)).to(u.percent / hundred_nm)

return SpectralGradient(S, wave=lambda_eff)

Expand Down Expand Up @@ -287,7 +284,7 @@ def renormalize(self, wave0):
"""

if self.wave0 is None:
raise ValueError('wave0 attribute must be defined.')
raise ValueError("wave0 attribute must be defined.")

delta = wave0 - self.wave0
S0 = 1 + self.to(delta.unit**-1) * delta
Expand All @@ -296,8 +293,9 @@ def renormalize(self, wave0):
return S


class SpectralModel():
class SpectralModel:
"""Range of spectral models"""

def haser():
"""Haser model
Expand All @@ -312,15 +310,15 @@ def reflectance():
"""Reflectance spectrum (asteroids)"""


class Spectrum():
class Spectrum:

def __init__(self, flux, dispersionaxis, unit):
self.flux = flux
self.axis = dispersionaxis
self.unit = unit

@classmethod
def read(cls, filename, columns='auto'):
def read(cls, filename, columns="auto"):
"""Read spectrum from file
Parameters
Expand All @@ -343,7 +341,7 @@ def read(cls, filename, columns='auto'):
"""

def write(self, filename, columns='all'):
def write(self, filename, columns="all"):
"""Write spectrum to file
Parameters
Expand Down Expand Up @@ -416,7 +414,7 @@ def slope(self, subtract=False):
"""

def integrated_flux(self, frequency, interval=1*u.km/u.s):
def integrated_flux(self, frequency, interval=1 * u.km / u.s):
"""
Calculate integrated flux of emission line.
Expand Down Expand Up @@ -446,7 +444,7 @@ def fit(self, spec):
Examples
--------
>>> spec_model = SpectralModel(type='Haser', molecule='H2O') # doctest: +SKIP
>>> spec_model = SpectralModel(type='Haser', molecule='H2O') # doctest: +SKIP
>>> spec.fit(spec_model) # doctest: +SKIP
>>> print(spec.fit_info) # doctest: +SKIP
Expand Down
1 change: 1 addition & 0 deletions sbpy/spectroscopy/tests/test_specgrad.py
Original file line number Diff line number Diff line change
Expand Up @@ -37,6 +37,7 @@ def test_new_wave_error(self):
0.21 * u.mag, 12 * u.percent / hundred_nm, 2),
((bandpass('SDSS g'), bandpass('SDSS r')),
-0.15 * u.mag, -10 * u.percent / hundred_nm, 0.5),
((500, 600) * u.nm, 3, 100 * u.percent / hundred_nm, 0.01),
))
def test_from_color(self, wfb, color, S0, atol):
"""Test from color to spectral gradient.
Expand Down

0 comments on commit e566f6d

Please sign in to comment.