File: periodic.py

package info (click to toggle)
astroplan 0.7-1
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 1,476 kB
  • sloc: python: 9,083; javascript: 161; makefile: 132; ansic: 88
file content (295 lines) | stat: -rw-r--r-- 10,665 bytes parent folder | download | duplicates (4)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
# Licensed under a 3-clause BSD style license - see LICENSE.rst
from __future__ import (absolute_import, division, print_function,
                        unicode_literals)
import numpy as np
import astropy.units as u
from astropy.time import Time

__all__ = ['PeriodicEvent', 'EclipsingSystem']


class PeriodicEvent(object):
    """
    A periodic event defined by an epoch and period.
    """
    @u.quantity_input(period=u.day, duration=u.day)
    def __init__(self, epoch, period, duration=None, name=None):
        """

        Parameters
        ----------
        epoch : `~astropy.time.Time`
            Time of event
        period : `~astropy.units.Quantity`
            Period of event
        duration : `~astropy.units.Quantity` (optional)
            Duration of event
        name : str (optional)
            Name of target/event
        """
        self.epoch = epoch
        self.period = period
        self.name = name
        self.duration = duration

    def phase(self, time):
        """
        Phase of periodic event, on interval [0, 1). For example, the phase
        could be an orbital phase for an eclipsing binary system.

        Parameters
        ----------
        time : `~astropy.time.Time`
            Evaluate the phase at this time or times

        Returns
        -------
        phase_array : `~numpy.ndarray`
            Phase at each ``time``, on range [0, 1)
        """
        return ((time - self.epoch).to(u.day).value %
                self.period.to(u.day).value) / self.period.to(u.day).value


