Integration with Pandas
py_outfit.pandas_pyoutfit
Pandas accessor for batch Initial Orbit Determination (IOD) with Outfit.
This module adds a DataFrame.outfit
accessor exposing a vectorized entry-point
to run Gauss IOD on a flat, columnar table of astrometric measurements.
Examples:
>>> df = pd.DataFrame({
... "tid": [0, 0, 0, 1, 1, 1],
... "mjd": [60000.0, 60000.01, 60000.02, 60000.0, 60000.01, 60000.02],
... "ra": [10.1, 10.2, 10.3, 33.4, 33.5, 33.6], # degrees by default
... "dec": [-5.0, -4.9, -4.8, 2.0, 2.1, 2.2], # degrees by default
... })
>>> env = PyOutfit(ephem="horizon:DE440", error_model="FCCT14")
>>> params = IODParams() # default settings
>>> obs = PyOutfit.get_observer_from_mpc_code("I41") # ZTF
>>> orb_pdf = pandas_traj.outfit.estimate_orbits(
... env, params, obs, ra_error=0.5, dec_error=0.5
... )
>>> orb_pdf.columns[:5]
Index(['object_id', 'variant', 'element_set', 'rms', ...], dtype='object')
Design
- Schema: maps your column names to the expected fields (
tid/mjd/ra/dec
). - Vector ingestion: we assemble a
TrajectorySet
directly from NumPy arrays. - Units:
units="degrees"
(default) interprets RA/DEC in degrees and uncertainties in arcseconds;"radians"
expects everything in radians.
Schema
dataclass
Column schema used to map a DataFrame to py_outfit's vector ingestion.
ATTRIBUTE | DESCRIPTION |
---|---|
tid |
Column containing object/trajectory IDs (int or str). Repeated per-row.
TYPE:
|
mjd |
Modified Julian Date (TT, days).
TYPE:
|
ra |
Right Ascension values (degrees if
TYPE:
|
dec |
Declination values (degrees if
TYPE:
|
Notes
- Only these four columns are required for the accessor.
- If your DataFrame uses different names, override the defaults:
schema=Schema(tid="object", mjd="epoch", ra="alpha", dec="delta")
.
OutfitAccessor
Pandas accessor for running Gauss IOD from a DataFrame.
Use via the attribute accessor DataFrame.outfit
. It exposes
:meth:estimate_orbits
to run a vectorized Initial Orbit Determination over
a flat table of astrometric measurements.
Examples:
Basic usage with degrees and arcseconds
Radians workflow
>>> out = df.outfit.estimate_orbits(
... env, params, observer,
... ra_error=1e-6, dec_error=1e-6,
... units="radians",
... )
estimate_orbits
estimate_orbits(
env: PyOutfit,
params: IODParams,
observer: Observer,
ra_error: float,
dec_error: float,
*,
schema: Schema = Schema(),
units: Literal["degrees", "radians"] = "degrees",
rng_seed: Optional[int] = None
) -> DataFrame
Run Gauss IOD on a flat astrometry table and return a one-row-per-object summary.
The input DataFrame must at least provide a trajectory/object identifier,
Modified Julian Date in TT, and right ascension/declination angles. The
names of these columns are defined by schema
(defaults are tid
, mjd
,
ra
, dec
).
PARAMETER | DESCRIPTION |
---|---|
env
|
Configured computation engine, including ephemerides, the error model and the available observers.
TYPE:
|
params
|
IOD configuration controlling triplet search, Monte Carlo perturbations, filters, and tolerances.
TYPE:
|
observer
|
Default observer applied to all rows. This covers the common single-station use case.
TYPE:
|
ra_error
|
Uncertainty on right ascension. When
TYPE:
|
dec_error
|
Uncertainty on declination. Follows the same unit convention as
TYPE:
|
schema
|
Column mapping for the current DataFrame. Use this to adapt to
non-standard column names. The default expects |
units
|
Angle units for
TYPE:
|
rng_seed
|
Optional seed to make randomized internals deterministic.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
DataFrame
|
A summary DataFrame with one row per object containing the RMS value,
the detected orbital element set, the orbit variant, and the native
orbital elements returned by the engine. The |
RAISES | DESCRIPTION |
---|---|
ValueError
|
Raised when required columns are missing in the DataFrame, or when
|
Notes
Ingestion is vectorized and avoids Python-level grouping. For
units="degrees"
, ra_error
and dec_error
are converted from
arcseconds to radians using RADSEC
. If needed, multi-observer or
per-row observatory support can be added by extending the
TrajectorySet.from_numpy_*
builder to accept vectorized observers.
See Also
-
TrajectorySet.estimate_all_orbits
: Batch Gauss IOD over trajectories. -
GaussResult.to_dict
: Native orbital element names used in the output. -
Schema
: Column mapping used to interpret the input DataFrame.