Plotting Methods
- orbitize.plot.plot_corner(results, param_list=None, **corner_kwargs)[source]
Make a corner plot of posterior on orbit fit from any sampler
- Parameters
param_list (list of strings) –
each entry is a name of a parameter to include. Valid strings:
sma1: semimajor axis ecc1: eccentricity inc1: inclination aop1: argument of periastron pan1: position angle of nodes tau1: epoch of periastron passage, expressed as fraction of orbital period per1: period K1: stellar radial velocity semi-amplitude [repeat for 2, 3, 4, etc if multiple objects] plx: parallax pm_ra: RA proper motion pm_dec: Dec proper motion alpha0: primary offset from reported Hipparcos RA @ alphadec0_epoch (generally 1991.25) delta0: primary offset from reported Hipparcos Dec @ alphadec0_epoch (generally 1991.25) gamma: rv offset sigma: rv jitter mi: mass of individual body i, for i = 0, 1, 2, ... (only if fit_secondary_mass == True) mtot: total mass (only if fit_secondary_mass == False)
**corner_kwargs – any remaining keyword args are sent to
corner.corner
. See here. Note: default axis labels used unless overwritten by user input.
- Returns
corner plot
- Return type
matplotlib.pyplot.Figure
Note
Example: Use
param_list = ['sma1,ecc1,inc1,sma2,ecc2,inc2']
to only plot posteriors for semimajor axis, eccentricity and inclination of the first two companionsWritten: Henry Ngo, 2018
- orbitize.plot.plot_orbits(results, 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, plot_errorbars=True, rv_time_series2=False, primary_instrument_name=None, fontsize=20, fig=None)[source]
Plots one orbital period for a select number of fitted orbits for a given object, with line segments colored according to time
- Parameters
object_to_plot (int) – which object to plot (default: 1)
start_mjd (float) – MJD in which to start plotting orbits (default: 51544, the year 2000)
num_orbits_to_plot (int) – number of orbits to plot (default: 100)
num_epochs_to_plot (int) – number of points to plot per orbit (default: 100)
square_plot (Boolean) – Aspect ratio is always equal, but if square_plot is True (default), then the axes will be square, otherwise, white space padding is used
show_colorbar (Boolean) – Displays colorbar to the right of the plot [True]
cmap (matplotlib.cm.ColorMap) – color map to use for making orbit tracks (default: modified Purples_r)
sep_pa_color (string) – any valid matplotlib color string, used to set the color of the orbit tracks in the Sep/PA panels (default: ‘lightgrey’).
sep_pa_end_year (float) – decimal year specifying when to stop plotting orbit tracks in the Sep/PA panels (default: 2025.0).
cbar_param (string) – options are the following: ‘Epoch [year]’, ‘sma1’, ‘ecc1’, ‘inc1’, ‘aop1’, ‘pan1’, ‘tau1’, ‘plx. Number can be switched out. Default is Epoch [year].
mod180 (Bool) – if True, PA will be plotted in range [180, 540]. Useful for plotting short arcs with PAs that cross 360 deg during observations (default: False)
rv_time_series (Boolean) – if fitting for secondary mass using MCMC for rv fitting and want to display time series, set to True.
plot_astrometry (Boolean) – set to True by default. Plots the astrometric data.
plot_astrometry_insts (Boolean) – set to False by default. Plots the astrometric data by instruments.
fig (matplotlib.pyplot.Figure) – optionally include a predefined Figure object to plot the orbit on. Most users will not need this keyword.
- Returns
the orbit plot if input is valid,
None
otherwise- Return type
matplotlib.pyplot.Figure
(written): Henry Ngo, Sarah Blunt, 2018 Additions by Malena Rice, 2019 Additions by Dino Hsu, 2023
- orbitize.plot.plot_propermotion(results, system, 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, tight_layout=False)[source]
Plots the proper motion of a host star as induced by a companion for one orbital period for a select number of fitted orbits for a given object, with line segments colored according to a given parameter (most informative is usually mass of companion)
Important Note: These plotted trajectories aren’t what are fitting in the likelihood evaluation for the HGCA runs. The implementation forward models the Hip/Gaia measurements per epoch and infers the differential proper motions. This plot is given only for the purposes of an approximate visualization.
- Parameters
system (object) – orbitize.system object with a HGCALogProb passed to system.gaia
object_to_plot (int) – which object to plot (default: 1)
start_mjd (float) – MJD in which to start plotting orbits (default: 51544, the year 2000)
periods_to_plot (int) – number of periods to plot (default: 1)
end_year (float) – decimal year specifying when to stop plotting orbit tracks in the Sep/PA panels (default: 2025.0).
alpha (float) – transparency of lines (default: 0.05)
num_orbits_to_plot (int) – number of orbits to plot (default: 100)
num_epochs_to_plot (int) – number of points to plot per orbit (default: 100)
show_colorbar (Boolean) – Displays colorbar to the right of the plot [True]
cmap (matplotlib.cm.ColorMap) – color map to use for making orbit tracks (default: modified Purples_r)
cbar_param (string) – options are the following: ‘sma1’, ‘ecc1’, ‘inc1’, ‘aop1’, ‘pan1’, ‘tau1’, ‘plx’, ‘m0’, ‘m1’, etc. Number can be switched out. Default is None.
tight_layout (bool) – apply plt.tight_layout function?
fig (matplotlib.pyplot.Figure) – optionally include a predefined Figure object to plot the orbit on. Most users will not need this keyword.
- Returns
the orbit plot if input is valid,
None
otherwise- Return type
matplotlib.pyplot.Figure
(written): William Balmer (2023), based on plot_orbits by Sarah Blunt and Henry Ngo
- orbitize.plot.plot_residuals(my_results, object_to_plot=1, start_mjd=51544, num_orbits_to_plot=100, num_epochs_to_plot=100, sep_pa_color='lightgrey', sep_pa_end_year=2025.0, cbar_param='Epoch [year]', mod180=False)[source]
Plots sep/PA residuals for a set of orbits
- Parameters
my_results (orbitiez.results.Results) – results to plot
object_to_plot (int) – which object to plot (default: 1)
start_mjd (float) – MJD in which to start plotting orbits (default: 51544, the year 2000)
num_orbits_to_plot (int) – number of orbits to plot (default: 100)
num_epochs_to_plot (int) – number of points to plot per orbit (default: 100)
sep_pa_color (string) – any valid matplotlib color string, used to set the color of the orbit tracks in the Sep/PA panels (default: ‘lightgrey’).
sep_pa_end_year (float) – decimal year specifying when to stop plotting orbit tracks in the Sep/PA panels (default: 2025.0).
cbar_param (string) – options are the following: ‘Epoch [year]’, ‘sma1’, ‘ecc1’, ‘inc1’, ‘aop1’, ‘pan1’, ‘tau1’, ‘plx. Number can be switched out. Default is Epoch [year].
mod180 (Bool) – if True, PA will be plotted in range [180, 540]. Useful for plotting short arcs with PAs that cross 360 deg during observations (default: False)
- Returns
the residual plots
- Return type
matplotlib.pyplot.Figure