Skip to content
105 changes: 86 additions & 19 deletions astroplan/constraints.py
Original file line number Diff line number Diff line change
Expand Up @@ -28,7 +28,8 @@
"is_observable", "is_always_observable", "time_grid_from_range",
"SunSeparationConstraint", "MoonSeparationConstraint",
"MoonIlluminationConstraint", "LocalTimeConstraint", "Constraint",
"TimeConstraint", "observability_table", "months_observable"]
"TimeConstraint", "observability_table", "months_observable",
"max_best_rescale", "min_best_rescale"]


def _get_altaz(times, observer, targets,
Expand Down Expand Up @@ -286,7 +287,7 @@ def compute_constraint(self, times, observer, targets):
uppermask = alt <= self.max
return lowermask & uppermask
else:
return _rescale_minmax(alt, self.min, self.max)
return max_best_rescale(alt, self.min, self.max)


class AirmassConstraint(AltitudeConstraint):
Expand Down Expand Up @@ -344,8 +345,8 @@ def compute_constraint(self, times, observer, targets):
mx = self.max

mi = 1 if self.min is None else self.min
# we reverse order so that airmass close to 1/min is good
return _rescale_airmass(secz, mi, mx)
# values below 1 should be disregarded
return min_best_rescale(secz, mi, mx, less_than_min=0)


class AtNightConstraint(Constraint):
Expand Down Expand Up @@ -997,24 +998,90 @@ def observability_table(constraints, observer, targets, times=None,
return tab


def _rescale_minmax(vals, min_val, max_val):
""" Rescale altitude into an observability score."""
rescaled = (vals - min_val) / (max_val - min_val)
below = rescaled < 0
above = rescaled > 1
rescaled[below] = 0
rescaled[above] = 1
def min_best_rescale(vals, min_val, max_val, less_than_min=1):
"""
rescales an input array ``vals`` to be a score (between zero and one),
where the ``min_val`` goes to one, and the ``max_val`` goes to zero.

Parameters
----------
vals : array-like
the values that need to be rescaled to be between 0 and 1
min_val : float
worst acceptable value (rescales to 0)
max_val : float
best value cared about (rescales to 1)
less_than_min : 0 or 1
what is returned for ``vals`` below ``min_val``. (in some cases
anything less than ``min_val`` should also return one,
in some cases it should return zero)

Returns
-------
array of floats between 0 and 1 inclusive rescaled so that
``vals`` equal to ``max_val`` equal 0 and those equal to
``min_val`` equal 1

Examples
--------
rescale airmasses to between 0 and 1, with the best (1)
and worst (2.25). All values outside the range should
return 0.
>>> from astroplan.constraints import min_best_rescale
>>> import numpy as np
>>> airmasses = np.array([1, 1.5, 2, 3, 0])
>>> min_best_rescale(airmasses, 1, 2.25, less_than_min = 0)
array([ 1. , 0.6, 0.2, 0. , 0. ])
"""
rescaled = (vals - max_val) / (min_val - max_val)
below = vals < min_val
above = vals > max_val
rescaled[below] = less_than_min
rescaled[above] = 0

return rescaled


def _rescale_airmass(vals, min_val, max_val):
""" Rescale airmass into an observability score."""
def max_best_rescale(vals, min_val, max_val, greater_than_max=1):
"""
rescales an input array ``vals`` to be a score (between zero and one),
where the ``max_val`` goes to one, and the ``min_val`` goes to zero.

Parameters
----------
vals : array-like
the values that need to be rescaled to be between 0 and 1
min_val : float
worst acceptable value (rescales to 0)
max_val : float
best value cared about (rescales to 1)
greater_than_max : 0 or 1
what is returned for ``vals`` above ``max_val``. (in some cases
anything higher than ``max_val`` should also return one,
in some cases it should return zero)

Returns
-------
array of floats between 0 and 1 inclusive rescaled so that
``vals`` equal to ``min_val`` equal 0 and those equal to
``max_val`` equal 1

Examples
--------
rescale an array of altitudes to be between 0 and 1,
with the best (60) going to 1 and worst (35) going to
0. For values outside the range, the rescale should
return 0 below 35 and 1 above 60.
>>> from astroplan.constraints import max_best_rescale
>>> import numpy as np
>>> altitudes = np.array([20, 30, 40, 45, 55, 70])
>>> max_best_rescale(altitudes, 35, 60)
array([ 0. , 0. , 0.2, 0.4, 0.8, 1. ])
"""
rescaled = (vals - min_val) / (max_val - min_val)
below = rescaled < 0
above = rescaled > 1
# In both cases, we want out-of-range airmasses to return a 0 score
rescaled[below] = 1
rescaled[above] = 1
below = vals < min_val
above = vals > max_val
rescaled[below] = 0
rescaled[above] = greater_than_max

return 1 - rescaled
return rescaled
13 changes: 12 additions & 1 deletion astroplan/tests/test_constraints.py
Original file line number Diff line number Diff line change
Expand Up @@ -15,7 +15,8 @@
is_observable, is_always_observable, observability_table,
time_grid_from_range, SunSeparationConstraint,
MoonSeparationConstraint, MoonIlluminationConstraint,
TimeConstraint, LocalTimeConstraint, months_observable)
TimeConstraint, LocalTimeConstraint, months_observable,
max_best_rescale, min_best_rescale)

APY_LT104 = not minversion('astropy','1.0.4')

Expand Down Expand Up @@ -369,6 +370,16 @@ def test_months_observable():
assert months == should_be


