From db1188137ff0aed1c50a6b91e8457cc61ae6cf88 Mon Sep 17 00:00:00 2001 From: ramonefoster Date: Fri, 9 Aug 2024 15:57:08 -0300 Subject: [PATCH] add near meridian constraint --- astroplan/constraints.py | 42 ++++++++++++++++++++++++++++- astroplan/tests/test_constraints.py | 17 ++++++++++-- 2 files changed, 56 insertions(+), 3 deletions(-) diff --git a/astroplan/constraints.py b/astroplan/constraints.py index f9857fe9..30df8145 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"] + "min_best_rescale", "PhaseConstraint", "is_event_observable", "NearMeridianConstraint"] _current_year = time.localtime().tm_year # needed for backward compatibility _current_year_time_range = Time( # needed for backward compatibility @@ -938,6 +938,46 @@ def compute_constraint(self, times, observer=None, targets=None): return mask +class NearMeridianConstraint(Constraint): + """ + Constraint near the Meridian. + """ + def __init__(self, min=None, boolean_constraint=True): + """ + Parameters + ---------- + min : `~astropy.units.Quantity` or `None`, optional + Minimum acceptable distance to meridian. + `None` indicates no limit. + + boolean_constraint : bool + + Examples + -------- + Constrain observations to targets that are 3 degrees away from the meridian. + >>> import astropy.units as u + >>> constraint = NearMeridianConstraint(min=3*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.boolean_constraint = boolean_constraint + + def compute_constraint(self, times, observer, targets): + lst = observer.local_sidereal_time(times) + meridian = SkyCoord(ra=lst, dec=targets.dec) + + meridian_separation = meridian.separation(targets) + + if self.boolean_constraint: + mask = (self.min < meridian_separation) + return mask + else: + rescale = min_best_rescale(meridian_separation, self.min, less_than_min=0) + return rescale + + def is_always_observable(constraints, observer, targets, times=None, time_range=None, time_grid_resolution=0.5*u.hour): """ diff --git a/astroplan/tests/test_constraints.py b/astroplan/tests/test_constraints.py index c0f66f39..4d1e5ebc 100644 --- a/astroplan/tests/test_constraints.py +++ b/astroplan/tests/test_constraints.py @@ -5,7 +5,7 @@ import numpy as np import astropy.units as u from astropy.time import Time -from astropy.coordinates import Galactic, SkyCoord, get_sun, get_body +from astropy.coordinates import Galactic, SkyCoord, EarthLocation, get_sun, get_body from astropy.utils import minversion import pytest @@ -20,7 +20,7 @@ TimeConstraint, LocalTimeConstraint, months_observable, max_best_rescale, min_best_rescale, PhaseConstraint, PrimaryEclipseConstraint, SecondaryEclipseConstraint, - is_event_observable) + is_event_observable, NearMeridianConstraint) from ..periodic import EclipsingSystem from ..exceptions import MissingConstraintWarning @@ -200,6 +200,18 @@ 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") + + # Pico dos Dias Observatory (Brazil) + opd = Observer(location=EarthLocation(lat=-22.53, lon=-45.58, height=1864)) + constraint = NearMeridianConstraint(min=3*u.deg) + + results = constraint(opd, target, times=time_grid_from_range(time_range)) + assert np.all(results == [True, True, False, True, True]) + + def test_moon_separation(): time = Time('2003-04-05 06:07:08') apo = Observer.at_site("APO") @@ -419,6 +431,7 @@ def test_rescale_minmax(): AtNightConstraint(), SunSeparationConstraint(min=90*u.deg), MoonSeparationConstraint(min=20*u.deg), + NearMeridianConstraint(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"])) ]