lsst.pipe.tasks gcf790cdeb6+0604939b8f
Loading...
Searching...
No Matches
lsst.pipe.tasks.ssp.util Namespace Reference

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
 

Function Documentation

◆ argjoin()

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'])

Definition at line 439 of file util.py.

◆ assoc_validate()

lsst.pipe.tasks.ssp.util.assoc_validate ( dia,
assoc )

Definition at line 22 of file util.py.

◆ earthlocation_from_obscode()

EarthLocation lsst.pipe.tasks.ssp.util.earthlocation_from_obscode ( str obscode)
Convert an MPC observatory code (e.g. 'X05') to an EarthLocation,
using MPC.get_observatory_codes() columns:
  Code, Longitude, cos, sin, Name.

Definition at line 224 of file util.py.

◆ group_by()

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

Definition at line 124 of file util.py.

◆ observatory_barycentric_posvel()

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.

Definition at line 257 of file util.py.

◆ packed_ascii_to_uint64_le()

lsst.pipe.tasks.ssp.util.packed_ascii_to_uint64_le ( mpc_packed)
Convert a pandas string[pyarrow] column of ASCII strings (<= 8 bytes)
to little-endian uint64 by left-padding with ASCII spaces to 8 chars.

Definition at line 54 of file util.py.

◆ solar_elongation_ndarray()

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.

Definition at line 84 of file util.py.

◆ struct_to_parquet()

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 )
Write a large NumPy structured array to a Parquet file using PyArrow.

Designed for dtypes like dia_dtype / orbit_dtype with up to ~1e8 rows.

Definition at line 303 of file util.py.

◆ tisserand_jupiter()

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.

Definition at line 392 of file util.py.

◆ unpack()

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.

Definition at line 416 of file util.py.

◆ values_grouped()

bool lsst.pipe.tasks.ssp.util.values_grouped ( np.ndarray a)
Return True if each distinct value in 1D array `a`
appears in a single contiguous block (all duplicates grouped).

Definition at line 202 of file util.py.

Variable Documentation

◆ A_JUP

float lsst.pipe.tasks.ssp.util.A_JUP = 5.2044

Definition at line 389 of file util.py.