Division added.
This commit is contained in:
parent
f0cba7ee92
commit
fbbfa287d4
@ -272,6 +272,46 @@ class SpectralQty:
|
||||
|
||||
__rmul__ = __mul__
|
||||
|
||||
def __truediv__(self, other: Union[int, float, u.Quantity, "SpectralQty", Callable[[u.Quantity], u.Quantity]]) ->\
|
||||
"SpectralQty":
|
||||
"""
|
||||
Calculate the quotient with another object
|
||||
|
||||
Parameters
|
||||
----------
|
||||
other : Union[int, float, u.Quantity, "SpectralQty", Callable]
|
||||
Divisor for this object. If the binning of the object on the right hand side differs
|
||||
from the binning of the left object, the object on the right hand side will be rebinned.
|
||||
|
||||
Returns
|
||||
-------
|
||||
sum : SpectralQty
|
||||
The quotient of both objects
|
||||
"""
|
||||
# Factor is of type int, float or Quantity, just multiply
|
||||
if isinstance(other, int) or isinstance(other, float) or isinstance(other, u.Quantity):
|
||||
return SpectralQty(self.wl, self.qty / other)
|
||||
# Factor is of type lambda
|
||||
elif isLambda(other):
|
||||
return SpectralQty(self.wl, self.qty / [other(wl).value for wl in self.wl] * other(self.wl[0]).unit)
|
||||
# Factor is of type SpectralQty
|
||||
else:
|
||||
if other.wl.unit.is_equivalent(self.wl.unit):
|
||||
# Wavelengths are matching, just multiply the quantities
|
||||
if len(self.wl) == len(other.wl) and (self.wl == other.wl).all():
|
||||
return SpectralQty(self.wl, self.qty / other.qty)
|
||||
# Wavelengths are not matching, rebinning needed
|
||||
else:
|
||||
# Rebin factor
|
||||
other_rebinned = other.rebin(self.wl)
|
||||
if len(self.wl) == len(other_rebinned.wl) and (self.wl == other_rebinned.wl).all():
|
||||
return SpectralQty(self.wl, self.qty / other_rebinned.qty)
|
||||
else:
|
||||
# Wavelengths are still not matching as extrapolation is disabled, rebin this spectral quantity
|
||||
return SpectralQty(other_rebinned.wl, self.rebin(other_rebinned.wl).qty / other_rebinned.qty)
|
||||
else:
|
||||
error("Units are not matching for division.")
|
||||
|
||||
def rebin(self, wl: u.Quantity) -> "SpectralQty":
|
||||
"""
|
||||
Resample the spectral quantity sqty(wl) over the new grid wl, rebinning if necessary, otherwise interpolates.
|
||||
|
@ -43,6 +43,29 @@ class TestSpectralQty(TestCase):
|
||||
self.assertEqual(self.sqty * (lambda wl: 0.7 * u.dimensionless_unscaled),
|
||||
SpectralQty(self.wl, self.qty * 0.7))
|
||||
|
||||
def test___truediv__(self):
|
||||
# Integer
|
||||
self.assertEqual(self.sqty / 2, SpectralQty(np.arange(200, 204, 1) << u.nm,
|
||||
np.arange(5.5e-1, 7.5e-1, 5e-2) << u.W / (u.m ** 2 * u.nm)))
|
||||
# Float
|
||||
self.assertEqual(self.sqty / 2., SpectralQty(np.arange(200, 204, 1) << u.nm,
|
||||
np.arange(5.5e-1, 7.5e-1, 5e-2) << u.W / (u.m ** 2 * u.nm)))
|
||||
# SpectralQty
|
||||
self.assertEqual(self.sqty / SpectralQty(self.wl, np.arange(1, 5, 1) << u.m),
|
||||
SpectralQty(self.wl, np.array([1.1, 0.6, 1.3 / 3, 0.35]) << u.W / (u.m ** 3 * u.nm)))
|
||||
# rebin without extrapolation and without reduction
|
||||
self.assertEqual(
|
||||
self.sqty / SpectralQty(np.arange(199.5, 204.5, 1) << u.nm, np.arange(1, 6, 1) << u.m),
|
||||
SpectralQty(self.wl, [1.1 / 1.5, 1.2 / 2.5, 1.3 / 3.5, 1.4 / 4.5] * u.W / (u.m ** 3 * u.nm)))
|
||||
# rebin without extrapolation and with reduction
|
||||
self.assertEqual(
|
||||
self.sqty / SpectralQty(np.arange(200.5, 204.5, 1) << u.nm, np.arange(1, 5, 1) << u.m, fill_value=False),
|
||||
SpectralQty(np.arange(201, 204) << u.nm, np.array([1.2 / 1.5, 1.3 / 2.5, 1.4 / 3.5]) << u.W /
|
||||
(u.m ** 3 * u.nm)))
|
||||
# lambda
|
||||
self.assertEqual(self.sqty / (lambda wl: 0.7 * u.dimensionless_unscaled),
|
||||
SpectralQty(self.wl, self.qty / 0.7))
|
||||
|
||||
def test___sub__(self):
|
||||
# Quantity
|
||||
self.assertEqual(self.sqty - 0.1 * u.W / (u.m ** 2 * u.nm),
|
||||
|
Loading…
Reference in New Issue
Block a user