API documentation for pykima
kima itself is written in C++, to leverage the increased performance. The Python bindings created with nanobind expose most of the functionality of the data (RVData) and model classes (e.g. RVmodel, GPmodel, etc.), as well as the sampler (kima.run).
To analyse the results and create figures, the pykima
sub-package was created.
It can be used in an IPython shell, a Jupyter notebook, or the standard Python
interpreter.
The typical interface with the results of a run is through the
[KimaResults
][pykima.results.KimaResults] class created by calling the
[kima.load_results
][pykima.results.load_results] function
KimaResults
class
A class to hold, analyse, and display the results from kima
Attributes:
Name | Type | Description |
---|---|---|
model |
str
|
The type of kima model |
priors |
dict
|
A dictionary with the priors used in the model |
ESS |
int
|
Effective sample size |
evidence |
float
|
The log-evidence (\(\ln Z\)) of the model |
information |
float
|
The Kullback-Leibler divergence between prior and posterior |
data |
data_holder
|
The data |
posteriors |
posterior_holder
|
The marginal posterior samples |
parameter_priors
property
A list of priors which can be indexed using self.indices
load(filename=None, diagnostic=False, **kwargs)
classmethod
Load a KimaResults object from the current directory, a pickle file, or a zip file.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
filename
|
str
|
If given, load the model from this file. Can be a zip or pickle file. Defaults to None. |
None
|
diagnostic
|
bool
|
Whether to plot the DNest4 diagnotics. Defaults to False. |
False
|
**kwargs
|
Extra keyword arguments passed to |
{}
|
Returns:
Name | Type | Description |
---|---|---|
res |
KimaResults
|
An object holding the results |
save_pickle(filename=None, postfix=None, verbose=True)
Pickle this KimaResults object into a file.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
filename
|
str
|
The name of the file where to save the model. If not given, a unique name will be generated from the properties of the model. |
None
|
postfix
|
str
|
A string to add to the filename, after the timestamp. |
None
|
verbose
|
bool
|
Print a message. Defaults to True. |
True
|
save_zip(filename, verbose=True)
Save this KimaResults object and the text files into a zip.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
filename
|
str
|
The name of the file where to save the model |
required |
verbose
|
bool
|
Print a message. Defaults to True. |
True
|
get_marginals()
Get the marginal posteriors from the posterior_sample matrix. They go into self.T, self.A, self.E, etc
get_medians()
return the median values of all the parameters
log_prior(sample, debug=False)
map_sample(Np=None, mask=None, printit=True, from_posterior=False)
cached
Get the maximum a posteriori (MAP) sample.
Note
This requires recalculation of the prior for all samples, so it can be a bit slow, depending on the number of posterior samples.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
Np
|
int
|
If given, select only samples with that number of planets. |
None
|
printit
|
bool
|
Whether to print the sample. |
True
|
from_posterior
|
bool
|
If True, return the highest likelihood sample from those that represent the posterior. |
False
|
maximum_likelihood_sample(Np=None, printit=True, mask=None, from_posterior=False)
Get the maximum likelihood sample.
By default, this is the highest likelihood sample found by DNest4.
Note
If from_posterior=True
, the returned sample may change, due to
random choices, between different calls to load_results
.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
Np
|
int
|
If given, select only samples with that number of planets. |
None
|
printit
|
bool
|
Whether to print the sample |
True
|
from_posterior
|
bool
|
If True, return the highest likelihood sample from those that represent the posterior. |
False
|
median_sample(Np=None, printit=True)
Get the median posterior sample.
Warning
Unless all posteriors are Gaussian or otherwise well-behaved, the median sample is usually not the appropriate choice for plots, etc.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
Np
|
int
|
If given, select only samples with that number of planets. |
None
|
printit
|
bool
|
Whether to print the sample |
True
|
eval_model(sample, t=None, include_planets=True, include_known_object=True, include_transiting_planet=True, include_indicator_correlations=True, include_trend=True, single_planet=None, except_planet=None)
Evaluate the deterministic part of the model at one posterior sample
.
If t
is None, use the observed times.
Note
Instrument offsets are only added if t
is None, but the systemic
velocity is always added.
Note
This function does not evaluate the GP component of the model.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
sample
|
array
|
One posterior sample, with shape (npar,) |
required |
t
|
array
|
Times at which to evaluate the model, or None to use observed times |
None
|
include_planets
|
bool
|
Whether to include the contribution from the planets |
True
|
include_known_object
|
bool
|
Whether to include the contribution from the known object planet(s) |
True
|
include_transiting_planet
|
bool
|
Whether to include the contribution from the transiting planet(s) |
True
|
include_indicator_correlations
|
bool
|
Whether to include the indicator correlation model |
True
|
include_trend
|
bool
|
Whether to include the contribution from the trend |
True
|
single_planet
|
int
|
Index of a single planet to include in the model, starting at 1. Use positive values (1, 2, ...) for the Np planets and negative values (-1, -2, ...) for the known object and transiting planets. |
None
|
except_planet
|
Union[int, List]
|
Index (or list of indices) of a single planet to exclude from the model, starting at 1. Use positive values (1, 2, ...) for the Np planets and negative values (-1, -2, ...) for the known object and transiting planets. |
None
|
planet_model(sample, t=None, include_known_object=True, include_transiting_planet=True, single_planet=None, except_planet=None)
Evaluate the planet part of the model at one posterior sample
.
If t
is None, use the observed times.
Note
This function does not evaluate the GP component of the model nor the systemic velocity and instrument offsets.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
sample
|
array
|
One posterior sample, with shape (npar,) |
required |
t
|
array
|
Times at which to evaluate the model, or None to use observed times |
None
|
include_known_object
|
bool
|
Whether to include the contribution from the known object planet(s) |
True
|
include_transiting_planet
|
bool
|
Whether to include the contribution from the transiting planet(s) |
True
|
single_planet
|
int
|
Index of a single planet to include in the model, starting at 1. Use positive values (1, 2, ...) for the Np planets and negative values (-1, -2, ...) for the known object and transiting planets. |
None
|
except_planet
|
Union[int, List]
|
Index (or list of indices) of a single planet to exclude from the model, starting at 1. Use positive values (1, 2, ...) for the Np planets and negative values (-1, -2, ...) for the known object and transiting planets. |
None
|
Tip
To evaluate at all posterior samples, consider using
Examples:
To get the Keplerian contribution from the first planet in a
posterior sample p
use:
For, e.g., the second known object in the model, use:
or to get the contributions from all planets except that one
stochastic_model(sample, t=None, return_std=False, derivative=False, include_jitters=True, **kwargs)
Evaluate the stochastic part of the model (GP) at one posterior sample.
If t
is None, use the observed times.
Note
Instrument offsets are only added if t
is None, but the systemic
velocity is always added.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
sample
|
array
|
One posterior sample, with shape (npar,) |
required |
t
|
ndarray
|
Times at which to evaluate the model, or |
None
|
return_std
|
bool
|
Whether to return the standard deviation of the predictive. Default is False. |
False
|
derivative
|
bool
|
Return the first time derivative of the GP prediction instead |
False
|
include_jitters
|
bool
|
Whether to include the jitter values in |
True
|
full_model(sample, t=None, **kwargs)
Evaluate the full model at one posterior sample, including the GP. If
t
is None
, use the observed times. Instrument offsets are only added
if t
is None
, but the systemic velocity is always added.
To evaluate at all posterior samples, consider using
Parameters:
Name | Type | Description | Default |
---|---|---|---|
sample
|
array
|
One posterior sample, with shape (npar,) |
required |
t
|
ndarray
|
Times at which to evaluate the model, or |
None
|
**kwargs
|
Keyword arguments passed directly to |
{}
|
burst_model(sample, t=None, v=None)
For models with multiple instruments, this function "bursts" the
computed RV into n_instruments
individual arrays. This is mostly
useful for plotting the RV model together with the original data.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
sample
|
array
|
One posterior sample, with shape (npar,) |
required |
t
|
array
|
Times at which to evaluate the model |
None
|
v
|
array
|
Pre-computed RVs. If |
None
|
from_prior(n=1)
Generate n
samples from the priors for all parameters.
print_results(show_prior=False)
Print a summary of the results, showing the posterior estimates for each parameter.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
show_prior
|
bool
|
Whether to show the prior distribution. |
False
|
Note
This function is heavily inspired by a similar implementation in UltraNest.
kima.load_results
API
Load the results from a kima run
Parameters:
Name | Type | Description | Default |
---|---|---|---|
model_or_file
|
str or Model
|
If a string, load results from a pickle or zip file. If a model
(e.g. |
required |
data
|
RVData
|
An instance of |
None
|
diagnostic
|
bool
|
Show the diagnostic plots |
False
|
verbose
|
bool
|
Print some information about the results |
True
|
moreSamples
|
int
|
The total number of posterior samples will be equal to ESS * moreSamples |
1
|
Raises:
Type | Description |
---|---|
FileNotFoundError
|
If |
Returns:
Name | Type | Description |
---|---|---|
res |
KimaResults
|
An instance of |
A number of convenience plots are redirected from the
pykima.display
module into any
KimaResults
instance:
Note
Several of the plotting functions are aliased and numbered, as in
This can be convenient, but it's mostly for historical reasons and the API might change at any time. Therefore, these funcions are not documented.
kima.pykima.analysis
API
np_most_probable(results)
Return the value of Np with the highest posterior probability.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
results
|
KimaResults
|
A results instance |
required |
np_bayes_factor_threshold(results, threshold=150)
Return the value of Np supported by the data considering a posterior ratio (Bayes factor) threshold.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
results
|
KimaResults
|
A results instance |
required |
threshold
|
float
|
Posterior ratio threshold. |
150
|
np_posterior_threshold(results, threshold=0.9)
Return the value of Np supported by the data considering an absolute posterior probability threshold.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
results
|
KimaResults
|
A results instance |
required |
threshold
|
float
|
Posterior probability threshold |
0.9
|
get_planet_mass(P, K, e, star_mass=1.0, full_output=False)
Calculate the planet (minimum) mass, \(M_p \sin i\), given orbital period P
,
semi-amplitude K
, eccentricity e
, and stellar mass. If star_mass
is a
tuple with (estimate, uncertainty), this (Gaussian) uncertainty will be
taken into account in the calculation.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
P
|
Union[float, ndarray]
|
orbital period [days] |
required |
K
|
Union[float, ndarray]
|
semi-amplitude [m/s] |
required |
e
|
Union[float, ndarray]
|
orbital eccentricity |
required |
star_mass
|
Union[float, Tuple]
|
stellar mass, or (mass, uncertainty) [Msun] |
1.0
|
This function returns different results depending on the inputs.
If P
, K
, and e
are floats and star_mass
is a float
Returns:
Name | Type | Description |
---|---|---|
Msini |
float
|
planet mass, in \(M_{\rm Jup}\) |
Msini |
float
|
planet mass, in \(M_{\rm Earth}\) |
If P
, K
, and e
are floats and star_mass
is a tuple
Returns:
Name | Type | Description |
---|---|---|
Msini |
tuple
|
planet mass and uncertainty, in \(M_{\rm Jup}\) |
Msini |
tuple
|
planet mass and uncertainty, in \(M_{\rm Earth}\) |
If P
, K
, and e
are arrays and full_output=True
Returns:
Name | Type | Description |
---|---|---|
m_Msini |
float
|
posterior mean for the planet mass, in \(M_{\rm Jup}\) |
s_Msini |
float
|
posterior standard deviation for the planet mass, in \(M_{\rm Jup}\) |
Msini |
array
|
posterior samples for the planet mass, in \(M_{\rm Jup}\) |
If P
, K
, and e
are arrays and full_output=False
Returns:
Name | Type | Description |
---|---|---|
m_Msini |
float
|
posterior mean for the planet mass, in \(M_{\rm Jup}\) |
s_Msini |
float
|
posterior standard deviation for the planet mass, in \(M_{\rm Jup}\) |
m_Msini |
float
|
posterior mean for the planet mass, in \(M_{\rm Earth}\) |
s_Msini |
float
|
posterior standard deviation for the planet mass, in \(M_{\rm Earth}\) |
get_planet_semimajor_axis(P, K, star_mass=1.0, full_output=False)
Calculate the semi-major axis of the planet's orbit given orbital period
P
, semi-amplitude K
, and stellar mass.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
P
|
Union[float, ndarray]
|
orbital period [days] |
required |
K
|
Union[float, ndarray]
|
semi-amplitude [m/s] |
required |
star_mass
|
Union[float, Tuple]
|
stellar mass, or (mass, uncertainty) [Msun] |
1.0
|
This function returns different results depending on the inputs.
If P
and K
are floats and star_mass
is a float
Returns:
Name | Type | Description |
---|---|---|
a |
float
|
planet semi-major axis, in AU |
If P
and K
are floats and star_mass
is a tuple
Returns:
Name | Type | Description |
---|---|---|
a |
tuple
|
semi-major axis and uncertainty, in AU |
If P
and K
are arrays and full_output=True
Returns:
Name | Type | Description |
---|---|---|
m_a |
float
|
posterior mean for the semi-major axis, in AU |
s_a |
float
|
posterior standard deviation for the semi-major axis, in AU |
a |
array
|
posterior samples for the semi-major axis, in AU |
If P
and K
are arrays and full_output=False
Returns:
Name | Type | Description |
---|---|---|
m_a |
float
|
posterior mean for the semi-major axis, in AU |
s_a |
float
|
posterior standard deviation for the semi-major axis, in AU |
get_planet_mass_and_semimajor_axis(P, K, e, star_mass=1.0, full_output=False, verbose=False)
Calculate the planet (minimum) mass Msini and the semi-major axis given
orbital period P
, semi-amplitude K
, eccentricity e
, and stellar mass.
If star_mass is a tuple with (estimate, uncertainty), this (Gaussian)
uncertainty will be taken into account in the calculation.
Units
P [days] K [m/s] e [] star_mass [Msun]
Returns: (M, A) where M is the output of get_planet_mass A is the output of get_planet_semimajor_axis
FIP(results, plot=True, show_peaks=True, oversampling=1)
Calculate (and plot) the True and False Inclusion Probability (TIP/FIP)
Parameters:
Name | Type | Description | Default |
---|---|---|---|
res
|
KimaResults
|
The |
required |
plot
|
bool
|
Plot the TIP and FIP. Defaults to True. |
True
|
show_peaks
|
bool
|
Identify and show prominent TIP peaks. Defaults to True. |
True
|
oversampling
|
int
|
Oversampling factor for the period binning. Defaults to 1. |
1
|
Returns:
Name | Type | Description |
---|---|---|
bins |
ndarray
|
The period bins |
FIP |
ndarray
|
The False Inclusion Probability |
find_outliers(results, sample, threshold=10, full_output=False)
Estimate which observations are outliers, for a model with a Student t
likelihood. This function first calculates the residuals given the
parameters in sample
. Then it computes the relative probability of each
residual point given a Student-t (Td) and a Gaussian (Gd) likelihoods. If
the probability Td is larger than Gd (by a factor of threshold
), the point
is flagged as an outlier. The function returns an "outlier mask". If
full_output
is True, it also returns the ratio Td/Gd.
detection_limits(results, star_mass=1.0, Np=None, bins=200, plot=True, ax=None, semi_amplitude=False, semi_major_axis=False, logX=True, return_mask=False, remove_nan=True, show_eccentricity=False, K_lines=(5, 3, 1), smooth=False, smooth_degree=3)
Calculate detection limits using samples with more than Np
planets. By
default, this function uses the value of Np
which passes the posterior
probability threshold.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
star_mass
|
Union[float, Tuple]
|
Stellar mass and optionally its uncertainty [in solar masses]. |
1.0
|
Np
|
int
|
Consider only posterior samples with more than |
None
|
bins
|
int
|
Number of bins at which to calculate the detection limits. The period ranges from the minimum to the maximum orbital period in the posterior. |
200
|
plot
|
bool
|
Whether to plot the detection limits |
True
|
semi_amplitude
|
bool
|
Show the detection limits for semi-amplitude, instead of planet mass |
False
|
semi_major_axis
|
bool
|
Show semi-major axis in the x axis, instead of orbital period |
False
|
logX
|
bool
|
Should X-axis bins be logarithmic? |
True
|
return_mask
|
bool
|
|
False
|
remove_nan
|
bool
|
remove bins with no counts (no posterior samples) |
True
|
smooth
|
bool
|
Smooth the binned maximum with a polynomial |
False
|
smooth_degree
|
int
|
Degree of the polynomial used for smoothing |
3
|
Returns:
Name | Type | Description |
---|---|---|
P,K,E,M (ndarray):
Orbital periods, semi-amplitudes, eccentricities, and planet masses
used in the calculation of the detection limits. These correspond to
all the posterior samples with more than |
||
s |
namedtuple
|
Detection limits result, with attributes |
hpd_grid(sample, alpha=0.05, bw_method=None, roundto=2)
Calculate highest posterior density (HPD) of array for given alpha. The HPD is the minimum width Bayesian credible interval (BCI). The function works for multimodal distributions, returning more than one mode.
Warning
This function is experimental.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
sample
|
array
|
An array containing MCMC samples |
required |
alpha
|
float
|
Desired probability of type I error (defaults to 0.05) |
0.05
|
bw_method
|
float
|
The method used to estimate the KDE bandwidth. See
|
None
|
roundto
|
integer
|
Number of digits after the decimal point for the results |
2
|
kima.pykima.display
API
plot_posterior_np(res, ax=None, errors=False, show_ESS=True, show_detected=True, show_probabilities=False, show_title=True, **kwargs)
Plot the histogram of the posterior for Np
Parameters:
Name | Type | Description | Default |
---|---|---|---|
res
|
KimaResults
|
The |
required |
ax
|
Axes
|
An existing matplotlib axes where to draw the plot |
None
|
errors
|
bool
|
Whether to estimate and display errors on the Np posterior. |
False
|
show_ESS
|
bool
|
Display the effective sample size on the plot. |
True
|
show_detected
|
bool
|
Highlight the detected Np in the plot (from
|
True
|
show_probabilities
|
bool
|
Display the probabilities on top of the histogram bars. |
False
|
show_title
|
bool
|
Display the title on the plot |
True
|
**kwargs
|
Keyword arguments to pass to the |
{}
|
Returns:
Name | Type | Description |
---|---|---|
fig |
Figure
|
The matplotlib figure with the plot |
plot_posterior_period(res, nbins=100, bins=None, plims=None, logx=True, kde=False, kde_bw=None, show_peaks=False, show_prior=False, show_year=True, show_timespan=True, show_aliases=False, include_known_object=False, include_transiting_planet=False, separate_colors=False, return_bins=False, mark_periods=None, mark_periods_text=False, **kwargs)
Plot the histogram (or the kde) of the posterior for the orbital period(s)
Parameters:
Name | Type | Description | Default |
---|---|---|---|
res
|
KimaResults
|
The |
required |
nbins
|
int
|
Number of bins in the histogram |
100
|
bins
|
array
|
Histogram bins |
None
|
plims
|
tuple
|
Period limits, as (pmin, pmax) |
None
|
logx
|
bool
|
Plot the x axis in lograithmic scale |
True
|
kde
|
bool
|
Show a kernel density estimation (KDE) instead of the histogram |
False
|
kde_bw
|
float
|
Bandwith for the KDE |
None
|
show_peaks
|
bool
|
Try to identify prominent peaks in the posterior |
False
|
show_prior
|
bool
|
Plot the prior together with the posterior |
False
|
show_year
|
bool
|
Show a vertical line at 1 year |
True
|
show_timespan
|
bool
|
Show a vertical line at the timespan of the data |
True
|
show_aliases
|
bool or int
|
Show daily and yearly aliases for top peak(s) |
False
|
separate_colors
|
bool
|
Show different Keplerians as different colors |
False
|
return_bins
|
bool
|
Return the bins used for the histogram |
False
|
mark_periods
|
(list, tuple, array)
|
Mark specific periods in the plot |
None
|
mark_periods_text
|
bool
|
Write the period values next to the marker |
False
|
Returns:
Name | Type | Description |
---|---|---|
fig |
Figure
|
The matplotlib figure with the plot |
plot_PKE(res, mask=None, include_known_object=False, include_transiting_planet=False, show_prior=False, show_aliases=None, reorder_P=False, sort_by_increasing_P=False, points=True, colors_np=True, gridsize=50, **kwargs)
Plot the 2d histograms of the posteriors for semi-amplitude and orbital
period and for eccentricity and orbital period. If points
is True, plot
each posterior sample, else plot hexbins
plot_gp(res, Np=None, ranges=None, show_prior=False, fig=None, axs=None, **hist_kwargs)
Plot histograms for the GP hyperparameters. If Np is not None, highlight the samples with Np Keplerians.
plot_gp_rvfwhm(res, Np=None, ranges=None, show_prior=False, fig=None, **hist_kwargs)
Plot histograms for the GP hyperparameters. If Np is not None, highlight the samples with Np Keplerians.
plot_gp_corner(res, include_jitters=False, ranges=None)
Corner plot for the GP hyperparameters
corner_planet_parameters(res, fig=None, Np=None, true_values=None, period_ranges=None, include_known_object=False, include_transiting_planet=False, KO_Np=None, TR_Np=None, show_prior=False, wrap_M0=False, replace_angles_with_mass=False, star_mass=1.0, a_factor=1.0)
Corner plots of the posterior samples for the planet parameters
hist_vsys(res, show_offsets=True, specific=None, show_prior=False, **kwargs)
Plot the histogram of the posterior for the systemic velocity and for the between-instrument offsets. (if the model has multiple instruments).
Parameters:
Name | Type | Description | Default |
---|---|---|---|
res
|
KimaResults
|
The |
required |
show_offsets
|
bool
|
Whether to plot the histograms for the between-instrument offsets |
True
|
specific
|
tuple
|
If not None, it should be a tuple with two instrument names
(matching |
None
|
show_prior
|
bool
|
Whether to plot the histogram of the prior distribution |
False
|
hist_jitter(res, show_prior=False, show_stats=False, show_title=True, ax=None, **kwargs)
Plot the histogram of the posterior for the additional white noise
hist_correlations(res, show_prior=False)
Plot the histogram of the posterior for the activity correlations
hist_trend(res, per_year=True, ax=None, show_prior=False, show_title=True)
Plot the histogram of the posterior for the coefficients of the trend
hist_MA(res)
Plot the histogram of the posterior for the MA parameters
hist_nu(res, show_prior=False, **kwargs)
Plot the histogram of the posterior for the Student-t degrees of freedom
plot_RVData(data, **kwargs)
Simple plot of RV data
phase_plot(res, sample, phase_axs=None, sort_by_increasing_P=False, sort_by_decreasing_K=True, highlight=None, highlight_points=None, only=None, add_titles=True, sharey=False, show_gls_residuals=False, **kwargs)
Plot the planet phase curves, GP, and residuals, for a given sample
.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
res
|
KimaResults
|
The |
required |
sample
|
array
|
Array with one posterior sample |
required |
phase_axs
|
list[Axes]
|
One or more axes for the phase plot(s) |
None
|
sort_by_increasing_P
|
bool
|
Sort the planets by increasing period |
False
|
sort_by_decreasing_K
|
bool
|
Sort the planets by decreasing semi-amplitude |
True
|
highlight
|
list
|
Highlight all data points from a specific instrument |
None
|
highlight_points
|
list
|
Highlight specific data points by index |
None
|
only
|
list
|
Only show data from specific instrument(s) |
None
|
add_titles
|
bool
|
Add titles to each phase plot |
True
|
sharey
|
bool
|
Share the y-axis of the phase plots |
False
|
show_gls_residuals
|
bool
|
Add a panel with the Lomb-Scargle periodogram of the residuals |
False
|
**kwargs
|
dict
|
Keyword arguments passed to |
{}
|
Warning
This is probably the most complicated function in the whole package! For one, the layout of the axes in the figure may not always be optimal.
plot_random_samples_multiseries(res, ncurves=50, samples=None, over=0.1, ntt=10000, show_vsys=False, isolate_known_object=True, include_jitters=True, just_rvs=False, full_plot=False, highest_likelihood=False, **kwargs)
Display the RV data together with curves from the posterior predictive.
A total of ncurves
random samples are chosen, and the Keplerian
curves are calculated covering 100 + over
% of the data timespan.
If the model has a GP component, the prediction is calculated using the
GP hyperparameters for each of the random samples.
kima.pykima.utils
API
Small (but sometimes important) utility functions
Fixed(value)
Not a real distribution, just a parameter fixed at a given value
logpdf(x)
Logarithm of the probability density function
pdf(x)
Probability density function: 1.0 if x==value, otherwise 0.0
rvs()
Random sample (always returns value
)
ZeroDist()
A dummy probability distribution which always returns 0.0
apply_argsort(arr1, arr2, axis=-1)
Apply arr1
.argsort() on arr2
, along axis
.
chdir(dir, create=True)
A simple context manager to switch directories temporarily
clipped_mean(arr, min, max)
Mean of arr
between min
and max
clipped_std(arr, min, max)
Standard deviation of arr
between min
and max
date_date_to_jd(year, month, day)
Convert a date (year, month, day) to Julian Day.
Algorithm from 'Practical Astronomy with your Calculator or Spreadsheet', 4th ed., Duffet-Smith and Zwart, 2011.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
year
|
int
|
Year as integer. Years preceding 1 A.D. should be 0 or negative. The year before 1 A.D. is 0, 10 B.C. is year -9. |
required |
month
|
int
|
Month as integer, Jan = 1, Feb. = 2, etc. |
required |
day
|
float
|
Day, may contain fractional part. |
required |
Returns:
Name | Type | Description |
---|---|---|
jd |
float
|
Julian Day |
Examples:
Convert 6 a.m., February 17, 1985 to Julian Day
date_datetime_to_jd(date)
Convert a datetime.datetime
object to Julian day.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
date
|
datetime
|
the date to be converted to Julian day |
required |
Returns:
Name | Type | Description |
---|---|---|
jd |
float
|
Julian day |
Examples:
date_days_to_hmsm(days)
Convert fractional days to hours, minutes, seconds, and microseconds. Precision beyond microseconds is rounded to the nearest microsecond.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
days
|
float
|
A fractional number of days. Must be less than 1. |
required |
Returns:
Name | Type | Description |
---|---|---|
hour |
int
|
Hour number |
min |
int
|
Minute number |
sec |
int
|
Second number |
micro |
int
|
Microsecond number |
Raises:
Type | Description |
---|---|
ValueError
|
If |
Examples:
date_hmsm_to_days(hour=0, min=0, sec=0, micro=0)
Convert hours, minutes, seconds, and microseconds to fractional days.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
hour
|
int
|
Hour number. Defaults to 0. |
0
|
min
|
int
|
Minute number. Defaults to 0. |
0
|
sec
|
int
|
Second number. Defaults to 0. |
0
|
micro
|
int
|
Microsecond number. Defaults to 0. |
0
|
Returns:
Name | Type | Description |
---|---|---|
days |
float
|
Fractional days |
Examples:
date_jd_to_date(jd)
Convert Julian Day to date.
Algorithm from 'Practical Astronomy with your Calculator or Spreadsheet', 4th ed., Duffet-Smith and Zwart, 2011.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
jd
|
float
|
Julian Day |
required |
Returns:
Name | Type | Description |
---|---|---|
year |
int
|
Year as integer. Years preceding 1 A.D. should be 0 or negative. The year before 1 A.D. is 0, 10 B.C. is year -9. |
month |
int
|
Month as integer, Jan = 1, Feb. = 2, etc. |
day |
float
|
Day, may contain fractional part. |
Examples:
Convert Julian Day 2446113.75 to year, month, and day.
date_jd_to_mjd(jd)
Convert Julian Day to Modified Julian Day
Parameters:
Name | Type | Description | Default |
---|---|---|---|
jd
|
float
|
Julian Day |
required |
Returns mjd (float): Modified Julian Day
date_mjd_to_jd(mjd)
Convert Modified Julian Day to Julian Day.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
mjd
|
float
|
Modified Julian Day |
required |
Returns jd (float): Julian Day
distribution_rvs(dist, size=1)
Generate random samples from a distribution.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
dist
|
Distribution or rv_continuous
|
The distribution to sample from |
required |
size
|
int
|
Number of samples to generate |
1
|
distribution_support(dist)
Return the support of a distribution
find_constrained_invgamma(low, upp, low_tail=0.01, upp_tail=0.01)
Find the parameters of the inverse gamma distribution that is constrained
between low
and upp
with tail probabilities low_tail
and upp_tail
.
find_prior_limits(prior)
Find lower and upper limits of a prior from the kima_model_setup.txt file.
find_prior_parameters(prior)
Find parameters of a prior from the kima_model_setup.txt file.
get_gaussian_prior_vsys(data, use_ptp=True, use_std=False, use_value=None)
Get an informative Gaussian prior for the systemic velocity using the data.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
data
|
RVData
|
The RV dataset |
required |
use_ptp
|
bool
|
Width of the prior is the range of the data. |
True
|
use_std
|
bool
|
Width of the prior is the standard deviation of the data. |
False
|
use_value
|
float
|
Width of the prior is provided directly. |
None
|
Returns:
Name | Type | Description |
---|---|---|
prior |
Gaussian
|
A Gaussian prior for the systemic velocity |
Raises:
Type | Description |
---|---|
ValueError
|
If both |
get_gaussian_priors_individual_offsets(data, use_ptp=True, use_std=False, use_value=None)
Get a informative Gaussian priors for the between-instrument offsets using the data.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
data
|
RVData
|
The RV dataset (must have multiple instruments) |
required |
use_ptp
|
bool
|
Width of the prior is the range of the data. |
True
|
use_std
|
bool
|
Width of the prior is the standard deviation of the data. |
False
|
use_value
|
float
|
Width of the prior is provided directly. |
None
|
Returns:
Name | Type | Description |
---|---|---|
prior |
Gaussian
|
Gaussian priors appropriate for the between-instrument offsets |
Raises:
Type | Description |
---|---|
ValueError
|
|
get_instrument_name(data_file)
Try to find the instrument name from a filename
get_mean_longitude(M0, w, fold=False, deg=False)
Return the mean longitude λ for a given mean anomaly M0 and longitude of periastron w (given by λ = M0 + w).
Parameters:
Name | Type | Description | Default |
---|---|---|---|
M0
|
float or ndarray
|
Mean anomaly (M0=0 at periastron) [rad] |
required |
w
|
float or ndarray
|
Longitude of periastron [rad] |
required |
fold
|
bool
|
Fold the mean longitude into the range [-π, π) [default: False] |
False
|
deg
|
bool
|
Return mean longitude in degrees [default: False] |
False
|
Returns:
Name | Type | Description |
---|---|---|
λ |
float or ndarray
|
Mean longitude [rad or deg] |
get_planet_teq(Tstar=5777, Rstar=1, a=1, A=0, f=1)
Calculate the planet's equlibrium temperature
Parameters:
Name | Type | Description | Default |
---|---|---|---|
Tstar
|
float
|
Stellar effective temperature [K]. Defaults to 5777. |
5777
|
Rstar
|
float
|
Stellar radius [Rsun]. Defaults to 1.0. |
1
|
a
|
float
|
Planet semi-major axis [AU]. Defaults to 1.0. |
1
|
A
|
float
|
Bond albedo. Defaults to 0. |
0
|
f
|
float
|
Redistribution factor. Defaults to 1. |
1
|
Returns:
Name | Type | Description |
---|---|---|
Teq |
float
|
Planet equiolibrium temperature [K] |
get_prior(prior)
Return a scipt.stats-like prior from a kima_model_setup.txt prior
get_star_name(data_file)
Find star name (usually works for approx. standard filenames)
get_transit_probability(Rstar=1.0, a=1.0)
Transit probability, simple. Eq 6.1 in Exoplanet Handbook Rstar: stellar radius [Rsun] a: semi-major axis [au]
hide_stdout()
lighten_color(color, amount=0.5)
mean_anomaly_from_epoch(t, P, M0, epoch)
Calculate mean anomaly for times t, given period P, mean anomaly at the epoch M0, and epoch.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t
|
(float, ndarray)
|
Times at which to calculate mean anomaly [days] |
required |
P
|
float
|
Orbital period [days] |
required |
M0
|
float
|
Mean anomaly at the epoch [rad] |
required |
epoch
|
float
|
Epoch [days] |
required |
Returns:
Name | Type | Description |
---|---|---|
M |
(float, ndarray)
|
Mean anomaly [rad] |
mean_anomaly_from_time_periastron(t, P, Tp, deg=False)
Calculate mean anomaly for times t, given period P and time of periastron Tp.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
t
|
(float, ndarray)
|
Times at which to calculate mean anomaly [days] |
required |
P
|
float
|
Orbital period [days] |
required |
Tp
|
float
|
Time of periastron [days] |
required |
deg
|
bool
|
Return mean anomaly in degrees [default: False] |
False
|
Returns:
Name | Type | Description |
---|---|---|
M |
(float, ndarray)
|
Mean anomaly [rad or deg] |
percentile68_ranges(a, min=None, max=None)
Calculate the 16th and 84th percentiles of values in a
, clipped between
min
and max
.
minus median plus
-------==========|==========------
16% | 68% | 16%
Returns:
Name | Type | Description |
---|---|---|
median |
float
|
Median value of |
plus |
float
|
The 84th percentile minus the median. |
minus |
float
|
The median minus the 16th percentile. |
percentile68_ranges_latex(a, min=None, max=None, collapse=True, include_dollar=True)
Return a LaTeX-formatted string of the 68% range of values in a
, clipped
between a
and b
, in the form
$ median ^{+ plus} _{- minus} $
Parameters:
Name | Type | Description | Default |
---|---|---|---|
collapse
|
bool
|
If True and plus=minus, return $ median \pm plus $ |
True
|
include_dollar
|
bool
|
Whether to include dollar signs in the output, so it's a valid LaTeX math mode string |
True
|
percentile_ranges(a, percentile=68, min=None, max=None)
Calculate the given percentile range of values in a
, clipped between min
and max
.
Returns:
Name | Type | Description |
---|---|---|
median |
float
|
Median value of |
plus |
float
|
The 84th percentile minus the median. |
minus |
float
|
The median minus the 16th percentile. |
read_datafile(datafile, skip)
Read data from datafile
for multiple instruments.
Can be str, in which case the last column is assumed to contain an integer
identifier of the instrument.
Or list, in which case each element will be one different filename
containing three columns each.
read_datafile_rvfwhm(datafile, skip)
Read data from datafile
for multiple instruments and RV-FWHM data.
read_datafile_rvfwhmrhk(datafile, skip)
Read data from datafile
for multiple instruments and RV-FWHM-RHK data.
rms(array)
Root mean square of array
show_tips()
Show a few tips on how to use kima (but only sometimes)
wrms(array, weights)
Weighted root mean square of array, given weights