Skip to content

Commit 235db81

Browse files
author
Brett M. Morris
authored
Merge pull request #213 from bmorris3/fix-issue-211
Removing location kwarg from moon illumination/phase computations
2 parents 7657ab7 + 70f0274 commit 235db81

File tree

3 files changed

+34
-50
lines changed

3 files changed

+34
-50
lines changed

astroplan/constraints.py

Lines changed: 16 additions & 11 deletions
Original file line numberDiff line numberDiff line change
@@ -83,10 +83,10 @@ def _get_altaz(times, observer, targets,
8383
return observer._altaz_cache[aakey]
8484

8585

86-
def _get_moon_data(times, observer,
87-
force_zero_pressure=False):
86+
def _get_moon_data(times, observer, force_zero_pressure=False):
8887
"""
89-
Calculate moon altitude az and illumination for an array of times for ``observer``.
88+
Calculate moon altitude az and illumination for an array of times for
89+
``observer``.
9090
9191
Cache the result on the ``observer`` object.
9292
@@ -119,8 +119,7 @@ def _get_moon_data(times, observer,
119119
observer.pressure = 0
120120

121121
altaz = observer.moon_altaz(times)
122-
illumination = np.array(moon_illumination(times,
123-
observer.location))
122+
illumination = np.array(moon_illumination(times))
124123
observer._moon_cache[aakey] = dict(times=times,
125124
illum=illumination,
126125
altaz=altaz)
@@ -133,7 +132,8 @@ def _get_moon_data(times, observer,
133132

134133
def _get_meridian_transit_times(times, observer, targets):
135134
"""
136-
Calculate next meridian transit for an array of times for ``targets`` and ``observer``.
135+
Calculate next meridian transit for an array of times for ``targets`` and
136+
``observer``.
137137
138138
Cache the result on the ``observer`` object.
139139
@@ -162,7 +162,8 @@ def _get_meridian_transit_times(times, observer, targets):
162162

163163
if aakey not in observer._meridian_transit_cache:
164164
meridian_transit_times = Time([observer.target_meridian_transit_time(
165-
time, target, which='next') for target in targets
165+
time, target, which='next')
166+
for target in targets
166167
for time in times])
167168
observer._meridian_transit_cache[aakey] = dict(times=meridian_transit_times)
168169

@@ -735,15 +736,18 @@ def __init__(self, min=None, max=None):
735736
self.max = max
736737

737738
if self.min is None and self.max is None:
738-
raise ValueError("You must at least supply either a minimum or a maximum time.")
739+
raise ValueError("You must at least supply either a minimum or a "
740+
"maximum time.")
739741

740742
if self.min is not None:
741743
if not isinstance(self.min, Time):
742-
raise TypeError("Time limits must be specified as astropy.time.Time objects.")
744+
raise TypeError("Time limits must be specified as "
745+
"astropy.time.Time objects.")
743746

744747
if self.max is not None:
745748
if not isinstance(self.max, Time):
746-
raise TypeError("Time limits must be specified as astropy.time.Time objects.")
749+
raise TypeError("Time limits must be specified as "
750+
"astropy.time.Time objects.")
747751

