Orbits¶
If they represent binary or triple star systems,
vespa.stars.StarPopulation objects are created with a large
population of randomized orbits. This is done using
the OrbitPopulation and TripleOrbitPopulation
objects.
The engine that makes it possible to initialize large numbers of
random orbital positions nearly instantaneously is
the kepler.Efn() function (as used
by utils.orbit_posvel()), which uses a precomputed grid to
interpolate the solutions to Kepler’s equation for a given mean
anomaly and eccentricity (or arrays thereof).
The final coordinate system of these populations is
“observer-oriented,” with the z axis along the line of sight, and
the x-y plane being the plane of the sky. Practically, this is
accomplished by first simulating all the random orbits in the x-y
plane, and then “observing” them from lines of sight randomly oriented
on the unit sphere, and projecting appropriately.
Coordinates are handled using astropy.coordinates.SkyCoord
objects.
Orbit Populations¶
-
class
vespa.orbits.populations.OrbitPopulation(M1, M2, P, ecc=0, n=None, mean_anomaly=None, obsx=None, obsy=None, obsz=None, obspos=None)[source]¶ Population of orbits.
Parameters: - M1,M2 – Primary and secondary masses (if not
Quantity, assumed to be in solar masses). Can befloat, array-like orQuantity. - P (
float, array-like orQuantity.) – Orbital period(s) (if notQuantity, assumed to be in days) - ecc – (
floator array-like, optional) Eccentricities. - n – (optional)
Number of instances to simulate. If not provided, then this number
will be the length of
M2(orP) provided. - mean_anomaly – (optional)
Mean anomalies of orbits. Usually this will just be set randomly,
but can be provided to initialize a particular state (e.g., when
restoring an
OrbitPopulationobject from a saved state). - obsy, obsz (obsx,) – (optional) “Observer” positions to define coordinates. Will be set randomly if not provided.
- obspos (
astropy.coordinates.SkyCoord) – (optional) “Observer” positions may be set with aSkyCoordobject (replaces obsx, obsy, obsz)
-
RV¶ Relative radial velocities of two stars
-
RV_com1¶ RVs of star 1 relative to center-of-mass
-
RV_com2¶ RVs of star 2 relative to center-of-mass
-
RV_timeseries(ts, recalc=False)[source]¶ Radial Velocity time series for star 1 at given times ts.
Parameters: - ts (array-like or
Quantity) – Times. If notQuantity, assumed to be in days. - recalc – (optional)
If
False, then if called with the exact sametsas last call, it will return cached calculation.
- ts (array-like or
-
Rsky¶ Sky separation of stars, in projected AU.
-
dRV(dt, com=False)[source]¶ Change in RV of star 1 for time separation dt (default=days)
Parameters: - dt (float, array-like, or
Quantity) – Time separation for which to compute RV change. If not aQuantity, then assumed to be in days. - com – (
bool, optional) IfTrue, then return dRV of star 1 in center-of-mass frame.
Return dRV: Change in radial velocity over time
dt.- dt (float, array-like, or
-
dataframe¶ Summary DataFrame of OrbitPopulation
Used to save/restore state.
-
classmethod
from_df(df)[source]¶ Creates an OrbitPopulation from a DataFrame.
Parameters: df – pandas.DataFrameobject. Must contain the following columns:['M1','M2','P','ecc','mean_anomaly','obsx','obsy','obsz'], i.e., as what is accessed viaOrbitPopulation.dataframe.Returns: OrbitPopulation.
-
classmethod
load_hdf(filename, path='')[source]¶ Loads OrbitPopulation from HDF file.
Parameters: - filename – HDF file
- path – Path within HDF file store where
OrbitPopulationis saved.
-
scatterplot(fig=None, figsize=(7, 7), ms=0.5, rmax=None, log=False, **kwargs)[source]¶ Makes a scatter plot of projected X-Y sky separation
Parameters: - fig – (optional)
Passed to
plotutils.setfig() - figsize – (optional) Size of figure (in).
- ms – (optional) Marker size
- rmax – (optional) Maximum projected radius to plot.
- log – (optional) Whether to plot with log scale.
- **kwargs –
Additional keyword arguments passed to
plt.plot.
- fig – (optional)
Passed to
- M1,M2 – Primary and secondary masses (if not
-
class
vespa.orbits.populations.TripleOrbitPopulation(M1, M2, M3, Plong, Pshort, ecclong=0, eccshort=0, n=None, mean_anomaly_long=None, obsx_long=None, obsy_long=None, obsz_long=None, obspos_long=None, mean_anomaly_short=None, obsx_short=None, obsy_short=None, obsz_short=None, obspos_short=None)[source]¶ Stars 2 and 3 orbit each other (short orbit), far from star 1 (long orbit)
This object defines the orbits of a triple star system, with orbits calculated assuming the “long” orbit does not perturb the “short” orbit, which will not be true in the long run, but should be true over short timescales as long as
Plong >> Pshort.A
TripleOrbitPopulationis essentially made up of twoOrbitPopulationobjects: one for the “long” orbit and one for the “short.”Parameters: - M1,M2,M3 – Masses of stars. Stars 2 and 3 are in a short orbit, far away from star 1.
If not
astropy.units.Quantityobjects, then assumed to be in solar mass units. May be single value or array-like. - Plong,Pshort – Orbital Periods. Plong is orbital period of 2+3 and 1; Pshort is orbital
period of 2 and 3. If not
astropy.units.Quantityobjects, assumed to be in days. Can be single value or array-like. N.B. If any item in Pshort happens to be longer than the corresponding item in Plong, they will be switched. - ecclong,eccshort – (optional) Eccentricities. Same story (long vs. short). Default=0 (circular). Can be single value or array-like.
- n – (optional)
Number of systems to simulate (if
M1,M2,M3aren’t arrays of size > 1 already). - mean_anomaly_short,mean_anomaly_long – (optional)
Mean anomalies. This is only passed if you need to restore a
particular specific configuration (i.e., a particular saved simulation),
e.g., as done by
TripleOrbitPopulation.from_df(). If not provided, then randomized on (0, 2pi). - obsx_short,obsy_short,obsz_short – (optional) “Observer” positions for the short orbit. Also only passed for purposes of restoring configuration.
- obsx_long,obsy_long,obsz_long – (optional) “Observer” positions for long orbit. Also only passed for purposes of restoring configuration.
- obspos_short,obspos_long – (optional)
“Observer positions for short and long orbits, provided
as
astropy.SkyCoordobjects. These will replace obsx_short/long, obsy_short/long, obsz_short/long parameters if present.
-
RV¶ Instantaneous RV of star 1 with respect to system center-of-mass
-
RV_1¶ Instantaneous RV of star 1 with respect to system center-of-mass
-
RV_2¶ Instantaneous RV of star 2 with respect to system center-of-mass
-
RV_3¶ Instantaneous RV of star 3 with respect to system center-of-mass
-
Rsky¶ Projected separation of star 2+3 pair from star 1 [projected AU]
-
dRV(dt)[source]¶ Returns difference in RVs (separated by time dt) of star 1.
Parameters: dt – Time separation for which to compute RV change. If not an astropy.units.Quantityobject, then assumed to be in days.
-
dRV_1(dt)[source]¶ Returns difference in RVs (separated by time dt) of star 1.
Parameters: dt – Time separation for which to compute RV change. If not an astropy.units.Quantityobject, then assumed to be in days.
-
dRV_2(dt)[source]¶ Returns difference in RVs (separated by time dt) of star 2.
Parameters: dt – Time separation for which to compute RV change. If not an astropy.units.Quantityobject, then assumed to be in days.
-
dRV_3(dt)[source]¶ Returns difference in RVs (separated by time dt) of star 3.
Parameters: dt – Time separation for which to compute RV change. If not an astropy.units.Quantityobject, then assumed to be in days.
-
classmethod
from_df(df_long, df_short)[source]¶ Builds TripleOrbitPopulation from DataFrame
DataFrameobjects must be of appropriate form to pass toOrbitPopulation.from_df().Parameters: df_short (df_long,) – pandas.DataFrameobjects to pass toOrbitPopulation.from_df().
- M1,M2,M3 – Masses of stars. Stars 2 and 3 are in a short orbit, far away from star 1.
If not
Utility Functions¶
The following functions are used in the creation
of OrbitPopulation objects. kepler.Efn() is used for
instanteous solution of Kepler’s equation (via interpolation),
and utils.orbit_posvel() does the projecting of random orbits
into 3-d Cartesian coordinates, assisted by
utils.orbitproject() and utils.random_spherepos().
-
vespa.orbits.kepler.Efn(Ms, eccs)[source]¶ Returns Eccentric anomaly, interpolated from pre-computed grid of M, ecc
Instantaneous solution of Kepler’s equation!
Works for
-2*np.pi < Ms < 2*np.piandeccs <= 0.97Parameters: - Ms – (
floator array-like) Mean anomaly - eccs – (
floator array-like)
- Ms – (
-
vespa.orbits.utils.orbit_posvel(Ms, eccs, semimajors, mreds, obspos=None)[source]¶ Returns positions in projected AU and velocities in km/s for given mean anomalies.
Returns 3-D positions and velocities as SkyCoord objects, in “observer” reference frame. Uses
kepler.Efn()to calculate eccentric anomalies using interpolation.Parameters: - Ms,eccs,semimajors,mreds – (
floator array-like) Mean anomalies, eccentricities, semimajor axes [AU], reduced masses [Msun]. - obspos – (
None,(x,y,z)tupleorSkyCoordobject) Locations of observers for which to return coordinates. IfNonethen populate randomly on sphere. If(x,y,z)orSkyCoordobject provided, then use those.
Returns pos,vel: SkyCoordObjects representing the positions and velocities, the coordinates of which areQuantityobjects that have units. Positions are in projected AU and velocities in km/s.- Ms,eccs,semimajors,mreds – (
-
vespa.orbits.utils.orbitproject(x, y, inc, phi=0, psi=0)[source]¶ Transform x,y planar coordinates into observer’s coordinate frame.
x,yare coordinates inz=0plane (plane of the orbit)observer is at
(inc, phi)on celestial sphere (angles in radians);psiis orientation of finalx-yaxes about the(inc,phi)vector.Returns
x,y,zvalues in observer’s coordinate frame, wherex,yare now plane-of-sky coordinates andzis along the line of sight.Parameters: - x,y – (
floator array-like) Coordinates to transform. - inc – (
floator array-like) Polar angle(s) of observer (whereinc=0corresponds to north pole of originalx-yplane). This angle is the same as standard “inclination.” - phi – (
floator array-like, optional) Azimuthal angle of observer aroundz-axis - psi – (
floator array-like, optional) Orientation of final observer coordinate frame (azimuthal around(inc,phi)vector.
Return x,y,z: (
ndarray) Coordinates in observers’ frames.x,yin “plane of sky” andzalong line of sight.- x,y – (