|
lsst.pipe.tasks gcf790cdeb6+e07a3617c0
|
Functions | |
| assoc_validate (dia, assoc) | |
| packed_ascii_to_uint64_le (mpc_packed) | |
| solar_elongation_ndarray (ra_deg, dec_deg, t) | |
| group_by (arrs, key, func, out=None, check_grouped=True) | |
| bool | values_grouped (np.ndarray a) |
| EarthLocation | earthlocation_from_obscode (str obscode) |
| observatory_barycentric_posvel (str obscode, Time obstime) | |
| None | struct_to_parquet (np.ndarray arr, str path, *, Optional[int] chunk_size=None, Optional[int] row_group_size=None) |
| tisserand_jupiter (a, e, inc_deg, a_j=A_JUP) | |
| unpack (df, to_numpy=True) | |
| argjoin (a, v) | |
Variables | |
| float | A_JUP = 5.2044 |
| lsst.pipe.tasks.ssp.util.argjoin | ( | a, | |
| v ) |
Perform an efficient inner join between two 1-D NumPy arrays, returning
the index pairs that match by value.
Parameters
----------
a : ndarray
The left-hand array to join on. Must be 1-dimensional.
v : ndarray
The right-hand array to join on. Must be 1-dimensional.
Returns
-------
aidx : ndarray (int)
Indices into `a` selecting the rows that participate in the join.
vidx : ndarray (int)
Indices into `v` selecting the corresponding matching rows.
After the join:
a[aidx] == v[vidx]
is guaranteed to be true for all elements.
Notes
-----
This function implements a pure NumPy equivalent of an SQL-style
INNER JOIN on the key columns `a` and `v`.
The algorithm:
1. Sort `a` to produce a permutation `i` so that `a[i]` is sorted.
2. Use `np.searchsorted(a[i], v)` to find, for each element of `v`,
the candidate matching location in the sorted array.
3. Map these positions back to the coordinates of the original array `a`
using the permutation `i`.
4. Filter out non-matches (values in `v` not present in `a`).
The remaining pairs form the inner join.
Complexity
----------
Sorting: O(len(a) log len(a))
Searching: O(len(v) log len(a))
Total: O(n log n)
This is optimal for join-like operations on unsorted arrays in NumPy.
Examples
--------
>>> a = np.array(["b", "a", "c", "b"])
>>> v = np.array(["a", "b", "x", "b"])
>>> aidx, vidx = argjoin(a, v)
>>> a[aidx]
array(['a', 'b', 'b'])
>>> v[vidx]
array(['a', 'b', 'b'])
| EarthLocation lsst.pipe.tasks.ssp.util.earthlocation_from_obscode | ( | str | obscode | ) |
| lsst.pipe.tasks.ssp.util.group_by | ( | arrs, | |
| key, | |||
| func, | |||
| out = None, | |||
| check_grouped = True ) |
Group multiple NumPy arrays by arrs[0][key], assuming the key column
is already grouped (equal keys are contiguous), but group blocks may
appear in any order.
If out is provided:
func(row_view, *subarrs)
If out is None:
func(*subarrs)
Parameters
----------
arrs : list/tuple of ndarray
Equal-length arrays. arrs[0] contains the grouping key.
key : str
Column name in arrs[0] to group by.
func : callable
Called either as func(row, *subarrs) or func(*subarrs).
out : ndarray or None
Optional preallocated structured output array.
check_grouped : bool
If True, verify that key values are grouped contiguously.
Returns
-------
ndarray or dict
| lsst.pipe.tasks.ssp.util.observatory_barycentric_posvel | ( | str | obscode, |
| Time | obstime ) |
Barycentric (ICRS) position and velocity of an observatory given an
MPC obscode, using JPL DE440 for the Earth ephemeris.
Returns
-------
r_bary : Quantity, shape (3, ...)
Barycentric position in AU.
v_bary : Quantity, shape (3, ...)
Barycentric velocity in AU/day.
| lsst.pipe.tasks.ssp.util.packed_ascii_to_uint64_le | ( | mpc_packed | ) |
| lsst.pipe.tasks.ssp.util.solar_elongation_ndarray | ( | ra_deg, | |
| dec_deg, | |||
| t ) |
Very fast computation of solar elongation (ICRS great-circle separation)
using astropy.coordinates.angular_separation.
Parameters
----------
ra_deg : ndarray
Target RA in degrees (ICRS).
dec_deg : ndarray
Target Dec in degrees (ICRS).
t : astropy.time.Time
Observation times.
Returns
-------
elong_deg : ndarray
Solar elongation in degrees.
| None lsst.pipe.tasks.ssp.util.struct_to_parquet | ( | np.ndarray | arr, |
| str | path, | ||
| * | , | ||
| Optional[int] | chunk_size = None, | ||
| Optional[int] | row_group_size = None ) |
| lsst.pipe.tasks.ssp.util.tisserand_jupiter | ( | a, | |
| e, | |||
| inc_deg, | |||
| a_j = A_JUP ) |
Compute Tisserand parameter with respect to Jupiter.
Parameters
----------
a : float or ndarray
Semimajor axis of the small body [AU].
e : float or ndarray
Eccentricity.
inc_deg : float or ndarray
Inclination [degrees], typically to the ecliptic.
a_j : float
Semimajor axis of Jupiter [AU]. Default ~5.2044.
Returns
-------
T_J : float or ndarray
Tisserand parameter with respect to Jupiter.
| lsst.pipe.tasks.ssp.util.unpack | ( | df, | |
| to_numpy = True ) |
Return all DataFrame columns as a tuple.
Parameters
----------
df : pandas.DataFrame
Input dataframe.
to_numpy : bool, default True
If True, return each column as a NumPy array.
If False, return each column as a pandas Series.
Returns
-------
tuple
Tuple of columns in the original order.
| bool lsst.pipe.tasks.ssp.util.values_grouped | ( | np.ndarray | a | ) |