Results
- class orbitize.results.Results(system=None, sampler_name=None, post=None, lnlike=None, version_number=None, curr_pos=None)[source]
A class to store accepted orbital configurations from the sampler
- Parameters
system (orbitize.system.System) – System object used to do the fit.
sampler_name (string) – name of sampler class that generated these results (default: None).
post (np.array of float) – MxN array of orbital parameters (posterior output from orbit-fitting process), where M is the number of orbits generated, and N is the number of varying orbital parameters in the fit (default: None).
lnlike (np.array of float) – M array of log-likelihoods corresponding to the orbits described in
post
(default: None).version_number (str) – version of orbitize that produced these results.
data (astropy.table.Table) – output from
orbitize.read_input.read_file()
curr_pos (np.array of float) – for MCMC only. A multi-D array of the current walker positions that is used for restarting a MCMC sampler.
Written: Henry Ngo, Sarah Blunt, 2018
API Update: Sarah Blunt, 2021
- add_samples(orbital_params, lnlikes, curr_pos=None)[source]
Add accepted orbits, their likelihoods, and the orbitize version number to the results
- Parameters
orbital_params (np.array) – add sets of orbital params (could be multiple) to results
lnlike (np.array) – add corresponding lnlike values to results
curr_pos (np.array of float) – for MCMC only. A multi-D array of the current walker positions
Written: Henry Ngo, 2018
API Update: Sarah Blunt, 2021
- load_results(filename, append=False)[source]
Populate the
results.Results
object with data from a datafile- Parameters
filename (string) – filepath where data is saved
append (boolean) – if True, then new data is added to existing object. If False (default), new data overwrites existing object
See the
save_results()
method in this module for information on how the data is structured.Written: Henry Ngo, 2018
API Update: Sarah Blunt, 2021
- plot_orbits(object_to_plot=1, start_mjd=51544.0, num_orbits_to_plot=100, num_epochs_to_plot=100, square_plot=True, show_colorbar=True, cmap=<matplotlib.colors.LinearSegmentedColormap object>, sep_pa_color='lightgrey', sep_pa_end_year=2025.0, cbar_param='Epoch [year]', mod180=False, rv_time_series=False, plot_astrometry=True, plot_astrometry_insts=False, fig=None)[source]
Wrapper for orbitize.plot.plot_orbits
- plot_propermotion(object_to_plot=1, start_mjd=44239.0, periods_to_plot=1, end_year=2030.0, alpha=0.05, num_orbits_to_plot=100, num_epochs_to_plot=100, show_colorbar=True, cmap=<matplotlib.colors.LinearSegmentedColormap object>, cbar_param=None)[source]
Wrapper for orbitize.plot.plot_propermotion
- save_results(filename)[source]
Save results.Results object to an hdf5 file
- Parameters
filename (string) – filepath to save to
Save attributes from the
results.Results
object.sampler_name
,tau_ref_epcoh
,version_number
are attributes of the root group.post
,lnlike
, andparameter_labels
are datasets that are members of the root group.Written: Henry Ngo, 2018
API Update: Sarah Blunt, 2021