class EclipsingSystem(PeriodicEvent):
    """
    Define parameters for an eclipsing system; useful for an eclipsing binary or
    transiting exoplanet.

    .. warning::
        There are currently two major caveats in the implementation of
        ``EclipsingSystem``. The secondary eclipse time approximation is
        only accurate when the orbital eccentricity is small, and the eclipse
        times are computed without any barycentric corrections. The current
        implementation should only be used forapproximate mid-eclipse times for
        low eccentricity orbits, with event durations longer than the
        barycentric correction error (<=16 minutes).
    """
    @u.quantity_input(period=u.day, duration=u.day)
    def __init__(self, primary_eclipse_time, orbital_period, duration=None,
                 name=None, eccentricity=None, argument_of_periapsis=None):
        """
        Parameters
        ----------
        primary_eclipse_time : `~astropy.time.Time`
            Time of primary eclipse
        orbital_period : `~astropy.units.Quantity`
            Orbital period of eclipsing system
        duration : `~astropy.units.Quantity` (optional)
            Duration of eclipse
        name : str (optional)
            Name of target/event
        eccentricity : float (optional)
            Orbital eccentricity. Default is `None`, which assumes circular
            orbit (e=0).
        argument_of_periapsis : float (optional)
            Argument of periapsis for the eclipsing system, in radians.
            Default is `None`, which assumes pi/2.
        """
        self.epoch = primary_eclipse_time
        self.period = orbital_period
        self.name = name
        self.duration = duration

        if eccentricity is None:
            eccentricity = 0
        self.eccentricity = eccentricity

        if argument_of_periapsis is None:
            argument_of_periapsis = np.pi/2
        self.argument_of_periapsis = argument_of_periapsis

    def in_primary_eclipse(self, time):
        """
        Returns `True` when ``time`` is during a primary eclipse.

        .. warning::
            Barycentric offsets are ignored in the current implementation.

        Parameters
        ----------
        time : `~astropy.time.Time`
            Time to evaluate

        Returns
        -------
        in_eclipse : `~numpy.ndarray` or bool
            `True` if ``time`` is during primary eclipse
        """
        phases = self.phase(time)
        return ((phases < float(self.duration/self.period)/2) |
                (phases > 1 - float(self.duration/self.period)/2))

    def in_secondary_eclipse(self, time):
        r"""
        Returns `True` when ``time`` is during a secondary eclipse

        If the eccentricity of the eclipsing system is non-zero, then we compute
        the secondary eclipse time approximated to first order in eccentricity,
        as described in Winn (2010) Equation 33 [1]_:

        The time between the primary eclipse and secondary eclipse :math:`\delta t_c`
        is given by :math:`\delta t_c \approx 0.5 \left (\frac{4}{\pi} e \cos{\omega \right)`,
        where :math:`e` is the orbital eccentricity and :math:`\omega` is the
        angle of periapsis.

        .. warning::
            This approximation for the secondary eclipse time is only accurate
            when the orbital eccentricity is small; and barycentric offsets
            are ignored in the current implementation.

        Parameters
        ----------
        time : `~astropy.time.Time`
            Time to evaluate

        Returns
        -------
        in_eclipse : `~numpy.ndarray` or bool
            `True` if ``time`` is during secondary eclipse

        References
        ----------
        .. [1] Winn (2010) https://arxiv.org/abs/1001.2010
        """
        if self.eccentricity < 1e-5:
            secondary_eclipse_phase = 0.5
        else:
            secondary_eclipse_phase = 0.5 * (1 + 4/np.pi * self.eccentricity *
                                             np.cos(self.argument_of_periapsis))
        phases = self.phase(time)
        return ((phases < secondary_eclipse_phase + float(self.duration/self.period)/2) &
                (phases > secondary_eclipse_phase - float(self.duration/self.period)/2))

    def out_of_eclipse(self, time):
        """
        Returns `True` when ``time`` is not during primary or secondary eclipse.

        .. warning::
            Barycentric offsets are ignored in the current implementation.

        Parameters
        ----------
        time : `~astropy.time.Time`
            Time to evaluate

        Returns
        -------
        in_eclipse : `~numpy.ndarray` or bool
            `True` if ``time`` is not during primary or secondary eclipse
        """
        return np.logical_not(np.logical_or(self.in_primary_eclipse(time),
                                            self.in_secondary_eclipse(time)))

    def next_primary_eclipse_time(self, time, n_eclipses=1):
        """
        Time of the next primary eclipse after ``time``.

        .. warning::
            Barycentric offsets are ignored in the current implementation.

        Parameters
        ----------
        time : `~astropy.time.Time`
            Find the next primary eclipse after ``time``
        n_eclipses : int (optional)
            Return the times of eclipse for the next ``n_eclipses`` after
            ``time``. Default is 1.

        Returns
        -------
        primary_eclipses : `~astropy.time.Time`
            Times of the next ``n_eclipses`` primary eclipses after ``time``
        """
        eclipse_times = ((1-self.phase(time)) * self.period + time +
                         np.arange(n_eclipses) * self.period)
        return eclipse_times

    def next_secondary_eclipse_time(self, time, n_eclipses=1):
        """
        Time of the next secondary eclipse after ``time``.

        .. warning::
            Barycentric offsets are ignored in the current implementation.

        Parameters
        ----------
        time : `~astropy.time.Time`
            Find the next secondary eclipse after ``time``
        n_eclipses : int (optional)
            Return the times of eclipse for the next ``n_eclipses`` after
            ``time``. Default is 1.

        Returns
        -------
        secondary_eclipses : `~astropy.time.Time`
            Times of the next ``n_eclipses`` secondary eclipses after ``time``
        """
        phase = self.phase(time)
        if phase >= 0.5:
            next_eclipse_phase = 1.5
        else:
            next_eclipse_phase = 0.5
        eclipse_times = ((next_eclipse_phase - phase) * self.period + time +
                         np.arange(n_eclipses) * self.period)
        return eclipse_times

    def next_primary_ingress_egress_time(self, time, n_eclipses=1):
        """
        Calculate the times of ingress and egress for the next ``n_eclipses``
        primary eclipses after ``time``

        .. warning::
            Barycentric offsets are ignored in the current implementation.

        Parameters
        ----------
        time : `~astropy.time.Time`
            Find the next primary ingress and egress after ``time``
        n_eclipses : int (optional)
            Return the times of eclipse for the next ``n_eclipses`` after
            ``time``. Default is 1.

        Returns
        -------
        primary_eclipses : `~astropy.time.Time` of shape (``n_eclipses``, 2)
            Times of ingress and egress for the next ``n_eclipses`` primary
            eclipses after ``time``
        """
        next_mid_eclipses = self.next_primary_eclipse_time(time, n_eclipses=n_eclipses)
        next_ingresses = next_mid_eclipses - self.duration/2
        next_egresses = next_mid_eclipses + self.duration/2

        ing_egr = np.vstack([next_ingresses.utc.jd, next_egresses.utc.jd]).T

        return Time(ing_egr, format='jd', scale='utc')

    def next_secondary_ingress_egress_time(self, time, n_eclipses=1):
        """
        Calculate the times of ingress and egress for the next ``n_eclipses``
        secondary eclipses after ``time``

        .. warning::
            Barycentric offsets are ignored in the current implementation.

        Parameters
        ----------
        time : `~astropy.time.Time`
            Find the next secondary ingress and egress after ``time``
        n_eclipses : int (optional)
            Return the times of eclipse for the next ``n_eclipses`` after
            ``time``. Default is 1.

        Returns
        -------
        secondary_eclipses : `~astropy.time.Time` of shape (``n_eclipses``, 2)
            Times of ingress and egress for the next ``n_eclipses`` secondary
            eclipses after ``time``.
        """
        next_mid_eclipses = self.next_secondary_eclipse_time(time, n_eclipses=n_eclipses)
        next_ingresses = next_mid_eclipses - self.duration/2
        next_egresses = next_mid_eclipses + self.duration/2

        ing_egr = np.vstack([next_ingresses.utc.jd, next_egresses.utc.jd]).T

        return Time(ing_egr, format='jd', scale='utc')