Hipparcos API Module
- class orbitize.hipparcos.HipparcosLogProb(path_to_iad_file, hip_num, num_secondary_bodies, alphadec0_epoch=1991.25, renormalize_errors=False)[source]
Class to compute the log probability of an orbit with respect to the Hipparcos Intermediate Astrometric Data (IAD). If using a DVD file, queries Vizier for all metadata relevant to the IAD, and reads in the IAD themselves from a specified location. Follows Nielsen+ 2020 (studying the orbit of beta Pic b).
Fitting the Hipparcos IAD requires fitting for the following five parameters. They are added to the vector of fitting parameters in system.py, but are described here for completeness. See Nielsen+ 2020 for more detail.
- alpha0: RA offset from the reported Hipparcos position at a particular
epoch (usually 1991.25) [mas]
- delta0: Dec offset from the reported Hipparcos position at a particular
epoch (usually 1991.25) [mas]
pm_ra: RA proper motion [mas/yr]
pm_dec: Dec proper motion [mas/yr]
plx: parallax [mas]
Note
In orbitize, it is possible to perform a fit to just the Hipparcos IAD, but not to just the Gaia astrometric data.
- Parameters
path_to_iad_file (str) – location of IAD file to be used in your fit. See the Hipparcos tutorial for a walkthrough of how to download these files.
hip_num (str) – Hipparcos ID of star. Available on Simbad. Should have zeros in the prefix if number is <100,000. (i.e. 27321 should be passed in as ‘027321’).
num_secondary_bodies (int) – number of companions in the system
alphadec0_epoch (float) – epoch (in decimal year) that the fitting parameters alpha0 and delta0 are defined relative to (see above).
renormalize_errors (bool) – if True, normalize the scan errors to get chisq_red = 1, following Nielsen+ 2020 (eq 10). In general, this should be False, but it’s helpful for testing. Check out orbitize.hipparcos.nielsen_iad_refitting_test() for an example using this renormalization.
Written: Sarah Blunt & Rob de Rosa, 2021
- compute_lnlike(raoff_model, deoff_model, samples, param_idx)[source]
Computes the log likelihood of an orbit model with respect to the Hipparcos IAD. This is added to the likelihoods calculated with respect to other data types in
sampler._logl()
.- Parameters
raoff_model (np.array of float) – M-dimensional array of primary RA offsets from the barycenter incurred from orbital motion of companions (i.e. not from parallactic motion), where M is the number of epochs of IAD scan data.
deoff_model (np.array of float) – M-dimensional array of primary RA offsets from the barycenter incurred from orbital motion of companions (i.e. not from parallactic motion), where M is the number of epochs of IAD scan data.
samples (np.array of float) – R-dimensional array of fitting parameters, where R is the number of parameters being fit. Must be in the same order documented in
System
.param_idx – a dictionary matching fitting parameter labels to their indices in an array of fitting parameters (generally set to System.basis.param_idx).
- Returns
- array of length M, where M is the number of input
orbits, representing the log likelihood of each orbit with respect to the Hipparcos IAD.
- Return type
np.array of float
- class orbitize.hipparcos.PMPlx_Motion(epochs_mjd, alpha0, delta0, alphadec0_epoch=1991.25)[source]
Class to compute the predicted position at an array of epochs given a parallax and proper motion model (NO orbital motion is added in this class).
- Parameters
times_mjd (np.array of float) – times (in mjd) at which we have absolute astrometric measurements
alpha0 (float) – measured RA position (in degrees) of the object at alphadec0_epoch (see below).
delta0 (float) – measured Dec position (in degrees) of the object at alphadec0_epoch (see below).
alphadec0_epoch (float) – a (fixed) reference time. For stars with Hipparcos data, this should generally be 1991.25, but you can define it however you want. Absolute astrometric data (passed in via an orbitize! data table) should be defined as offsets from the reported position of the object at this epoch (with propagated uncertainties). For example, if you have two absolute astrometric measurements taken with GRAVITY, as well as a Hipparcos-derived position (at epoch 1991.25), alphadec0_epoch should be 1991.25, and you should pass in absolute astrometry in terms of mas offset from the Hipparcos catalog position, with propagated errors of your measurement and the Hipparcos measurement.
- compute_astrometric_model(samples, param_idx, epochs=None)[source]
Compute the astrometric prediction at self.epochs_mjd from parallax and proper motion alone, given an array of model parameters (no orbital motion is added here).
- Parameters
samples (np.array of float) – Length R array of fitting parameters, where R is the number of parameters being fit. Must be in the same order documented in
System
.param_idx – a dictionary matching fitting parameter labels to their indices in an array of fitting parameters (generally set to System.basis.param_idx).
epochs – if None, use self.epochs for astrometric predictions. Otherwise, use this array passed in [in decimalyear].
- Returns
- float: predicted RA*cos(delta0) position offsets from the measured
position at alphadec0_epoch, calculated for each input epoch [mas]
- float: predicted Dec position offsets from the measured position
at alphadec0_epoch, calculated for each input epoch [mas]
- Return type
tuple of
- orbitize.hipparcos.nielsen_iad_refitting_test(iad_file, hip_num='027321', saveplot='bPic_IADrefit.png', burn_steps=100, mcmc_steps=5000)[source]
Reproduce the IAD refitting test from Nielsen+ 2020 (end of Section 3.1). The default MCMC parameters are what you’d want to run before using the IAD for a new system. This fit uses 100 walkers.
- Parameters
iad_loc (str) – path to the IAD file.
hip_num (str) – Hipparcos ID of star. Available on Simbad. Should have zeros in the prefix if number is <100,000. (i.e. 27321 should be passed in as ‘027321’).
saveplot (str) – what to save the summary plot as. If None, don’t make a plot
burn_steps (int) – number of MCMC burn-in steps to run.
mcmc_steps (int) – number of MCMC production steps to run.
- Returns
numpy.array of float: n_steps x 5 array of posterior samples
- orbitize.hipparcos.HipparcosLogProb: the object storing relevant
metadata for the performed Hipparcos IAD fit
- Return type
tuple