API/Reference

Index

Light Curves

Transits.PolynomialLimbDarkType
PolynomialLimbDark(u::AbstractVector)

Polynomial limb darkening using analytical integrals. The length of the u vector is equivalent to the order of polynomial used; e.g., [0.2, 0.3] corresponds to quadratic limb darkening.

Mathematical form

\[I(\mu) \propto 1 - u_1(1-\mu) - u_2(1-\mu)^2 - \dots - u_N(1-\mu)^N\]

which is equivalent to the series

\[I(\mu) \propto -\sum_{i=0}^N{u_i(1-\mu)^i}\]

with the definition $u_0 \equiv -1$.

Examples

u = [0.4, 0.26] # quadratic and below is 100% analytical
ld = PolynomialLimbDark(u)
ld(0.1, 0.01)

# output
0.9998787880717668
u2 = vcat(u, ones(12) ./ 12)
ld2 = PolynomialLimbDark(u2)
ld2(0.1, 0.01)

# output
0.9998740059086433

References

Agol, Luger, Foreman-Mackey (2020)

"Analytic Planetary Transit Light Curves and Derivatives for Stars with Polynomial Limb Darkening"

Luger et al. (2019)

"starry: Analytic Occultation Light Curves"

source
Transits.QuadLimbDarkType
QuadLimbDark(u::AbstractVector)

A specialized implementation of PolynomialLimbDark with a maximum of two terms (quadratic form). This has a completely closed-form solution without any numerical integration. This means there are no intermediate allocations and reduced numerical error.

Mathematical form

\[I(\mu) \propto 1 - u_1(1-\mu) - u_2(1-\mu)^2\]

Higher-order terms

Higher-order terms will be ignored; no error will be thrown

Examples

ld = QuadLimbDark(Float64[]) # constant term only

b = [0, 1, 2] # impact parameter
r = 0.01 # radius ratio
ld.(b, r)

# output
3-element Vector{Float64}:
 0.9999
 0.9999501061035608
 1.0
ld = QuadLimbDark([0.4, 0.26]) # max two terms
ld.(b, r)

# output
3-element Vector{Float64}:
 0.9998785437247428
 0.999974726693709
 1.0

References

See references for PolynomialLimbDark

source
Transits.IntegratedLimbDarkType
IntegratedLimbDark(limbdark; N=21, basis=:legendre)
IntegratedLimbDark(u; kwargs...)

Computes the time-averaged flux in the middle of an exposure by wrapping a limb darkening law limbdark with a quadrature scheme. For each time step t, N extra points are super-sampled from t-texp/2 to t+texp/2and the time-averaged flux is calculated via quadrature.

If a set of limb darkening coefficients, u, is provided, a PolynomialLimbDark law will be used by default.

Mathematical form

