funcs

Functions relating observable properties of binary stars and exoplanet systems to their fundamental properties, and vice versa. Also functions related to Keplerian orbits.

Parameters

Functions are defined in terms of the following parameters. 1

  • a - orbital semi-major axis in solar radii = a_1 + a_2

  • P - orbital period in mean solar days

  • Mass - total system mass in solar masses, Mass = m_1 + m_2

  • ecc - orbital eccentricity

  • omdeg - longitude of periastron of star’s orbit, omega, in _degrees_

  • sini - sine of the orbital inclination

  • K - 2.pi.a.sini/(P.sqrt(1-e^2)) = K_1 + K_2

  • K_1, K_2 - orbital semi-amplitudes in km/s

  • q - mass ratio = m_2/m_1 = K_1/K_2 = a_1/a_2

  • f_m - mass function = m_2^3.sini^3/(m_1+m_2)^2 in solar masses

    = K_1^3.P/(2.pi.G).(1-e^2)^(3/2)

  • r_1 - radius of star 1 in units of the semi-major axis, r_1 = R_*/a

  • rhostar - mean stellar density = 3.pi/(GP^2(1+q)r_1^3)

  • rstar - host star radius/semi-major axis, rstar = R_*/a

  • k - planet/star radius ratio, k = R_planet/R_star

  • tzero - time of mid-transit (minimum on-sky star-planet separation).

  • b - impact parameter, b = a.cos(i)/R_star

1

Hilditch, R.W., An Introduction to Close Binary Stars, CUP 2001.

Functions

pycheops.funcs.a_rsun(P, Mass)

Semi-major axis in solar radii

Parameters
  • P – orbital period in mean solar days

  • Mass – total mass in solar masses, M

Returns

a = (G.M.P^2/(4.pi^2))^(1/3) in solar radii

pycheops.funcs.f_m(P, K, ecc=0)

Mass function in solar masses

Parameters
  • P – orbital period in mean solar days

  • K – semi-amplitude of the spectroscopic orbit in km/s

  • ecc – orbital eccentricity

Returns

f_m = m_2^3.sini^3/(m_1+m_2)^2 in solar masses

pycheops.funcs.m1sin3i(P, K_1, K_2, ecc=0)

Reduced mass of star 1 in solar masses

Parameters
  • K_1 – semi-amplitude of star 1 in km/s

  • K_2 – semi-amplitude of star 2 in km/s

  • P – orbital period in mean solar days

  • ecc – orbital eccentricity

Returns

m_1.sini^3 in solar masses

pycheops.funcs.m2sin3i(P, K_1, K_2, ecc=0)

Reduced mass of star 2 in solar masses

Parameters
  • K_1 – semi-amplitude of star 1 in km/s

  • K_2 – semi-amplitude of star 2 in km/s

  • P – orbital period in mean solar days

  • ecc – orbital eccentricity

Returns

m_2.sini^3 in solar masses

pycheops.funcs.asini(K, P, ecc=0)

a.sini in solar radii

Parameters
  • K – semi-amplitude of the spectroscopic orbit in km/s

  • P – orbital period in mean solar days

Returns

a.sin(i) in solar radii

pycheops.funcs.rhostar(r_1, P, q=0)

Mean stellar density from scaled stellar radius.

Parameters
  • r_1 – radius of star in units of the semi-major axis, r_1 = R_*/a

  • P – orbital period in mean solar days

  • q – mass ratio, m_2/m_1

Returns

Mean stellar density in solar units

pycheops.funcs.K_kms(m_1, m_2, P, sini, ecc)
Semi-amplitudes of the spectroscopic orbits in km/s
  • K = 2.pi.a.sini/(P.sqrt(1-ecc^2))

  • K_1 = K * m_2/(m_1+m_2)

  • K_2 = K * m_1/(m_1+m_2)

Parameters
  • m_1 – mass of star 1 in solar masses

  • m_2 – mass of star 2 in solar masses

  • P – orbital period in mean solar days

  • sini – sine of the orbital inclination

  • ecc – orbital eccentrcity

Returns

K_1, K_2 – semi-amplitudes in km/s

pycheops.funcs.m_comp(f_m, m_1, sini)

Companion mass in solar masses given mass function and stellar mass

Parameters
  • f_m – = K_1^3.P/(2.pi.G).(1-ecc^2)^(3/2) in solar masses

  • m_1 – mass of star 1 in solar masses

  • sini – sine of orbital inclination

