diff --git a/astroplan/constraints.py b/astroplan/constraints.py index 30df8145..9244ad05 100644 --- a/astroplan/constraints.py +++ b/astroplan/constraints.py @@ -35,7 +35,7 @@ "LocalTimeConstraint", "PrimaryEclipseConstraint", "SecondaryEclipseConstraint", "Constraint", "TimeConstraint", "observability_table", "months_observable", "max_best_rescale", - "min_best_rescale", "PhaseConstraint", "is_event_observable", "NearMeridianConstraint"] + "min_best_rescale", "PhaseConstraint", "is_event_observable", "MeridianSeparationConstraint"] _current_year = time.localtime().tm_year # needed for backward compatibility _current_year_time_range = Time( # needed for backward compatibility @@ -938,30 +938,36 @@ def compute_constraint(self, times, observer=None, targets=None): return mask -class NearMeridianConstraint(Constraint): +class MeridianSeparationConstraint(Constraint): """ - Constraint near the Meridian. + Constraint on angular separation from the meridian. """ - def __init__(self, min=None, boolean_constraint=True): + def __init__(self, min=None, max=None, boolean_constraint=True): """ Parameters ---------- min : `~astropy.units.Quantity` or `None`, optional Minimum acceptable distance to meridian. - `None` indicates no limit. + `None` indicates no lower limit. + + max : `~astropy.units.Quantity` or `None`, optional + Maximum acceptable distance to meridian. + `None` indicates no upper limit. boolean_constraint : bool Examples -------- - Constrain observations to targets that are 3 degrees away from the meridian. + Constrain observations to targets that are between 3 and 35 degrees + away from the meridian. >>> import astropy.units as u - >>> constraint = NearMeridianConstraint(min=3*u.deg) + >>> constraint = MeridianSeparationConstraint(min=3*u.deg, max=35*u.deg) This can be useful for observations using German-Equatorial Mounts, to avoid flipping the side of the pier during exposures. """ - self.min = min if min is not None else 0*u.hourangle + self.min = min if min is not None else 0*u.deg + self.max = max if max is not None else 180*u.deg self.boolean_constraint = boolean_constraint def compute_constraint(self, times, observer, targets): @@ -971,10 +977,10 @@ def compute_constraint(self, times, observer, targets): meridian_separation = meridian.separation(targets) if self.boolean_constraint: - mask = (self.min < meridian_separation) + mask = (self.min < meridian_separation) & (meridian_separation < self.max) return mask else: - rescale = min_best_rescale(meridian_separation, self.min, less_than_min=0) + rescale = min_best_rescale(meridian_separation, self.min, self.max, less_than_min=0) return rescale diff --git a/astroplan/tests/test_constraints.py b/astroplan/tests/test_constraints.py index 4d1e5ebc..ee528885 100644 --- a/astroplan/tests/test_constraints.py +++ b/astroplan/tests/test_constraints.py @@ -20,7 +20,7 @@ TimeConstraint, LocalTimeConstraint, months_observable, max_best_rescale, min_best_rescale, PhaseConstraint, PrimaryEclipseConstraint, SecondaryEclipseConstraint, - is_event_observable, NearMeridianConstraint) + is_event_observable, MeridianSeparationConstraint) from ..periodic import EclipsingSystem from ..exceptions import MissingConstraintWarning @@ -200,16 +200,16 @@ def test_sun_separation(): assert np.all(is_constraint_met == [False, True, True]) -def test_near_meridian(): - time_range = Time(["2024-10-08 20:28", "2024-10-08 22:30"]) - target = FixedTarget(coord=SkyCoord(ra=19.75*u.hour, dec=-22.05*u.deg), name="name") +def test_meridian_separation(): + time_range = Time(["2024-10-08 21:00", "2024-10-08 23:00"]) + target = FixedTarget(coord=SkyCoord(ra=19.75*u.hour, dec=-22.50*u.deg), name="name") # Pico dos Dias Observatory (Brazil) opd = Observer(location=EarthLocation(lat=-22.53, lon=-45.58, height=1864)) - constraint = NearMeridianConstraint(min=3*u.deg) + constraint = MeridianSeparationConstraint(min=3*u.deg, max=35*u.deg) results = constraint(opd, target, times=time_grid_from_range(time_range)) - assert np.all(results == [True, True, False, True, True]) + assert np.all(results == [True, False, True, True, True]) def test_moon_separation(): @@ -431,7 +431,7 @@ def test_rescale_minmax(): AtNightConstraint(), SunSeparationConstraint(min=90*u.deg), MoonSeparationConstraint(min=20*u.deg), - NearMeridianConstraint(min=3*u.deg), + MeridianSeparationConstraint(min=3*u.deg), LocalTimeConstraint(min=dt.time(23, 50), max=dt.time(4, 8)), TimeConstraint(*Time(["2015-08-28 03:30", "2015-09-05 10:30"])) ]