748752
def compute_constraint(self, times, observer, targets):
749753
with warnings.catch_warnings():
@@ -979,7 +983,8 @@ def observability_table(constraints, observer, targets, times=None,
979983
always_obs = np.all(contraint_arr, axis=1)
980984
frac_obs = np.sum(contraint_arr, axis=1) / contraint_arr.shape[1]
981985

982-
tab = table.Table(names=colnames, data=[target_names, ever_obs, always_obs, frac_obs])
986+
tab = table.Table(names=colnames, data=[target_names, ever_obs, always_obs,
987+
frac_obs])
983988

984989
if times is None and time_range is not None:
985990
times = time_grid_from_range(time_range,

astroplan/moon.py

Lines changed: 16 additions & 37 deletions
Original file line numberDiff line numberDiff line change
@@ -1,23 +1,21 @@
11
# Licensed under a 3-clause BSD style license - see LICENSE.rst
22
"""
3-
This version of the `moon` module contains temporary solutions to calculating
4-
the position of the moon using astropy.coordinates.get_moon, awaiting solutions
5-
to the astropy issue [#5069](https://github.com/astropy/astropy/issues/5069).
3+
This version of the `moon` module calculates lunar phase angle for a geocentric
64
"""
75

86
from __future__ import (absolute_import, division, print_function,
97
unicode_literals)
108

119
# Third-party
1210
import numpy as np
13-
from astropy.coordinates import (get_moon, get_sun, AltAz, Angle)
11+
from astropy.coordinates import get_moon, get_sun
1412

1513
__all__ = ["moon_phase_angle", "moon_illumination"]
1614

1715

18-
def moon_phase_angle(time, location, ephemeris=None):
16+
def moon_phase_angle(time, ephemeris=None):
1917
"""
20-
Calculate lunar orbital phase [radians].
18+
Calculate lunar orbital phase in radians.
2119
2220
Parameters
2321
----------
@@ -29,7 +27,7 @@ def moon_phase_angle(time, location, ephemeris=None):
2927
3028
ephemeris : str, optional
3129
Ephemeris to use. If not given, use the one set with
32-
``astropy.coordinates.solar_system_ephemeris.set`` (which is
30+
`~astropy.coordinates.solar_system_ephemeris` (which is
3331
set to 'builtin' by default).
3432
3533
Returns
@@ -39,35 +37,16 @@ def moon_phase_angle(time, location, ephemeris=None):
3937
"""
4038
# TODO: cache these sun/moon SkyCoord objects
4139

42-
# TODO: when astropy/astropy#5069 is resolved, replace this workaround which
43-
# handles scalar and non-scalar time inputs differently
44-
45-
if time.isscalar:
46-
altaz_frame = AltAz(location=location, obstime=time)
47-
sun = get_sun(time).transform_to(altaz_frame)
48-
49-
moon = get_moon(time, location=location, ephemeris=ephemeris).transform_to(altaz_frame)
50-
elongation = sun.separation(moon)
51-
return np.arctan2(sun.distance*np.sin(elongation),
52-
moon.distance - sun.distance*np.cos(elongation))
53-
54-
else:
55-
phase_angles = []
56-
for t in time:
57-
altaz_frame = AltAz(location=location, obstime=t)
58-
moon_coord = get_moon(t, location=location, ephemeris=ephemeris).transform_to(altaz_frame)
59-
sun_coord = get_sun(t).transform_to(altaz_frame)
60-
elongation = sun_coord.separation(moon_coord)
61-
phase_angle = np.arctan2(sun_coord.distance*np.sin(elongation),
62-
moon_coord.distance -
63-
sun_coord.distance*np.cos(elongation))
64-
phase_angles.append(phase_angle)
65-
return Angle(phase_angles)
66-
67-
68-
def moon_illumination(time, location, ephemeris=None):
40+
sun = get_sun(time)
41+
moon = get_moon(time, ephemeris=ephemeris)
42+
elongation = sun.separation(moon)
43+
return np.arctan2(sun.distance*np.sin(elongation),
44+
moon.distance - sun.distance*np.cos(elongation))
45+
46+
47+
def moon_illumination(time, ephemeris=None):
6948
"""
70-
Calculate fraction of the moon illuminated
49+
Calculate fraction of the moon illuminated.
7150
7251
Parameters
7352
----------
@@ -79,14 +58,14 @@ def moon_illumination(time, location, ephemeris=None):
7958
8059
ephemeris : str, optional
8160
Ephemeris to use. If not given, use the one set with
82-
``astropy.coordinates.solar_system_ephemeris.set`` (which is
61+
`~astropy.coordinates.solar_system_ephemeris` (which is
8362
set to 'builtin' by default).
8463
8564
Returns
8665
-------
8766
k : float
8867
Fraction of moon illuminated
8968
"""
90-
i = moon_phase_angle(time, location, ephemeris=ephemeris)
69+
i = moon_phase_angle(time, ephemeris=ephemeris)
9170
k = (1 + np.cos(i))/2.0
9271
return k.value

astroplan/observer.py

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -1398,7 +1398,7 @@ def moon_illumination(self, time):
13981398
if not isinstance(time, Time):
13991399
time = Time(time)
14001400

1401-
return moon_illumination(time, self.location)
1401+
return moon_illumination(time)
14021402

14031403
def moon_phase(self, time=None):
14041404
"""
@@ -1435,7 +1435,7 @@ def moon_phase(self, time=None):
14351435
if time is not None and not isinstance(time, Time):
14361436
time = Time(time)
14371437

1438-
return moon_phase_angle(time, self.location)
1438+
return moon_phase_angle(time)
14391439

14401440
def moon_altaz(self, time, ephemeris=None):
14411441
"""

0 commit comments

Comments
 (0)