Skip to content

Commit

Permalink
Fix test failures with older Numpy/Astropy.
Browse files Browse the repository at this point in the history
Make integrate methods in models.py consistent.

Add tests.
  • Loading branch information
pllim committed Mar 11, 2020
1 parent 987b6eb commit 4a29896
Show file tree
Hide file tree
Showing 3 changed files with 26 additions and 9 deletions.
25 changes: 17 additions & 8 deletions synphot/models.py
Original file line number Diff line number Diff line change
Expand Up @@ -136,8 +136,10 @@ def evaluate(x, temperature):
return bbflux.value

def integrate(self, *args):
return (const.sigma_sb * (u.Quantity(self.temperature, u.K) ** 4) /
math.pi) # per steradian
with u.add_enabled_equivalencies(u.temperature()):
t = u.Quantity(self.temperature, u.K)

return (const.sigma_sb * t ** 4 / math.pi) # per steradian


class BlackBodyNorm1D(BlackBody1D):
Expand Down Expand Up @@ -225,7 +227,10 @@ def sampleset(self, step=0.01, minimal=False):

def integrate(self, *args):
# TODO: Remove unit hardcoding when we use model with units natively.
return self.amplitude * (self.width * u.AA)
with u.add_enabled_equivalencies(u.spectral()):
w = u.Quantity(self.width, u.AA)

return self.amplitude * w


class ConstFlux1D(_models.Const1D):
Expand Down Expand Up @@ -290,8 +295,9 @@ def integrate(self, x):
wav_unit = u.Hz
with u.add_enabled_equivalencies(u.spectral()):
x = u.Quantity(x, wav_unit)
amp = u.Quantity(self.amplitude, self._flux_unit)

return (max(x) - min(x)) * (self.amplitude * self._flux_unit)
return (max(x) - min(x)) * amp


class Empirical1D(Tabular1D):
Expand Down Expand Up @@ -474,9 +480,12 @@ class Gaussian1D(BaseGaussian1D):
``sampleset`` defined.
"""
# TODO: Remove unit hardcoding when we use model with units natively.
def integrate(self, *args):
return self.amplitude * (self.stddev * u.AA) * self._sqrt_2_pi
# TODO: Remove unit hardcoding when we use model with units natively.
with u.add_enabled_equivalencies(u.spectral()):
stddev = u.Quantity(self.stddev, u.AA)

return self.amplitude * stddev * self._sqrt_2_pi


# TODO: Deprecate this?
Expand Down Expand Up @@ -559,8 +568,8 @@ def __init__(self, *args, **kwargs):
self.meta['expr'] = 'em({0:g}, {1:g}, {2:g}, {3})'.format(
self.mean.value, fwhm, total_flux, u_str)

# TODO: Remove unit hardcoding when we use model with units natively.
def integrate(self, *args):
# TODO: Remove unit hardcoding when we use model with units natively.
return super(GaussianFlux1D, self).integrate(*args) * units.PHOTLAM


Expand Down Expand Up @@ -592,8 +601,8 @@ def sampleset(self, factor_step=0.05, **kwargs):

return np.asarray(w)

# TODO: Remove unit hardcoding when we use model with units natively.
def integrate(self, x):
# TODO: Remove unit hardcoding when we use model with units natively.
with u.add_enabled_equivalencies(u.spectral()):
x = u.Quantity(x, u.AA)
x_0 = u.Quantity(self.x_0, u.AA)
Expand Down
2 changes: 1 addition & 1 deletion synphot/spectrum.py
Original file line number Diff line number Diff line change
Expand Up @@ -1539,7 +1539,7 @@ def check_overlap(self, other, wavelengths=None, threshold=0.01):
other.model.is_tapered() or
not isinstance(other.model,
(Empirical1D, CompoundModel))) and
np.allclose(other(x1[::x1.size - 1]), 0)):
np.allclose(other(x1[::x1.size - 1]).value, 0)):
result = 'full'

# Check if the lack of overlap is significant.
Expand Down
8 changes: 8 additions & 0 deletions synphot/tests/test_integrator.py
Original file line number Diff line number Diff line change
Expand Up @@ -157,6 +157,14 @@ def test_bandpass_RickerWavelet1D():
assert_quantity_allclose(bp.equivwidth(integration_type='analytical'),
ans, rtol=5e-3)

with pytest.raises(NotImplementedError, match='Partial analytic'):
bp.equivwidth(integration_type='analytical',
wavelengths=np.arange(4900, 5200) * u.AA)

with pytest.raises(NotImplementedError, match='Partial analytic'):
bp.equivwidth(integration_type='analytical',
wavelengths=np.arange(4800, 5100) * u.AA)


def test_source_RickerWavelet1D():
pytest.xfail(
Expand Down

0 comments on commit 4a29896

Please # to comment.