Utility Functions
- radvel.utils.Msini(K: float | ndarray[tuple[Any, ...], dtype[float]], P: float | ndarray[tuple[Any, ...], dtype[float]], Mstar: float | ndarray[tuple[Any, ...], dtype[float]], e: float | ndarray[tuple[Any, ...], dtype[float]], Msini_units: str = 'earth') float | ndarray[tuple[Any, ...], dtype[float]][source]
Calculate Msini
Calculate Msini for a given K, P, stellar mass, and e
- Parameters:
array (K (float or) – Doppler semi-amplitude [m/s]
P (float or array) – Orbital period [days]
Mstar (float or array) – Mass of star [Msun]
e (float or array) – eccentricity
Msini_units (Optional[str]) – Units of Msini {‘earth’,’jupiter’} default: ‘earth’
- Returns:
Msini [units = Msini_units]
- Return type:
float or array
- radvel.utils.bintels(t: ndarray[tuple[Any, ...], dtype[float]], vel: ndarray[float], err: ndarray[float], telvec: ndarray[str], binsize: float = 0.5) tuple[ndarray[tuple[Any, ...], dtype[float]], ndarray[tuple[Any, ...], dtype[float]], ndarray[tuple[Any, ...], dtype[float]], ndarray[tuple[Any, ...], dtype[str]]][source]
Bin velocities by instrument
Bin RV data with bins of with binsize in the units of t. Will not bin data from different telescopes together since there may be offsets between them.
- Parameters:
t (array) – array of timestamps
vel (array) – array of velocities
err (array) – array of velocity uncertainties
telvec (array) – array of strings corresponding to the instrument name for each velocity
binsize (float) – (optional) width of bin in units of t (default=1/2.)
- Returns:
(bin centers, binned measurements, binned uncertainties, binned instrument codes)
- Return type:
tuple
- radvel.utils.date2jd(date: datetime) float[source]
Convert datetime object to JD”
- Parameters:
date (datetime.datetime) – date to convert
- Returns:
Julian date
- Return type:
float
- radvel.utils.density(mass: float, radius: float, MR_units: str = 'earth') float[source]
Compute density from mass and radius
- Parameters:
mass (float) – mass [MR_units]
radius (float) – radius [MR_units]
MR_units (string) – (optional) units of mass and radius. Must be ‘earth’, or ‘jupiter’ (default ‘earth’).
- Returns:
density in g/cc
- Return type:
float
- radvel.utils.draw_models_from_chain(mod: ~radvel.model.GeneralRVModel, chain: <Mock name='mock.DataFrame' id='135641087121584'>, t: ~numpy.ndarray[tuple[~typing.Any, ...], ~numpy.dtype[float]], nsamples: int = 50) ndarray[tuple[Any, ...], dtype[float]][source]
Draw Models from Chain
Given an MCMC chain of parameters, draw representative parameters and synthesize models.
- Parameters:
mod (radvel.RVmodel) – RV model
chain (DataFrame) – pandas DataFrame with different values from MCMC chain
t (array) – time range over which to synthesize models
nsamples (int) – number of draws
- Returns:
2D array with the different models as different rows
- Return type:
array
- radvel.utils.fastbin(x: ndarray[tuple[Any, ...], dtype[float]], y: ndarray[tuple[Any, ...], dtype[float]], nbins: int = 30) tuple[ndarray[tuple[Any, ...], dtype[float]], ndarray[tuple[Any, ...], dtype[float]], ndarray[tuple[Any, ...], dtype[float]]][source]
Fast binning
Fast binning function for equally spaced data
- Parameters:
x (array) – independent variable
y (array) – dependent variable
nbins (int) – number of bins
- Returns:
(bin centers, binned measurements, binned uncertainties)
- Return type:
tuple
- radvel.utils.geterr(vec: ndarray[tuple[Any, ...], dtype[float]], angular: bool = False) tuple[float, float, float][source]
Calculate median, 15.9, and 84.1 percentile values for a given vector.
- Parameters:
vec (array) – vector, usually an MCMC chain for one parameter
angular (bool [optioanl]) – Is this an angular parameter? if True vec should be in radians. This will perform some checks to ensure proper boundary wrapping.
- Returns:
50, 15.9 and 84.1 percentiles
- Return type:
tuple
- radvel.utils.initialize_posterior(config_file: str, decorr: bool = False) tuple[ModuleType, Posterior][source]
Initialize Posterior object from a Python setup file.
Parse a setup file and initialize the RVModel, Likelihood, Posterior and priors.
- Parameters:
config_file (string) – path to config file
decorr (bool) – (optional) decorrelate RVs against columns defined in the decorr_vars list
- Returns:
(object representation of config file, radvel.Posterior object)
- Return type:
tuple
- radvel.utils.initialize_posterior_from_dict(config, decorr=False)[source]
Initialize Posterior object from an in-memory dict.
Mirrors the namespace of a Python setup file but accepts JSON-friendly types so callers (e.g. an HTTP API) can build a posterior without writing a file to disk first. Each top-level field may either contain a ready-made radvel object (e.g. a
radvel.Parametersinstance forparams) or a JSON-style representation that this function translates via the helpers below.- Parameters:
config (dict) – mapping with keys
starname,nplanets,instnames,fitting_basis,params,data,priors, optionallybjd0,time_base,planet_letters,stellar,planet,decorr_vars,hnames,kernel_name,any_basis,anybasis_params.decorr (bool) – same meaning as
initialize_posterior().
- Returns:
(SimpleNamespace mirroring the legacy
Pmodule, radvel.Posterior)- Return type:
tuple
- radvel.utils.jd2date(jd: float) datetime[source]
Convert JD to datetime.datetime object
- Parameters:
jd (float) – Julian date
- Returns:
calendar date
- Return type:
datetime.datetime
- radvel.utils.load_module_from_file(module_name: str, module_path: str) ModuleType[source]
Loads a python module from the path of the corresponding file.
- Parameters:
module_name (str) – namespace where the python module will be loaded, e.g.
foo.barmodule_path (str) – path of the python file containing the module
- Returns:
A valid module object
- Raises:
ImportError – when the module can’t be loaded
FileNotFoundError – when module_path doesn’t exist
- radvel.utils.round_sig(x: float, sig: int = 2) float[source]
Round by significant figures :param x: number to be rounded :type x: float :param sig: (optional) number of significant figures to retain :type sig: int
- Returns:
x rounded to sig significant figures
- Return type:
float
- radvel.utils.semi_amplitude(Msini: float, P: float, Mtotal: float, e: float, Msini_units: str = 'jupiter') float[source]
Compute Doppler semi-amplitude
- Parameters:
Msini (float) – mass of planet [Mjup]
P (float) – Orbital period [days]
Mtotal (float) – Mass of star + mass of planet [Msun]
e (float) – eccentricity
Msini_units (Optional[str]) – Units of Msini {‘earth’,’jupiter’} default: ‘jupiter’
- Returns:
Doppler semi-amplitude [m/s]
- radvel.utils.semi_major_axis(P: float, Mtotal: float) float[source]
Semi-major axis
Kepler’s third law
- Parameters:
P (float) – Orbital period [days]
Mtotal (float) – Mass [Msun]
- Returns:
semi-major axis in AU
- Return type:
float or array
- radvel.utils.sigfig(med: float, errlow: float, errhigh: float | None = None) tuple[float, float, float][source]
Format values with errors into an equal number of signficant figures.
- Parameters:
med (float) – median value
errlow (float) – lower errorbar
errhigh (float) – upper errorbar
- Returns:
(med,errlow,errhigh) rounded to the lowest number of significant figures
- Return type:
tuple
- radvel.utils.t_to_phase(params: Parameters, t: ndarray[tuple[Any, ...], dtype[float]], num_planet: int, cat: bool = False) ndarray[tuple[Any, ...], dtype[float]][source]
Time to phase
Convert JD to orbital phase
- Parameters:
params (radvel.params.RVParameters) – RV parameters object
t (array) – JD timestamps
num_planet (int) – Which planet’s ephemeris to phase fold on
cat (bool) – Concatenate/double the output phase array to extend from 0 to 2
- Returns:
orbital phase at each timestamp
- Return type:
array
- radvel.utils.time_print(tdiff: float) tuple[float, str][source]
Print time
Helper function to print time remaining in sensible units.
- Parameters:
tdiff (float) – time in seconds
- Returns:
(float time, string units)
- Return type:
tuple
- radvel.utils.timebin(time: ndarray[tuple[Any, ...], dtype[_ScalarT]], meas: ndarray[tuple[Any, ...], dtype[_ScalarT]], meas_err: ndarray[tuple[Any, ...], dtype[_ScalarT]], binsize: float) tuple[ndarray, ndarray, ndarray][source]
Bin in equal sized time bins
This routine bins a set of times, measurements, and measurement errors into time bins. All inputs and outputs should be floats or double. binsize should have the same units as the time array. (from Andrew Howard, ported to Python by BJ Fulton)
- Parameters:
time (array) – array of times
meas (array) – array of measurements to be comined
meas_err (array) – array of measurement uncertainties
binsize (float) – width of bins in same units as time array
- Returns:
(bin centers, binned measurements, binned uncertainties)
- Return type:
tuple