System

class orbitize.system.System(num_secondary_bodies, data_table, stellar_or_system_mass, plx, mass_err=0, plx_err=0, restrict_angle_ranges=False, tau_ref_epoch=58849, fit_secondary_mass=False, hipparcos_IAD=None, gaia=None, fitting_basis='Standard', use_rebound=False)[source]

A class to store information about a system (data & priors) and calculate model predictions given a set of orbital parameters.

Parameters
  • num_secondary_bodies (int) – number of secondary bodies in the system. Should be at least 1.

  • data_table (astropy.table.Table) – output from orbitize.read_input.read_file()

  • stellar_or_system_mass (float) – mass of the primary star (if fitting for dynamical masses of both components, for example when you have both astrometry and RVs) or total system mass (if fitting for total system mass only, as in the case of a vanilla 2-body fit using relative astrometry only ) [M_sol]

  • plx (float) – mean parallax of the system, in mas

  • mass_err (float, optional) – uncertainty on stellar_or_system_mass, in M_sol

  • plx_err (float, optional) – uncertainty on plx, in mas

  • restrict_angle_ranges (bool, optional) – if True, restrict the ranges of the position angle of nodes to [0,180) to get rid of symmetric double-peaks for imaging-only datasets.

  • tau_ref_epoch (float, optional) – reference epoch for defining tau (MJD). Default is 58849 (Jan 1, 2020).

  • fit_secondary_mass (bool, optional) – if True, include the dynamical mass of the orbiting body as a fitted parameter. If this is set to False, stellar_or_system_mass is taken to be the total mass of the system. (default: False)

  • hipparcos_IAD (orbitize.hipparcos.HipparcosLogProb) – an object containing information & precomputed values relevant to Hipparcos IAD fitting. See hipparcos.py for more details.

  • gaia (orbitize.gaia.GaiaLogProb) – an object containing information & precomputed values relevant to Gaia astrometrry fitting. See gaia.py for more details.

  • fitting_basis (str) – the name of the class corresponding to the fitting basis to be used. See basis.py for a list of implemented fitting bases.

  • use_rebound (bool) – if True, use an n-body backend solver instead of a Keplerian solver.

Priors are initialized as a list of orbitize.priors.Prior objects and stored in the variable System.sys_priors. You should initialize this class, then overwrite priors you wish to customize. You can use the System.param_idx attribute to figure out which indices correspond to which fitting parameters. See the “changing priors” tutorial for more detail.

Written: Sarah Blunt, Henry Ngo, Jason Wang, 2018

compute_all_orbits(params_arr, epochs=None, comp_rebound=False)[source]

Calls orbitize.kepler.calc_orbit and optionally accounts for multi-body interactions. Also computes total quantities like RV (without jitter/gamma)

Parameters
  • params_arr (np.array of float) – RxM array of fitting parameters, where R is the number of parameters being fit, and M is the number of orbits we need model predictions for. Must be in the same order documented in System() above. If M=1, this can be a 1d array.

  • epochs (np.array of float) – epochs (in mjd) at which to compute orbit predictions.

  • comp_rebound (bool, optional) – A secondary optional input for use of N-body solver Rebound; by default, this will be set to false and a Kepler solver will be used instead.

Returns

raoff (np.array of float): N_epochs x N_bodies x N_orbits array of

RA offsets from barycenter at each epoch.

decoff (np.array of float): N_epochs x N_bodies x N_orbits array of

Dec offsets from barycenter at each epoch.

vz (np.array of float): N_epochs x N_bodies x N_orbits array of

radial velocities at each epoch.

Return type

tuple

compute_model(params_arr, use_rebound=False)[source]

Compute model predictions for an array of fitting parameters. Calls the above compute_all_orbits() function, adds jitter/gamma to RV measurements, and propagates these predictions to a model array that can be subtracted from a data array to compute chi2.

Parameters
  • params_arr (np.array of float) – RxM array of fitting parameters, where R is the number of parameters being fit, and M is the number of orbits we need model predictions for. Must be in the same order documented in System() above. If M=1, this can be a 1d array.

  • use_rebound (bool, optional) – A secondary optional input for use of N-body solver Rebound; by default, this will be set to false and a Kepler solver will be used instead.

Returns

np.array of float: Nobsx2xM array model predictions. If M=1, this is

a 2d array, otherwise it is a 3d array.

np.array of float: Nobsx2xM array jitter predictions. If M=1, this is

a 2d array, otherwise it is a 3d array.

Return type

tuple of

convert_data_table_radec2seppa(body_num=1)[source]

Converts rows of self.data_table given in radec to seppa. Note that self.input_table remains unchanged.

Parameters

body_num (int) – which object to convert (1 = first planet)

save(hf)[source]

Saves the current object to an hdf5 file

Parameters

hf (h5py._hl.files.File) – a currently open hdf5 file in which to save the object.

orbitize.system.generate_synthetic_data(orbit_frac, mtot, plx, ecc=0.5, inc=0.7853981633974483, argp=0.7853981633974483, lan=0.7853981633974483, tau=0.8, num_obs=4, unc=2)[source]

Generate an orbitize-table of synethic data

Parameters
  • orbit_frac (float) – percentage of orbit covered by synthetic data

  • mtot (float) – total mass of the system [M_sol]

  • plx (float) – parallax of system [mas]

  • num_obs (int) – number of observations to generate

  • unc (float) – uncertainty on all simulated RA & Dec measurements [mas]

Returns

  • astropy.table.Table: data table of generated synthetic data

  • float: the semimajor axis of the generated data

Return type

2-tuple

orbitize.system.radec2seppa(ra, dec, mod180=False)[source]

Convenience function for converting from right ascension/declination to separation/ position angle.

Parameters
  • ra (np.array of float) – array of RA values, in mas

  • dec (np.array of float) – array of Dec values, in mas

  • mod180 (Bool) – if True, output PA values will be given in range [180, 540) (useful for plotting short arcs with PAs that cross 360 during observations) (default: False)

Returns

(separation [mas], position angle [deg])

Return type

tuple of float

orbitize.system.seppa2radec(sep, pa)[source]

Convenience function to convert sep/pa to ra/dec

Parameters
  • sep (np.array of float) – array of separation in mas

  • pa (np.array of float) – array of position angles in degrees

Returns

(ra [mas], dec [mas])

Return type

tuple

orbitize.system.transform_errors(x1, x2, x1_err, x2_err, x12_corr, transform_func, nsamps=100000)[source]

Transform errors and covariances from one basis to another using a Monte Carlo apporach

Parameters
  • x1 (float) – planet location in first coordinate (e.g., RA, sep) before transformation

  • x2 (float) – planet location in the second coordinate (e.g., Dec, PA) before transformation)

  • x1_err (float) – error in x1

  • x2_err (float) – error in x2

  • x12_corr (float) – correlation between x1 and x2

  • transform_func (function) – function that transforms between (x1, x2) and (x1p, x2p) (the transformed coordinates). The function signature should look like: x1p, x2p = transform_func(x1, x2)

  • nsamps – number of samples to draw more the Monte Carlo approach. More is slower but more accurate.