Skip to content

Commit

Permalink
Merge pull request #413 from lpsinger/galactic-plane-separation-const…
Browse files Browse the repository at this point in the history
…raint

Add GalacticLatitudeConstraint
  • Loading branch information
bmorris3 committed Jul 23, 2019
2 parents 6e3daf1 + e8aef24 commit 7a894ac
Show file tree
Hide file tree
Showing 3 changed files with 70 additions and 9 deletions.
3 changes: 3 additions & 0 deletions CHANGES.rst
Original file line number Diff line number Diff line change
@@ -1,6 +1,9 @@
0.6 (unreleased)
----------------

- Add ``GalacticLatitudeConstraint`` to constrain the galactic latitudes of
targets. This can be useful for planning surveys for which crowding due to
Galactic point sources is an issue. [#413]

0.5 (2019-07-08)
----------------
Expand Down
48 changes: 41 additions & 7 deletions astroplan/constraints.py
Original file line number Diff line number Diff line change
Expand Up @@ -15,7 +15,7 @@
# Third-party
from astropy.time import Time
import astropy.units as u
from astropy.coordinates import get_body, get_sun, get_moon, SkyCoord
from astropy.coordinates import get_body, get_sun, get_moon, Galactic, SkyCoord
from astropy import table

import numpy as np
Expand All @@ -28,12 +28,12 @@

__all__ = ["AltitudeConstraint", "AirmassConstraint", "AtNightConstraint",
"is_observable", "is_always_observable", "time_grid_from_range",
"SunSeparationConstraint", "MoonSeparationConstraint",
"MoonIlluminationConstraint", "LocalTimeConstraint",
"PrimaryEclipseConstraint", "SecondaryEclipseConstraint",
"Constraint", "TimeConstraint", "observability_table",
"months_observable", "max_best_rescale", "min_best_rescale",
"PhaseConstraint", "is_event_observable"]
"GalacticLatitudeConstraint", "SunSeparationConstraint",
"MoonSeparationConstraint", "MoonIlluminationConstraint",
"LocalTimeConstraint", "PrimaryEclipseConstraint",
"SecondaryEclipseConstraint", "Constraint", "TimeConstraint",
"observability_table", "months_observable", "max_best_rescale",
"min_best_rescale", "PhaseConstraint", "is_event_observable"]


def _make_cache_key(times, targets):
Expand Down Expand Up @@ -482,6 +482,40 @@ def compute_constraint(self, times, observer, targets):
return mask


class GalacticLatitudeConstraint(Constraint):
"""
Constrain the distance between the Galactic plane and some targets.
"""

def __init__(self, min=None, max=None):
"""
Parameters
----------
min : `~astropy.units.Quantity` or `None` (optional)
Minimum acceptable Galactic latitude of target (inclusive).
`None` indicates no limit.
max : `~astropy.units.Quantity` or `None` (optional)
Minimum acceptable Galactic latitude of target (inclusive).
`None` indicates no limit.
"""
self.min = min
self.max = max

def compute_constraint(self, times, observer, targets):
separation = abs(targets.transform_to(Galactic).b)

if self.min is None and self.max is not None:
mask = self.max >= separation
elif self.max is None and self.min is not None:
mask = self.min <= separation
elif self.min is not None and self.max is not None:
mask = ((self.min <= separation) & (separation <= self.max))
else:
raise ValueError("No max and/or min specified in "
"GalacticLatitudeConstraint.")
return mask


class SunSeparationConstraint(Constraint):
"""
Constrain the distance between the Sun and some targets.
Expand Down
28 changes: 26 additions & 2 deletions astroplan/tests/test_constraints.py
Original file line number Diff line number Diff line change
Expand Up @@ -5,15 +5,17 @@
import numpy as np
import astropy.units as u
from astropy.time import Time
from astropy.coordinates import SkyCoord, get_sun, get_moon
from astropy.coordinates import Galactic, SkyCoord, get_sun, get_moon
from astropy.utils import minversion
import pytest

from ..observer import Observer
from ..target import FixedTarget, get_skycoord
from ..constraints import (AltitudeConstraint, AirmassConstraint, AtNightConstraint,
is_observable, is_always_observable, observability_table,
time_grid_from_range, SunSeparationConstraint,
time_grid_from_range,
GalacticLatitudeConstraint,
SunSeparationConstraint,
MoonSeparationConstraint, MoonIlluminationConstraint,
TimeConstraint, LocalTimeConstraint, months_observable,
max_best_rescale, min_best_rescale, PhaseConstraint,
Expand Down Expand Up @@ -136,6 +138,28 @@ def test_compare_airmass_constraint_and_observer():
assert all(always_from_observer == always_from_constraint)


def test_galactic_plane_separation():
time = Time('2003-04-05 06:07:08')
apo = Observer.at_site("APO")
one_deg_away = SkyCoord(0*u.deg, 1*u.deg, frame=Galactic)
five_deg_away = SkyCoord(0*u.deg, 5*u.deg, frame=Galactic)
twenty_deg_away = SkyCoord(0*u.deg, 20*u.deg, frame=Galactic)

constraint = GalacticLatitudeConstraint(min=2*u.deg, max=10*u.deg)
is_constraint_met = constraint(apo, [one_deg_away, five_deg_away,
twenty_deg_away], times=time)
assert np.all(is_constraint_met == [False, True, False])

constraint = GalacticLatitudeConstraint(max=10*u.deg)
is_constraint_met = constraint(apo, [one_deg_away, five_deg_away,
twenty_deg_away], times=time)
assert np.all(is_constraint_met == [True, True, False])

constraint = GalacticLatitudeConstraint(min=2*u.deg)
is_constraint_met = constraint(apo, [one_deg_away, five_deg_away,
twenty_deg_away], times=time)
assert np.all(is_constraint_met == [False, True, True])

# in astropy before v1.0.4, a recursion error is triggered by this test
@pytest.mark.skipif('APY_LT104')
def test_sun_separation():
Expand Down

0 comments on commit 7a894ac

Please sign in to comment.