def test_rescale_minmax():
a = np.array([2])
rescaled = np.zeros(5)
rescaled[0] = (min_best_rescale(a, 1, 6))[0]
rescaled[1] = (max_best_rescale(a, 1, 6))[0]
rescaled[2] = (max_best_rescale(a, 0, 1))[0]
rescaled[3] = (min_best_rescale(a, 0, 1))[0]
rescaled[4] = (max_best_rescale(a, 0, 1, greater_than_max=0))[0]
assert all(np.array([0.8, 0.2, 1, 0, 0]) == rescaled)

constraint_tests = [
AltitudeConstraint(),
AirmassConstraint(2),
Expand Down
87 changes: 67 additions & 20 deletions docs/tutorials/constraints.rst
Original file line number Diff line number Diff line change
Expand Up @@ -242,16 +242,21 @@ be within some angular separation from Vega – we'll call it
``compute_constraints`` method when you check if a target is observable
using `~astroplan.is_observable` or `~astroplan.is_always_observable`.

* We also want to provide the option of having the constraint output
a non-boolean score. Where being closer to the minimum separation
returns a higher score than being closer to the maximum separation.

Here's our ``VegaSeparationConstraint`` implementation::

from astroplan import Constraint, is_observable
from astroplan import Constraint, is_observable, min_best_rescale
from astropy.coordinates import Angle
import astropy.units as u

class VegaSeparationConstraint(Constraint):
"""
Constraint the separation from Vega
"""
def __init__(self, min=None, max=None):
def __init__(self, min=None, max=None, boolean_constraint=True):
"""
min : `~astropy.units.Quantity` or `None` (optional)
Minimum acceptable separation between Vega and target. `None`
Expand All @@ -262,6 +267,7 @@ Here's our ``VegaSeparationConstraint`` implementation::
"""
self.min = min
self.max = max
self.boolean_constraint = boolean_constraint

def compute_constraint(self, times, observer, targets):

Expand All @@ -273,31 +279,51 @@ Here's our ``VegaSeparationConstraint`` implementation::
# Calculate separation between target and vega
vega_separation = Angle([vega.separation(target.coord)
for target in targets])
if self.boolean_constraint:
# If a maximum is specified but no minimum
if self.min is None and self.max is not None:
mask = vega_separation < self.max

# If a maximum is specified but no minimum
if self.min is None and self.max is not None:
mask = vega_separation < self.max
# If a minimum is specified but no maximum
elif self.max is None and self.min is not None:
mask = self.min < vega_separation

# If a minimum is specified but no maximum
elif self.max is None and self.min is not None:
mask = self.min < vega_separation
# If both a minimum and a maximum are specified
elif self.min is not None and self.max is not None:
mask = ((self.min < vega_separation) & (vega_separation < self.max))

# If both a minimum and a maximum are specified
elif self.min is not None and self.max is not None:
mask = ((self.min < vega_separation) & (vega_separation < self.max))
# Otherwise, raise an error
else:
raise ValueError("No max and/or min specified in "
"VegaSeparationConstraint.")

# Otherwise, raise an error
else:
raise ValueError("No max and/or min specified in "
"VegaSeparationConstraint.")

# Return an array that is True where the target is observable and
# False where it is not
# Must have shape (len(targets), len(times))

# Return an array that is True where the target is observable and
# False where it is not
# Must have shape (len(targets), len(times))
# currently mask has shape (len(targets), 1)
return np.tile(mask, len(times))

# if we want to return a non-boolean score
else:
# no min and no max still should error
if self.min is None and self.max is None:
raise ValueError("No max and/or min specified in "
"VegaSeparationConstraint.")
if self.min is None:
# if no minimum is given, set it at 0 degrees
self.min = 0*u.deg
if self.max is None:
# if no maximum is given, set it to 180 degrees
self.max = 180*u.deg

# rescale the vega_separation values so that they become
# scores between zero and one
rescale = min_best_rescale(vega_separation, self.min,
self.max, less_than_min=0)
return np.tile(rescale, len(times))

# currently mask has shape (len(targets), 1)
return np.tile(mask, len(times))

Then as in the earlier example, we can call our constraint::

Expand All @@ -310,3 +336,24 @@ Then as in the earlier example, we can call our constraint::
The resulting list of booleans indicates that the only target separated by
5 and 30 degrees from Vega is Albireo. Following this pattern, you can design
arbitrarily complex criteria for constraints.

To see the (target x time) array for the non-boolean score::

>>> constraint = VegaSeparationConstraint(min=5*u.deg, max=30*u.deg,
... boolean_constraint=False)
>>> print(constraint(subaru, targets, time_range=time_range)
[[ 0. 0. 0. 0. 0. 0. 0.
0. 0. 0. 0. 0. ]
[ 0. 0. 0. 0. 0. 0. 0.
0. 0. 0. 0. 0. ]
[ 0.57748686 0.57748686 0.57748686 0.57748686 0.57748686 0.57748686
0.57748686 0.57748686 0.57748686 0.57748686 0.57748686 0.57748686]
[ 0. 0. 0. 0. 0. 0. 0.
0. 0. 0. 0. 0. ]
[ 0. 0. 0. 0. 0. 0. 0.
0. 0. 0. 0. 0. ]
[ 0. 0. 0. 0. 0. 0. 0.
0. 0. 0. 0. 0. ]]

The score of .5775 for Albireo indicates that it is slightly closer to
the 5 degree minimum than to the 30 degree maximum.