\[\bar{F}(t) = \frac{1}{\Delta t}\int_{t-\Delta t / 2}^{t+\Delta t / 2}{F(t')dt'}\]

where $F$ is the wrapped limb darkening law and $\Delta t$ is the exposure time.

Quadrature

The integration is approximated via Guassian quadrature

\[\frac{1}{\Delta t} \int{F(t')dt'} \approx \frac12\sum_i^N{w_i * F(\frac{\Delta t}{2}\xi_i + t)}\]

where the weights w_i and nodes ξ_i are defined by the given quadrature rule. The nodes are defined by evaluating orthogonal polynomials N times between -1 and 1. Notice the change of interval required to go from the natural bounds of the orthogonal polynomial basis, -1, 1, to the range defined by the exposure time.

The following bases are available from FastGaussQuadrature.jl. In addition, a function can be passed which calculates nodes, weights = f(N).

  • :legendre - Legendre polynomial base on the open (-1, 1)
  • :radau - Legendre polynomial base on the semi-open [-1, 1) interval
  • :lobatto - Legendre polynomial base on the closed [-1, 1] interval
source
Transits.SecondaryLimbDarkType
SecondaryLimbDark(primary::AbstractLimbDark,
                  secondary::AbstractLimbDark; 
                  brightness_ratio=1)
SecondaryLimbDark(u_p::AbstractVector, u_s=u_p; kwargs...)

Compose two limb darkening laws together to add a secondary eclipse. If vectors of coefficients are provided, laws will automatically be constructed using PolynomialLimbDark. The surface brightness ratio is given in terms of the host; e.g., if the companion is half as bright as the host, the ratio would be 0.5.

Interface

SecondaryLimbDark only works with an orbit, since the companion's reference frame needs to be calculated. This means you can't call it using an impact parameter like ld(b, r) directly.

Mathematical form

\[f(t, r) = \frac{2f_p(t, r) + \eta r^2 f_s(t', r')}{1 + f_p(t, r) + \eta r^2 f_s(t', r')}\]

where $f_p$ is to the primary flux, $f_s$ is to the secondary flux, and $\eta$ is the surface brightness ratio. $t'$ and $r'$ correspond to the time and radius ratio from the companion's reference frame.

Examples

# equal size and limb darkening
r = 1.0
u = [0.4, 0.26]
# companion is 1/10 as bright
brightness_ratio = 0.1
ld = SecondaryLimbDark(u; brightness_ratio)
orbit = SimpleOrbit(period=2, duration=0.5)
fp = ld(orbit, 0, r) # primary egress
fs = ld(orbit, 1, r) # secondary egress

fp ≈ brightness_ratio * fs

# output
true
source
Transits.computeFunction
compute(::AbstractLimbDark, b, r; kwargs...)

Compute the relative flux for the given impact parameter b and radius ratio r. The impact parameter is unitless. The radius ratio is given in terms of the host; e.g., if the companion is half the size of the host, r=0.5.

source
Transits.computeMethod
compute(::AbstractLimbDark, orbit::AbstractOrbit, t, r)

Compute the relative flux by calculating the impact parameter at time t from the given orbit. The time needs to be compatible with the period of the orbit, nominally in days.

Examples

julia> ld = PolynomialLimbDark([0.4, 0.26]);

julia> orbit = SimpleOrbit(period=3, duration=1);

julia> ld(orbit, 0, 0.1) # primary egress
0.9878664434953113

julia> ld(orbit, 0.1, 0.1) # 0.1 d
0.9879670695533511

this works effortlessly with libraries like Unitful.jl

julia> using Unitful

julia> orbit = SimpleOrbit(period=3u"d", duration=3u"hr");

julia> ld(orbit, 0u"d", 0.1)
0.9878664434953113
source

Gradients

Gradients and jacobians are integrated directly into ChainRules.jl via frules and rrules. For most users, this just means using AD libraries like ForwardDiff.jl and Zygote.jl is effortless and fast.

using Transits
using Zygote

lightcurve(X) = compute(PolynomialLimbDark(X[3:end]), X[1], X[2])
grad(X) = lightcurve'(X) # Zygote gradient
grad([0.1, 0.1, 0.4, 0.26])

# output
4-element Vector{Float64}:
  0.0004972185834858653
 -0.2419262730830416
 -0.0048107583897073185
 -0.0024501564976671724

To help demonstrate the logic behind these chain rules, here we derive a simple gradient function manually.

using ChainRules

u_n = [0.4, 0.26]
μ = 0.1
ror = 0.1
X0 = [μ, ror, u_n...]

function gradr(X)
    ld, ld_pullback = rrule(PolynomialLimbDark, X[3:end])
    f, f_pullback = rrule(compute, ld, X[1], X[2])

    f̄ = one(eltype(X))
    _, l̄d, b̄, r̄ = f_pullback(f̄)
    _, ū_n = ld_pullback(l̄d)
    return [b̄, r̄, ū_n...]
end

gradr([0.1, 0.1, 0.4, 0.26])

# output
4-element Vector{Float64}:
  0.0004972185834858653
 -0.2419262730830416
 -0.0048107583897073185
 -0.0024501564976671724

For the most granular support for gradients and jacobians, peer into the depths of polynomial/poly-grad.jl and polynomial/quad-grad.jl. These functions are not part of the public API and are not guaranteed any stability according to semantic versioning.

Orbits

Transits.Orbits.SimpleOrbitType
SimpleOrbit(; period, duration, t0=0, b=0.0)

Circular orbit parameterized by the basic observables of a transiting system.

Parameters

  • period - The orbital period of the planets, nominally in days
  • duration - The duration of the transit, similar units as period.
  • t0 - The midpoint time of the reference transit, similar units as period
  • b - The impact parameter of the orbit, unitless
source
Transits.Orbits.KeplerianOrbitType
KeplerianOrbit(; kwargs...)

Keplerian orbit parameterized by the basic observables of a transiting 2-body system.

Parameters

  • period/P – The orbital period of the planet [d].
  • t0/t_0 – The midpoint time of the reference transit [d].
  • tp/t_p – The time of periastron [d].
  • duration/τ/T – The transit duration [d].
  • a – The semi-major axis [R⊙].
  • aR_star/aRs – The ratio of the semi-major axis to star radius.
  • R_planet/Rp – The radius of the planet [R⊙].
  • R_star/Rs – The radius of the star [R⊙].
  • rho_star/ρ_star – The spherical star density [M⊙/R⊙³].
  • r/RpRs – The ratio of the planet radius to star radius.
  • b – The impact parameter, bounded between 0 ≤ b ≤ 1.
  • ecc/e – The eccentricity of the closed orbit, bounded between 0 ≤ ecc < 1.
  • incl – The inclination of the orbital plane relative to the axis perpendicular to the reference plane [rad]
  • omega/ω – The argument of periapsis [rad].
  • cos_omega/cos_ω – The cosine of the argument of periapsis.
  • sin_omega/sin_ω – The sine of the argument of periapsis.
  • Omega/Ω – The longitude of the ascending node [rad].
  • M_planet/Mp – The mass of the planet [M⊙].
  • M_star/Ms – The mass of the star [M⊙].

Valid combinations

The following flowchart can be used to determine which parameters can define a KeplerianOrbit:

  1. The period or a must be given. If both given, then neither M_star or rho_star can be defined because the stellar density is now implied.
  2. Only incl or b can be given.
  3. If ecc is given, then omega must also be given.
  4. If no stellar parameters are given, the central body is assumed to be the Sun. If only rho_star is given, then R_star is defined to be 1 solar radius. Otherwise, at most two of M_star, R_star, and rho_star can be given.
  5. Either t0 or tp must be given, but not both.
source
Transits.Orbits.relative_positionFunction
relative_position(::AbstractOrbit, t)

The relative position, [x, y, z], of the companion compared to the host at time t. In other words, this is the vector pointing from the host to the companion along the line of sight. Nominally, the units of this distance are relative to the host's radius. For example, a distance of 2 is 2 stellar radii.

source

Distributions

Transits.Kipping13Type
Kipping13()

A non-informative prior for two-parameter limb-darkening coefficients using triangular sampling (Kipping 2013).

Examples

julia> using Random; rng = Random.seed!(10);

julia> rand(rng, Kipping13())
2-element Vector{Float64}:
 0.24716310305467298
 0.08836997249298882

julia> rand(rng, Kipping13(), 5)
2×5 Matrix{Float64}:
 0.0664907  0.124817   1.00732    0.806902  0.74165
 0.520411   0.222718  -0.389412  -0.314755  0.0768429

References

Kipping (2013)

"Efficient, uninformative sampling of limb darkening coefficients for two-parameter laws"

source