Skip to content

Commit ba29e41

Browse files
committed
accepts different data types for MAST queries, tests, docs
1 parent 8a0a427 commit ba29e41

File tree

4 files changed

+157
-64
lines changed

4 files changed

+157
-64
lines changed

astroquery/mast/observations.py

Lines changed: 89 additions & 28 deletions
Original file line numberDiff line numberDiff line change
@@ -1121,6 +1121,63 @@ def service_request_async(self, service, params, *, pagesize=None, page=None, **
11211121

11221122
return self._portal_api_connection.service_request_async(service, params, pagesize, page, **kwargs)
11231123

1124+
def _normalize_filter_value(self, key: str, val) -> list:
1125+
"""
1126+
Normalize a filter value into a list suitable for MAST filters.
1127+
1128+
Parameters
1129+
----------
1130+
key : str
1131+
Parameter name (used for error messages).
1132+
val : any
1133+
Raw filter value.
1134+
1135+
Returns
1136+
-------
1137+
list
1138+
Normalized filter values.
1139+
"""
1140+
# Range filters must be dicts with 'min' and 'max'
1141+
if isinstance(val, dict):
1142+
if not {"min", "max"}.issubset(val.keys()):
1143+
raise InvalidQueryError(
1144+
f'Range filter for "{key}" must be a dictionary with "min" and "max" keys.'
1145+
)
1146+
return [val]
1147+
1148+
# Convert numpy arrays to lists
1149+
if isinstance(val, np.ndarray):
1150+
val = val.tolist()
1151+
1152+
# Convert numpy arrays, sets, or tuples to lists
1153+
if isinstance(val, (set, tuple)):
1154+
val = list(val)
1155+
1156+
# Wrap scalars into a list
1157+
return val if isinstance(val, list) else [val]
1158+
1159+
def _build_filters(self, service_params):
1160+
"""
1161+
Construct filters for filtered services.
1162+
1163+
Parameters
1164+
----------
1165+
service_params : dict
1166+
Parameters not classified as request/position keys.
1167+
1168+
Returns
1169+
-------
1170+
list of dict
1171+
Filters suitable for a MAST filtered query.
1172+
"""
1173+
filters = []
1174+
for key, val in service_params.items():
1175+
filters.append({
1176+
"paramName": key,
1177+
"values": self._normalize_filter_value(key, val)
1178+
})
1179+
return filters
1180+
11241181
def mast_query(self, service, columns=None, **kwargs):
11251182
"""
11261183
Given a Mashup service and parameters as keyword arguments, builds and excecutes a Mashup query.
@@ -1129,53 +1186,57 @@ def mast_query(self, service, columns=None, **kwargs):
11291186
----------
11301187
service : str
11311188
The Mashup service to query.
1132-
columns : str, optional
1189+
columns : str or list, optional
11331190
Specifies the columns to be returned as a comma-separated list, e.g. "ID, ra, dec".
11341191
**kwargs :
11351192
Service-specific parameters and MashupRequest properties. See the
11361193
`service documentation <https://mast.stsci.edu/api/v0/_services.html>`__ and the
11371194
`MashupRequest Class Reference <https://mast.stsci.edu/api/v0/class_mashup_1_1_mashup_request.html>`__
11381195
for valid keyword arguments.
11391196
1197+
For filtered services (i.e. those with "filtered" in the service name),
1198+
parameters that are not related to position or MashupRequest properties
1199+
are treated as filters. If the column has discrete values, the parameter value should be a
1200+
single value or list of values, and values will be matched exactly. If the column is continuous,
1201+
you can filter by a single value, a list of values, or a range of values. If filtering by a range of values,
1202+
the parameter value should be a dict in the form ``{'min': minVal, 'max': maxVal}``.
1203+
11401204
Returns
11411205
-------
11421206
response : `~astropy.table.Table`
11431207
"""
11441208
# Specific keywords related to positional and MashupRequest parameters.
1145-
position_keys = ['ra', 'dec', 'radius', 'position']
1146-
request_keys = ['format', 'data', 'filename', 'timeout', 'clearcache',
1147-
'removecache', 'removenullcolumns', 'page', 'pagesize']
1209+
position_keys = {'ra', 'dec', 'radius', 'position'}
1210+
request_keys = {'format', 'data', 'filename', 'timeout', 'clearcache',
1211+
'removecache', 'removenullcolumns', 'page', 'pagesize'}
11481212

1149-
# Explicit formatting for Mast's filtered services
1150-
if 'filtered' in service.lower():
1213+
# Split params into categories
1214+
position_params = {k: v for k, v in kwargs.items() if k.lower() in position_keys}
1215+
request_params = {k: v for k, v in kwargs.items() if k.lower() in request_keys}
1216+
service_params = {k: v for k, v in kwargs.items() if k.lower() not in position_keys | request_keys}
11511217

1152-
# Separating the filter params from the positional and service_request method params.
1153-
filters = [{'paramName': k, 'values': kwargs[k]} for k in kwargs
1154-
if k.lower() not in position_keys+request_keys]
1155-
position_params = {k: v for k, v in kwargs.items() if k.lower() in position_keys}
1156-
request_params = {k: v for k, v in kwargs.items() if k.lower() in request_keys}
1218+
# Handle filtered vs. non-filtered services
1219+
if 'filtered' in service.lower():
1220+
filters = self._build_filters(service_params)
11571221

1158-
# Mast's filtered services require at least one filter
1159-
if filters == []:
1160-
raise InvalidQueryError("Please provide at least one filter.")
1222+
if not filters:
1223+
raise InvalidQueryError('Please provide at least one filter.')
11611224

1162-
# Building 'params' for Mast.service_request
1163-
if columns is None:
1164-
columns = '*'
1225+
if columns is not None and isinstance(columns, list):
1226+
columns = ','.join(columns)
11651227

1166-
params = {'columns': columns,
1167-
'filters': filters,
1168-
**position_params
1169-
}
1228+
params = {
1229+
'columns': columns or '*',
1230+
'filters': filters,
1231+
**position_params,
1232+
}
11701233
else:
1171-
1172-
# Separating service specific params from service_request method params
1173-
params = {k: v for k, v in kwargs.items() if k.lower() not in request_keys}
1174-
request_params = {k: v for k, v in kwargs.items() if k.lower() in request_keys}
1175-
1176-
# Warning for wrong input
11771234
if columns is not None:
1178-
warnings.warn("'columns' parameter will not mask non-filtered services", InputWarning)
1235+
warnings.warn(
1236+
"'columns' parameter is ignored for non-filtered services.",
1237+
InputWarning
1238+
)
1239+
params = {**service_params, **position_params}
11791240

11801241
return self.service_request(service, params, **request_params)
11811242

astroquery/mast/tests/test_mast.py

Lines changed: 17 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -8,6 +8,7 @@
88
from unittest.mock import patch
99

1010
import pytest
11+
import numpy as np
1112

1213
from astropy.table import Table, unique
1314
from astropy.coordinates import SkyCoord
@@ -551,20 +552,30 @@ def test_mast_query(patch_post):
551552

552553
# filtered search
553554
result = mast.Mast.mast_query('Mast.Caom.Filtered',
554-
dataproduct_type=['image'],
555-
proposal_pi=['Osten, Rachel A.'],
556-
s_dec=[{'min': 43.5, 'max': 45.5}])
555+
dataproduct_type=['image', 'spectrum'],
556+
proposal_pi={'Osten, Rachel A.'},
557+
calib_level=np.asarray(3),
558+
s_dec={'min': 43.5, 'max': 45.5},
559+
columns=['proposal_pi', 's_dec', 'obs_id'])
557560
pp_list = result['proposal_pi']
558561
sd_list = result['s_dec']
559562
assert isinstance(result, Table)
560563
assert len(set(pp_list)) == 1
561564
assert max(sd_list) < 45.5
562565
assert min(sd_list) > 43.5
563566

