Phase curves¶
-
pytransit.utils.phasecurves.
doppler_boosting_alpha
(teff: float, flt)[source]¶ The photon weighted bandpass-integrated boosting factor.
Parameters: - teff – Effective temperature of the star [K]
- flt – Passband transmission
Returns: The photon weighted bandpass-integrated boosting factor.
Return type:
-
pytransit.utils.phasecurves.
doppler_boosting_amplitude
(mp: Union[float, numpy.ndarray], ms: Union[float, numpy.ndarray], period: Union[float, numpy.ndarray], alpha: Union[float, numpy.ndarray]) → Union[float, numpy.ndarray][source]¶ The amplitude of the doppler boosting signal.
Calculates the amplitude of the doppler boosting (beaming, reflex doppler effect) signal following the approach described by Loeb & Gaudi in [Loeb2003] . Note that you need to pre-calculate the photon-weighted bandpass-integrated boosting factor (alpha) [Bloemen2010] [Barclay2012] for the star and the instrument using
doppler_boosting_alpha
.Parameters: Returns: Doppler boosting signal amplitude
Return type: float or ndarray
References
[Loeb2003] Loeb, A. & Gaudi, B. S. Periodic Flux Variability of Stars due to the Reflex Doppler Effect Induced by Planetary Companions. Astrophys. J. 588, L117–L120 (2003). [Bloemen2010] Bloemen, S. et al. Kepler observations of the beaming binary KPD 1946+4340. MNRAS 410, (2010). [Barclay2012] Barclay, T. et al. PHOTOMETRICALLY DERIVED MASSES AND RADII OF THE PLANET AND STAR IN THE TrES-2 SYSTEM. AspJ 761, 53 (2012).
-
pytransit.utils.phasecurves.
ellipsoidal_variation_amplitude
(mp: Union[float, numpy.ndarray], ms: Union[float, numpy.ndarray], a: Union[float, numpy.ndarray], i: Union[float, numpy.ndarray], u: Union[float, numpy.ndarray], g: Union[float, numpy.ndarray]) → Union[float, numpy.ndarray][source]¶ The amplitude of the ellipsoidal variation signal.
Calculates the amplitude of the ellipsoidal variation signal following the approach described by Lillo-Box et al. in [Lillo-Box2014], page 11.
Parameters: - mp (float or ndarray) – Planetary mass [MJup]
- ms (float or ndarray) – Stellar mass [MSun]
- a (float or ndarray) – Semi-major axis of the orbit divided by the stellar radius
- i (float or ndarray) – Orbital inclination [rad]
- u (float or ndarray) – Linear limb darkening coefficient
- g (float or ndarray) – Gravity darkening coefficient
Returns: ev_amplitude – The amplitude of the ellipsoidal variation signal
Return type: float or ndarray
References
[Lillo-Box2014] Lillo-Box, J. et al. Kepler-91b: a planet at the end of its life. A&A 562, A109 (2014).
-
pytransit.utils.phasecurves.
ellipsoidal_variation_signal
(f: Union[float, numpy.ndarray], theta: Union[float, numpy.ndarray], e: float) → Union[float, numpy.ndarray][source]¶ Parameters: - f – True anomaly [rad]
- theta – Angle between the line-of-sight and the star-planet direction
- e – Eccentricity
-
pytransit.utils.phasecurves.
emission
(tp: Union[float, numpy.ndarray], tstar: Union[float, numpy.ndarray], k: Union[float, numpy.ndarray], flt) → Union[float, numpy.ndarray][source]¶ Thermal emission from the planet.
Parameters: Returns: - float or ndarray
- References
-
pytransit.utils.phasecurves.
equilibrium_temperature
(tstar: Union[float, numpy.ndarray], a: Union[float, numpy.ndarray], f: Union[float, numpy.ndarray], ab: Union[float, numpy.ndarray]) → Union[float, numpy.ndarray][source]¶ Planetary equilibrium temperature [K].
Parameters: - tstar – Effective stellar temperature [K]
- a – Scaled semi-major axis [Rsun]
- f – Redistribution factor
- ab – Bond albedo
Returns: Teq – Equilibrium temperature [K]
Return type: float or ndarray
-
pytransit.utils.phasecurves.
flux_ratio
(tstar: Union[float, numpy.ndarray], a: Union[float, numpy.ndarray], f: Union[float, numpy.ndarray], ab: Union[float, numpy.ndarray], l: Union[float, numpy.ndarray], r: Union[float, numpy.ndarray] = 1.5, ti: Union[float, numpy.ndarray] = 0) → Union[float, numpy.ndarray][source]¶ Total flux ratio per projected area element.
Parameters: - tstar – Effective stellar temperature [K]
- a – Scaled semi-major axis [Rs]
- f – Redistribution factor
- ab – Bond albedo
- l – Wavelength [m]
- r – Inverse of the phase integral
- ti – Temperature [K]
Returns: fr – Total flux ratio
Return type:
-
pytransit.utils.phasecurves.
planck
(t: Union[float, numpy.ndarray], l: Union[float, numpy.ndarray]) → Union[float, numpy.ndarray][source]¶ Radiance of a black body as a function of wavelength.
Parameters: - t – Black body temperature [K]
- l – Wavelength [m]
Returns: L – Back body radiance [W m^-2 sr^-1]
Return type: float or ndarray
-
pytransit.utils.phasecurves.
reflected_fr
(a: Union[float, numpy.ndarray], ab: Union[float, numpy.ndarray], r: Union[float, numpy.ndarray] = 1.5) → Union[float, numpy.ndarray][source]¶ Reflected flux ratio per projected area element.
Parameters: - a – Scaled semi-major axis [Rsun]
- ab – Bond albedo
- r – Inverse of the phase integral
Returns: fr – Reflected flux ratio
Return type:
-
pytransit.utils.phasecurves.
solve_ab
(fr, tstar, a, f, l, r=1.5, ti=0)[source]¶ Solve the Bond albedo.
Parameters: - fr – Flux ratio [-]
- tstar – Effective stellar temperature [K]
- a – Scaled semi-major axis [Rs]
- A – Bond albedo [-]
- l – Wavelength [m]
- r – Inverse of the phase integral [-]
- ti – Temperature [K]
Returns: A
Return type: Bond albedo
-
pytransit.utils.phasecurves.
solve_redistribution
(fr, tstar, a, ab, l)[source]¶ Solve the redistribution factor.
Parameters: - fr – Flux ratio [-]
- tstar – Effective stellar temperature [K]
- a – Scaled semi-major axis [Rs]
- ab – Bond albedo [-]
- l – Wavelength [m]
- r – Inverse of the phase integral [-]
- Ti (t) – Temperature [K]
Returns: f
Return type: Redistribution factor
-
pytransit.utils.phasecurves.
solve_teq
(fr, tstar, a, ab, l, r=1.5, ti=0)[source]¶ Solve the equilibrium temperature.
Parameters: - fr – Flux ratio
- tstar – Effective stellar temperature [K]
- a – Scaled semi-major axis [Rs]
- ab – Bond albedo
- l – Wavelength [m]
- r – Inverse of the phase integral
- ti – Temperature [K]
Returns: Teq – Equilibrium temperature
Return type: float or ndarray