EclipsingSystem

class astroplan.EclipsingSystem(primary_eclipse_time, orbital_period, duration=None, name=None, eccentricity=None, argument_of_periapsis=None)[source]

Bases: 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).

Parameters
primary_eclipse_timeTime

Time of primary eclipse

orbital_periodQuantity

Orbital period of eclipsing system

durationQuantity (optional)

Duration of eclipse

namestr (optional)

Name of target/event

eccentricityfloat (optional)

Orbital eccentricity. Default is None, which assumes circular orbit (e=0).

argument_of_periapsisfloat (optional)

Argument of periapsis for the eclipsing system, in radians. Default is None, which assumes pi/2.

Methods Summary

in_primary_eclipse(time)

Returns True when time is during a primary eclipse.

in_secondary_eclipse(time)

Returns True when time is during a secondary eclipse

next_primary_eclipse_time(time[, n_eclipses])

Time of the next primary eclipse after time.

next_primary_ingress_egress_time(time[, ...])

Calculate the times of ingress and egress for the next n_eclipses primary eclipses after time

next_secondary_eclipse_time(time[, n_eclipses])

Time of the next secondary eclipse after time.

next_secondary_ingress_egress_time(time[, ...])

Calculate the times of ingress and egress for the next n_eclipses secondary eclipses after time

out_of_eclipse(time)

Returns True when time is not during primary or secondary eclipse.

Methods Documentation

in_primary_eclipse(time)[source]

Returns True when time is during a primary eclipse.

Warning

Barycentric offsets are ignored in the current implementation.

Parameters
timeTime

Time to evaluate

Returns
in_eclipsendarray or bool

True if time is during primary eclipse

in_secondary_eclipse(time)[source]

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 \(\delta t_c\) is given by \(\delta t_c \approx 0.5 \left (\frac{4}{\pi} e \cos{\omega \right)\), where \(e\) is the orbital eccentricity and \(\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
timeTime

Time to evaluate

Returns
in_eclipsendarray or bool

True if time is during secondary eclipse

References

1

Winn (2010) https://arxiv.org/abs/1001.2010

next_primary_eclipse_time(time, n_eclipses=1)[source]

Time of the next primary eclipse after time.

Warning

Barycentric offsets are ignored in the current implementation.

Parameters
timeTime

Find the next primary eclipse after time

n_eclipsesint (optional)

Return the times of eclipse for the next n_eclipses after time. Default is 1.

Returns
primary_eclipsesTime

Times of the next n_eclipses primary eclipses after time

next_primary_ingress_egress_time(time, n_eclipses=1)[source]

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
timeTime

Find the next primary ingress and egress after time

n_eclipsesint (optional)

Return the times of eclipse for the next n_eclipses after time. Default is 1.

Returns
primary_eclipsesTime of shape (n_eclipses, 2)

Times of ingress and egress for the next n_eclipses primary eclipses after time

next_secondary_eclipse_time(time, n_eclipses=1)[source]

Time of the next secondary eclipse after time.

Warning

Barycentric offsets are ignored in the current implementation.

Parameters
timeTime

Find the next secondary eclipse after time

n_eclipsesint (optional)

Return the times of eclipse for the next n_eclipses after time. Default is 1.

Returns
secondary_eclipsesTime

Times of the next n_eclipses secondary eclipses after time

next_secondary_ingress_egress_time(time, n_eclipses=1)[source]

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
timeTime

Find the next secondary ingress and egress after time

n_eclipsesint (optional)

Return the times of eclipse for the next n_eclipses after time. Default is 1.

Returns
secondary_eclipsesTime of shape (n_eclipses, 2)

Times of ingress and egress for the next n_eclipses secondary eclipses after time.

out_of_eclipse(time)[source]

Returns True when time is not during primary or secondary eclipse.

Warning

Barycentric offsets are ignored in the current implementation.

Parameters
timeTime

Time to evaluate

Returns
in_eclipsendarray or bool

True if time is not during primary or secondary eclipse