ravest.model
Radial velocity models for planets, stars, and trends.
Provides classes for modelling radial velocity signals from planetary orbits and stellar trends (constant, linear, quadratic).
Classes
Planet defined by its orbital parameters. |
|
Instrument-specific parameters for RV observations. |
|
System-wide trend in the radial velocity of the star. |
|
Star with orbiting planet(s). |
Functions
|
Solve Kepler's equation for a single mean anomaly value. |
|
Compute cos(f) and sin(f) of the true anomaly from the eccentric anomaly. |
|
Compute radial velocity from the true anomaly of a single observation. |
|
Solve Kepler's equation and compute radial velocity for an array of times. |
|
Compute radial velocity from mean anomaly, dispatching by eccentricity. |
|
Calculate the minimum mass of the planet. |
|
Fold time series to orbital phase and return sorted arrays. |
Module Contents
- ravest.model._solve_kepler(Mi: float, e: float) tuple[float, float][source]
Solve Kepler’s equation for a single mean anomaly value.
Kepler’s equation \(E - e\sin E = M\) must be solved iteratively for the eccentric anomaly \(E\). Uses Halley’s method (cubic convergence) with the same tolerance and iteration limits as scipy’s newton solver.
- Parameters:
Mi (float) – Mean anomaly at a single time (radians).
e (float) – Eccentricity of the orbit, 0 < e < 1 (dimensionless).
- Returns:
cos_E (float) – Cosine of the eccentric anomaly.
sin_E (float) – Sine of the eccentric anomaly.
- ravest.model._true_anomaly(cos_E: float, sin_E: float, e: float, sqrt_1me2: float) tuple[float, float][source]
Compute cos(f) and sin(f) of the true anomaly from the eccentric anomaly.
Calculates the true anomaly directly from cos(E) and sin(E) rather than via the arctan formula, avoiding a trig call per element.
- Parameters:
cos_E (float) – Cosine of the eccentric anomaly.
sin_E (float) – Sine of the eccentric anomaly.
e (float) – Eccentricity of the orbit (dimensionless).
sqrt_1me2 (float) – Precomputed \(\sqrt{1 - e^2}\).
- Returns:
cos_f (float) – Cosine of the true anomaly.
sin_f (float) – Sine of the true anomaly.
Notes
Rather than computing \(f\) as an angle via the arctan formula \(f = 2\arctan\left(\sqrt{\frac{1+e}{1-e}} \tan\frac{E}{2}\right)\), we obtain \(\cos f\) and \(\sin f\) directly from \(\cos E\) and \(\sin E\):
\[ \begin{align}\begin{aligned}\cos f = \frac{\cos E - e}{1 - e\cos E}\\\sin f = \frac{\sqrt{1 - e^2} \sin E}{1 - e\cos E}\end{aligned}\end{align} \]The \(\cos f\) identity is a standard result from the geometry of the auxiliary circle (Perryman, Exoplanet Handbook, eq. 2.6). \(\sin f\) follows from \(\sin^2 f = 1 - \cos^2 f\), which simplifies to \((1 - e^2)\sin^2 E / (1 - e\cos E)^2\).
Since \(\cos E\) and \(\sin E\) are already available from the Kepler solver, this avoids any additional trig calls.
- ravest.model._radial_velocity_from_f(cos_f: float, sin_f: float, K: float, cos_w: float, sin_w: float, e_cos_w: float) float[source]
Compute radial velocity from the true anomaly of a single observation.
Uses the trig identity \(\cos(f + \omega) = \cos f \cos \omega - \sin f \sin \omega\) to avoid computing \(f\) explicitly.
- Parameters:
cos_f (float) – Cosine of the true anomaly.
sin_f (float) – Sine of the true anomaly.
K (float) – Semi-amplitude of the stellar radial velocity (m/s).
cos_w (float) – Precomputed cos(w).
sin_w (float) – Precomputed sin(w).
e_cos_w (float) – Precomputed e*cos(w).
- Returns:
Radial velocity (m/s).
- Return type:
float
Notes
The RV equation is:
\[v_r = K \left[ \cos(f + \omega) + e\cos\omega \right]\]Rather than computing the angle \(f\) and evaluating \(\cos(f + \omega)\), we expand via the addition identity \(\cos(f + \omega) = \cos f \cos \omega - \sin f \sin \omega\), using \(\cos f\) and \(\sin f\) from
_true_anomalydirectly. \(\cos\omega\), \(\sin\omega\), and \(e\cos\omega\) are constant across all observations and precomputed by the caller.
- ravest.model._njit_kepler_rv(M: numpy.ndarray, e: float, K: float, w: float) numpy.ndarray[source]
Solve Kepler’s equation and compute radial velocity for an array of times.
Loops over mean anomaly values element-by-element, calling the scalar Kepler solver, true anomaly, and RV functions for each. Numba inlines these calls at compile time, so performance is identical to a single monolithic loop while keeping each step readable.
- Parameters:
M (np.ndarray) – Mean anomaly at each time (radians).
e (float) – Eccentricity of the orbit, 0 < e < 1 (dimensionless).
K (float) – Semi-amplitude of the stellar radial velocity (m/s).
w (float) – Argument of periapsis of the star (radians).
- Returns:
Radial velocity of the star due to the planet at each time (m/s).
- Return type:
np.ndarray
- ravest.model._compute_rv(M: numpy.ndarray, e: float, K: float, w: float) numpy.ndarray[source]
Compute radial velocity from mean anomaly, dispatching by eccentricity.
For circular orbits (e=0) the RV simplifies to K*cos(M + w) and can be computed directly with numpy. For eccentric orbits, delegates to the numba-compiled solver.
- Parameters:
M (np.ndarray) – Mean anomaly at each time (radians).
e (float) – Eccentricity of the orbit, 0 <= e < 1 (dimensionless).
K (float) – Semi-amplitude of the stellar radial velocity (m/s).
w (float) – Argument of periapsis of the star (radians).
- Returns:
Radial velocity of the star due to the planet at each time (m/s).
- Return type:
np.ndarray
- class ravest.model.Planet(letter: str, parameterisation: ravest.param.Parameterisation, params: dict[str, float])[source]
Planet defined by its orbital parameters.
- Parameters:
letter (str) – The label of the planet, e.g. “b”, “c”. Must be a single letter.
parameterisation (parameterisation) – The set of planetary parameters used to define the planet.
params (dict) – The orbital parameters, matching the parameterisation.
- letter
- parameterisation
- params
- _rvparams
- _calculate_mean_motion(period: float) float[source]
Calculate mean motion (mean angular rate of orbit in radians/day).
This is the mean angular rate calculated by dividing the orbital period by 2*pi radians.
- Parameters:
period (float) – The orbital period in days.
- Returns:
The mean motion, the average angular rate of orbit (radians/day).
- Return type:
float
- _calculate_mean_anomaly(t: numpy.ndarray, n: float, time_peri: float) numpy.ndarray[source]
Calculate mean anomaly (radians).
For an eccentric orbit with period P, mean anomaly is a fictitious angle that increases linearly with time
tfor an angular raten, the mean motion. It is the angle from the point of periapsis passage that would have been swept in a circular orbit of periodPwith fixed angular raten.- Parameters:
n (float) – The mean motion (radians/day).
t (float) – The time to evaluate the mean anomaly at (day).
time_peri (float) – The time of periapsis of the star’s orbit (day).
- Returns:
The mean anomaly (radians).
- Return type:
float
- radial_velocity(t: numpy.ndarray) numpy.ndarray[source]
Calculate radial velocity of the star due to the planet, at time t.
Uses a numba-compiled Kepler solver that solves Kepler’s equation and computes the radial velocity in a single pass.
- Parameters:
t (np.ndarray) – The time to calculate the radial velocity at (day).
- Returns:
Radial velocity of the reflex motion of star due to the planet (m/s).
- Return type:
np.ndarray
- mpsini(mass_star: float, unit: str = 'kg') float[source]
Calculate the minimum mass of the planet.
- Parameters:
mass_star (float) – The mass of the star in solar masses.
unit (str) – The unit to return the planetary minimum mass in. Options are “kg”, “M_earth”, “M_jupiter”.
- Returns:
The minimum mass of the planet (solar masses).
- Return type:
float
- class ravest.model.Instrument(name: str, g: float, jit: float)[source]
Instrument-specific parameters for RV observations.
Represents the measurement characteristics of a specific instrument/telescope, including its velocity offset (gamma) and additional noise (jitter).
- Parameters:
name (str) – The name/identifier of the instrument (e.g., “HARPS”). Must match the labels used in the data’s instrument column.
g (float) – Gamma offset - the constant RV offset for this instrument [m/s].
jit (float) – Jitter - additional noise added in quadrature to uncertainties [m/s]. Must be >= 0.
Examples
>>> harps = Instrument("HARPS", g=5.0, jit=2.0) >>> print(harps) Instrument HARPS: γ=5.0 m/s, jitter=2.0 m/s
>>> hires = Instrument("HIRES", g=-3.6, jit=1.5) >>> hires.g -3.6
- name
- g
- jit
- class ravest.model.Trend(t0: float, params: dict[str, float])[source]
System-wide trend in the radial velocity of the star.
Represents long-term velocity trends that apply across all instruments, such as acceleration from a distant companion. The constant offset (gamma) is now handled per-instrument via the Instrument class.
- Parameters:
t0 (float) – The reference zero-point time for the linear and quadratic trend. Recommended to be the mean of the input times.
params (dict) – The parameters of the trend: linear and quadratic components. These must be named “gd” and “gdd” respectively (gamma-dot, gamma-dot-dot). These are in units of m/s/day and m/s/day^2 respectively.
- Returns:
The radial velocity of the star due to the trend (m/s).
- Return type:
float
Notes
The radial velocity contribution from the trend is calculated as the sum of the linear and quadratic components: gd*(t-t0) and gdd*((t-t0)**2).
The constant offset (gamma) is now instrument-specific and handled via the Instrument class. Use star.gamma_offsets(instrument) to get per-instrument offsets.
In general the trend is used to account for long-term effects such as acceleration from a distant companion. If you see a strong linear or quadratic trend in the data, it is worth investigating.
Examples
>>> trend = Trend(t0=2458000.0, params={"gd": 0.001, "gdd": 0.0}) >>> rv_trend = trend.radial_velocity(times)
- gammadot
- gammadotdot
- class ravest.model.Star(name: str, mass: float)[source]
Star with orbiting planet(s).
- Parameters:
name (str) – The name of the star system.
mass (float) – The mass of the star (solar masses)
- planets
The dict storing the Planet objects, the key is the planet.letter attribute.
- Type:
dict
- num_planets
The number of Planet objects stored in the Star object.
- Type:
int
- name
- mass
- planets
- instruments
- num_planets = 0
- add_planet(planet: Planet) None[source]
Store Planet object in planets dict with the key Planet.letter.
Planets cannot share letters; if two planets have the same letter then the second one will overwrite the first.
- Parameters:
planet (Planet) – A ravest.model.Planet object
- add_trend(trend: Trend) None[source]
Store Trend object in trend attribute.
- Parameters:
trend (Trend) – A ravest.model.Trend object
- add_instrument(instrument: Instrument) None[source]
Store Instrument object in instruments dict with key instrument.name.
Instruments cannot share names; if two instruments have the same name then the second one will overwrite the first.
- Parameters:
instrument (Instrument) – A ravest.model.Instrument object
- gamma_offsets(instrument: numpy.ndarray) numpy.ndarray[source]
Return gamma offset for each observation based on instrument column.
- Parameters:
instrument (np.ndarray) – Array of instrument names for each observation.
- Returns:
Array of gamma offsets (m/s) for each observation.
- Return type:
np.ndarray
- Raises:
ValueError – If an instrument in the data is not found in the Star’s instruments.
- jitter_values(instrument: numpy.ndarray) numpy.ndarray[source]
Return jitter value for each observation based on instrument column.
- Parameters:
instrument (np.ndarray) – Array of instrument names for each observation.
- Returns:
Array of jitter values (m/s) for each observation.
- Return type:
np.ndarray
- radial_velocity(t: numpy.ndarray) numpy.ndarray[source]
Calculate the radial velocity of the star at time
tdue to the planets and trend.Calculate the radial velocity of the star at time
tdue to the planets and the trend (constant velocity, linear, and quadratic.) This is a linear sum of the RVs of each of the Planet objects stored in the Star, and the RV of the Trend object.- Parameters:
t (float) – The time at which to calculate the radial velocity (day)
- Returns:
The radial velocity at time t (m/s)
- Return type:
float
- mpsini(planet_letter: str, unit: str = 'kg') float[source]
Calculate the minimum mass of a planet in the star system.
- Parameters:
planet_letter (str) – The letter of the planet to calculate the minimum mass of.
unit (str) – The unit to return the planetary minimum mass in. Options are “kg”, “M_earth”, “M_jupiter”.
- Returns:
The minimum mass of the planet.
- Return type:
float
- phase_plot(t: numpy.ndarray, ydata: numpy.ndarray, yerr: numpy.ndarray, instrument: numpy.ndarray) None[source]
Generate a phase plot for each planet.
Given RV
ydataat timetwith errorbarsyerr, generates a phase plot for each planet. Data is coloured by instrument and gamma offsets are subtracted before plotting.- Parameters:
t (np.ndarray) – The time of the observations
ydata(day)ydata (np.ndarray) – The observed radial velocity at time
t(m/s)yerr (np.ndarray) – The measurement error of the datapoint
ydata(m/s)instrument (np.ndarray) – The instrument/telescope name for each observation.
- ravest.model.calculate_mpsini(mass_star: float, period: float, semi_amplitude: float, eccentricity: float, unit: str = 'kg') float[source]
Calculate the minimum mass of the planet.
- Parameters:
mass_star (float) – The mass of the star in solar masses.
period (float) – The orbital period of the planet in days.
semi_amplitude (float) – The semi-amplitude of the radial velocity of the star in m/s.
eccentricity (float) – The eccentricity of the orbit, 0 <= e < 1 (dimensionless).
unit (str) – The unit to return the planetary minimum mass in. Options are “kg”, “M_earth”, “M_jupiter”.
- Returns:
The minimum mass mpsini of the planet.
- Return type:
float
- Raises:
ValueError – If the unit is not one of “kg”, “M_earth”, “M_jupiter”.
- ravest.model.fold_time_series(times: numpy.ndarray, period: float, t_ref: float) tuple[numpy.ndarray, numpy.ndarray][source]
Fold time series to orbital phase and return sorted arrays.
Converts times to orbital phases in the range [-0.5, 0.5] and returns both the sorted phases and the indices needed to sort other arrays consistently.
- Parameters:
times (np.ndarray) – Time values to fold
period (float) – Orbital period
t_ref (float) – Reference time (usually Tc - time of transit/conjunction)
- Returns:
phases_sorted (np.ndarray) – Phase-folded times in range [-0.5, 0.5], sorted ascending
sort_indices (np.ndarray) – Indices that sort the original times by phase
Examples
>>> times = np.array([0, 1, 2, 3, 4]) >>> period = 2.0 >>> t_ref = 0.5 # Reference time >>> phases, indices = fold_time_series(times, period, t_ref) >>> # phases will be sorted from -0.5 to 0.5 >>> # indices can be used to sort other arrays consistently