564-
# error handling
565-
with pytest.raises(InvalidQueryError) as invalid_query:
567+
# warn if columns provided for non-filtered query
568+
with pytest.warns(InputWarning, match="'columns' parameter is ignored"):
569+
mast.Mast.mast_query('Mast.Caom.Cone', ra=23.34086, dec=60.658, radius=0.2, columns=['obs_id', 's_ra'])
570+
571+
# error if no filters provided for filtered query
572+
with pytest.raises(InvalidQueryError, match="Please provide at least one filter."):
566573
mast.Mast.mast_query('Mast.Caom.Filtered')
567-
assert "Please provide at least one filter." in str(invalid_query.value)
574+
575+
# error if a full range if not provided for range filter
576+
with pytest.raises(InvalidQueryError,
577+
match='Range filter for "s_ra" must be a dictionary with "min" and "max" keys.'):
578+
mast.Mast.mast_query('Mast.Caom.Filtered', s_ra={'min': 10.0})
568579

569580

570581
def test_resolve_object_single(patch_post):

astroquery/mast/tests/test_mast_remote.py

Lines changed: 20 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -417,17 +417,30 @@ def test_mast_service_request(self):
417417
assert len(result) == 10
418418

419419
def test_mast_query(self):
420-
result = Mast.mast_query('Mast.Caom.Cone', ra=184.3, dec=54.5, radius=0.2)
421-
422-
# Is result in the right format
420+
# Cone search (unfiltered)
421+
result = Mast.mast_query('Mast.Caom.Cone', ra=184.3, dec=54.5, radius=0.005)
423422
assert isinstance(result, Table)
424-
425-
# Are the GALEX observations in the results table
426423
assert "GALEX" in result['obs_collection']
427-
428-
# Are the two GALEX observations with obs_id 6374399093149532160 in the results table
429424
assert len(result[np.where(result["obs_id"] == "6374399093149532160")]) == 2
430425