Returns

m_2 = mass of companion to star 1 in solar masses

pycheops.funcs.transit_width(r, k, b, P=1)

Total transit duration for a circular orbit.

See equation (3) from Seager and Malen-Ornelas, 2003ApJ…585.1038S.

Parameters
  • r – R_star/a

  • k – R_planet/R_star

  • b – impact parameter = a.cos(i)/R_star

  • P – orbital period (optional, default P=1)

Returns

Total transit duration in the same units as P.

pycheops.funcs.t2z(t, tzero, P, sini, rstar, ecc=0, omdeg=90, returnMask=False)

Calculate star-planet separation relative to scaled stellar radius, z

Optionally, return a flag/mask to indicate cases where the planet is further from the observer than the star, i.e., whether phases with z<1 are transits (mask==True) or eclipses (mask==False)

Parameters
  • t – time of observation (scalar or array)

  • tzero – time of inferior conjunction, i.e., mid-transit

  • P – orbital period

  • sini – sine of orbital inclination

  • rstar – scaled stellar radius, R_star/a

  • ecc – eccentricity (optional, default=0)

  • omdeg – longitude of periastron in degrees (optional, default=90)

  • returnFlag – return a flag to distinguish transits from eclipses.

N.B. omdeg is the longitude of periastron for the star’s orbit

Returns

z [, mask]

Example

>>> from pycheops.funcs import t2z
>>> from numpy import linspace
>>> import matplotlib.pyplot as plt
>>> t = linspace(0,1,1000)
>>> sini = 0.999
>>> rstar = 0.1
>>> plt.plot(t, t2z(t,0,1,sini,rstar))
>>> plt.xlim(0,1)
>>> plt.ylim(0,12)
>>> ecc = 0.1
>>> for omdeg in (0, 90, 180, 270):
>>>     plt.plot(t, t2z(t,0,1,sini,rstar,ecc,omdeg))
>>> plt.show()
pycheops.funcs.tzero2tperi(tzero, P, sini, ecc, omdeg)

Calculate time of periastron from time of mid-transit

Uses the method by Lacy, 1992AJ….104.2213L

Parameters
  • tzero – times of mid-transit

  • P – orbital period

  • sini – sine of orbital inclination

  • ecc – eccentricity

  • omdeg – longitude of periastron in degrees

Returns

time of periastron prior to tzero

Example
>>> from pycheops.funcs import tzero2tperi
>>> tzero = 54321.6789
>>> P = 1.23456
>>> sini = 0.987
>>> ecc = 0.123
>>> omdeg = 89.01
>>> print(tzero2tperi(tzero,P,sini,ecc,omdeg))
54321.6762764
pycheops.funcs.vrad(t, tzero, P, K, ecc=0, omdeg=90, sini=1, primary=True)

Calculate radial velocity, V_r, for body in a Keplerian orbit

Parameters
  • t – array of input times

  • tzero – time of inferior conjunction, i.e., mid-transit

  • P – orbital period

  • K – radial velocity semi-amplitude

  • ecc – eccentricity (optional, default=0)

  • omdeg – longitude of periastron in degrees (optional, default=90)

  • sini – sine of orbital inclination (to convert tzero to t_peri)

  • primary – if false calculate V_r for companion

Returns

V_r in same units as K relative to the barycentre of the binary

pycheops.funcs.xyz_planet(t, tzero, P, sini, ecc=0, omdeg=90)

Position of the planet in Cartesian coordinates.

The position of the ascending node is taken to be Omega=0 and the semi-major axis is taken to be a=1.

Parameters
  • t – time of observation (scalar or array)

  • tzero – time of inferior conjunction, i.e., mid-transit

  • P – orbital period

  • sini – sine of orbital inclination

  • ecc – eccentricity (optional, default=0)

  • omdeg – longitude of periastron in degrees (optional, default=90)

N.B. omdeg is the longitude of periastron for the star’s orbit

Returns

(x, y, z)

Example

>>> from pycheops.funcs import phase_angle
>>> from numpy import linspace
>>> import matplotlib.pyplot as plt
>>> t = linspace(0,1,1000)
>>> sini = 0.9
>>> ecc = 0.1
>>> omdeg = 90
>>> x, y, z = xyz_planet(t, 0, 1, sini, ecc, omdeg)
>>> plt.plot(x, y)
>>> plt.plot(x, z)
>>> plt.show()