426+
# Filtered query
427+
columns = ['target_name', 'obs_collection', 'calib_level', 'sequence_number', 't_min']
428+
result = Mast.mast_query('Mast.Caom.Filtered',
429+
target_name=375422201,
430+
obs_collection={'TESS'},
431+
calib_level=np.asarray(3),
432+
sequence_number=[15, 16],
433+
t_min={'min': 58710, 'max': 58720},
434+
columns=columns)
435+
assert isinstance(result, Table)
436+
assert all(result['target_name'] == '375422201')
437+
assert all(result['obs_collection'] == 'TESS')
438+
assert all(result['calib_level'] == 3)
439+
assert all((result['sequence_number'] == 15) | (result['sequence_number'] == 16))
440+
assert (result['t_min'] >= 58710).all() and (result['t_min'] <= 58720).all()
441+
assert all(c in list(result.columns.keys()) for c in columns)
442+
assert len(result.columns) == 5
443+
431444
def test_mast_session_info(self):
432445
sessionInfo = Mast.session_info(verbose=False)
433446
assert sessionInfo['ezid'] == 'anonymous'

docs/mast/mast_mastquery.rst

Lines changed: 31 additions & 23 deletions
Original file line numberDiff line numberDiff line change
@@ -7,7 +7,7 @@ The Mast class provides more direct access to the MAST interface. It requires
77
more knowledge of the inner workings of the MAST API, and should be rarely
88
needed. However in the case of new functionality not yet implemented in
99
astroquery, this class does allow access. See the
10-
`MAST api documentation <https://mast.stsci.edu/api/v0/>`__ for more
10+
`MAST API documentation <https://mast.stsci.edu/api/v0/>`__ for more
1111
information.
1212

1313
The basic MAST query function allows users to query through the following
@@ -18,12 +18,17 @@ their corresponding parameters and returns query results as an
1818
Filtered Mast Queries
1919
=====================
2020

21-
MAST's Filtered services use the parameters 'columns' and 'filters'. The 'columns'
22-
parameter is a required string that specifies the columns to be returned as a
23-
comma-separated list. The 'filters' parameter is a required list of filters to be
24-
applied. The `~astroquery.mast.MastClass.mast_query` method accepts that list of
25-
filters as keyword arguments paired with a list of values, similar to
26-
`~astroquery.mast.ObservationsClass.query_criteria`.
21+
MAST's filtered services (i.e. those with "filtered" in the service name) accept service-specific parameters, MashupRequest properties, and
22+
column filters as keyword arguments and return a table of matching observations. See the
23+
`service documentation <https://mast.stsci.edu/api/v0/_services.html>`__ and the
24+
`MashupRequest Class Reference <https://mast.stsci.edu/api/v0/class_mashup_1_1_mashup_request.html>`__ for valid keyword arguments.
25+
26+
Parameters that are not related to position or MashupRequest properties are treated as column filters.
27+
If the column has discrete values, the parameter value should be a single value or list of values, and values will be matched exactly.
28+
If the column is continuous, you can filter by a single value, a list of values, or a range of values. If filtering by a range of values,
29+
the parameter value should be a dict in the form ``{'min': minVal, 'max': maxVal}``.
30+
31+
The ``columns`` parameter specifies the columns to be returned in the response as a comma-separated string or list of strings.
2732

2833
The following example uses a JWST service with column names and filters specific to
2934
JWST services. For the full list of valid parameters view the
@@ -34,22 +39,25 @@ JWST services. For the full list of valid parameters view the
3439
>>> from astroquery.mast import Mast
3540
...
3641
>>> observations = Mast.mast_query('Mast.Jwst.Filtered.Nirspec',
37-
... columns='title, instrume, targname',
38-
... targoopp=['T'])
42+
... targoopp='T',
43+
... productLevel=['2a', '2b'],
44+
... duration={'min': 810, 'max': 820},
45+
... columns=['filename', 'targoopp', 'productLevel', 'duration'])
3946
>>> print(observations) # doctest: +IGNORE_OUTPUT
40-
title instrume targname
41-
------------------------------- -------- ----------------
42-
ToO Comet NIRSPEC ZTF (C/2022 E3)
43-
ToO Comet NIRSPEC ZTF (C/2022 E3)
44-
ToO Comet NIRSPEC ZTF (C/2022 E3)
45-
ToO Comet NIRSPEC ZTF (C/2022 E3)
46-
De-Mystifying SPRITEs with JWST NIRSPEC SPIRITS18nu
47-
ToO Comet NIRSPEC ZTF (C/2022 E3)
48-
... ... ...
49-
ToO Comet NIRSPEC ZTF (C/2022 E3)
50-
ToO Comet NIRSPEC ZTF (C/2022 E3)
51-
ToO Comet NIRSPEC ZTF (C/2022 E3)
52-
Length = 319 rows
47+
filename targoopp productLevel duration
48+
-------------------------------------------- -------- ------------ --------
49+
jw05324004001_03102_00004_nrs2_rate.fits t 2a 816.978
50+
jw05324004001_03102_00004_nrs2_rateints.fits t 2a 816.978
51+
jw05324004001_03102_00001_nrs2_rate.fits t 2a 816.978
52+
jw05324004001_03102_00001_nrs2_rateints.fits t 2a 816.978
53+
jw05324004001_03102_00005_nrs2_rate.fits t 2a 816.978
54+
... ... ... ...
55+
jw05324004001_03102_00003_nrs1_s2d.fits t 2b 816.978
56+
jw05324004001_03102_00003_nrs1_x1d.fits t 2b 816.978
57+
jw05324004001_03102_00002_nrs1_cal.fits t 2b 816.978
58+
jw05324004001_03102_00002_nrs1_s2d.fits t 2b 816.978
59+
jw05324004001_03102_00002_nrs1_x1d.fits t 2b 816.978
60+
Length = 25 rows
5361

5462

5563
TESS Queries
@@ -181,7 +189,7 @@ result in a warning.
181189
... ra=254.287,
182190
... dec=-4.09933,
183191
... radius=0.02) # doctest: +SHOW_WARNINGS
184-
InputWarning: 'columns' parameter will not mask non-filtered services
192+
InputWarning: 'columns' parameter is ignored for non-filtered services.
185193

186194
Advanced Service Request
187195
========================

0 commit comments

Comments
 (0)