Skip to content

Commit 71b00bb

Browse files
authored
Merge pull request #3076 from bsipocz/ENH_irsa_ssa
ENH: adding SSA method for IRSA
2 parents db6b28c + eecd8e1 commit 71b00bb

File tree

5 files changed

+138
-6
lines changed

5 files changed

+138
-6
lines changed

CHANGES.rst

+2
Original file line numberDiff line numberDiff line change
@@ -39,6 +39,8 @@ gaia
3939
ipac.irsa
4040
^^^^^^^^^
4141

42+
- Method to run Simple Spectral Access (SSA) VO queries, ``query_ssa``, is added. [#3076]
43+
4244
- Adding the "servicetype" kwarg to ``list_collections`` to be able to list SIA
4345
and SSA collections separately. [#3200]
4446

astroquery/ipac/irsa/__init__.py

+1
Original file line numberDiff line numberDiff line change
@@ -28,6 +28,7 @@ class Conf(_config.ConfigNamespace):
2828
60,
2929
'Time limit for connecting to the IRSA server.')
3030
sia_url = _config.ConfigItem('https://irsa.ipac.caltech.edu/SIA', 'IRSA SIA URL')
31+
ssa_url = _config.ConfigItem('https://irsa.ipac.caltech.edu/SSA', 'IRSA SSA URL')
3132
tap_url = _config.ConfigItem('https://irsa.ipac.caltech.edu/TAP', 'IRSA TAP URL')
3233

3334

astroquery/ipac/irsa/core.py

+45-4
Original file line numberDiff line numberDiff line change
@@ -12,10 +12,8 @@
1212
from astropy import units as u
1313
from astropy.utils.decorators import deprecated_renamed_argument
1414

15-
from pyvo.dal import TAPService
16-
17-
from pyvo.dal.sia2 import SIA2Service, SIA2_PARAMETERS_DESC
18-
15+
from pyvo.dal import TAPService, SIA2Service, SSAService
16+
from pyvo.dal.sia2 import SIA2_PARAMETERS_DESC
1917
from astroquery import log
2018
from astroquery.query import BaseVOQuery
2119
from astroquery.utils.commons import parse_coordinates
@@ -31,8 +29,10 @@ class IrsaClass(BaseVOQuery):
3129
def __init__(self):
3230
super().__init__()
3331
self.sia_url = conf.sia_url
32+
self.ssa_url = conf.ssa_url
3433
self.tap_url = conf.tap_url
3534
self._sia = None
35+
self._ssa = None
3636
self._tap = None
3737

3838
@property
@@ -41,6 +41,12 @@ def sia(self):
4141
self._sia = SIA2Service(baseurl=self.sia_url, session=self._session)
4242
return self._sia
4343

44+
@property
45+
def ssa(self):
46+
if not self._ssa:
47+
self._ssa = SSAService(baseurl=self.ssa_url, session=self._session)
48+
return self._ssa
49+
4450
@property
4551
def tap(self):
4652
if not self._tap:
@@ -122,6 +128,41 @@ def query_sia(self, *, pos=None, band=None, time=None, pol=None,
122128

123129
query_sia.__doc__ = query_sia.__doc__.replace('_SIA2_PARAMETERS', SIA2_PARAMETERS_DESC)
124130

131+
def query_ssa(self, *, pos=None, radius=None, band=None, time=None, collection=None):
132+
"""
133+
Use standard SSA attributes to query the IRSA SSA service.
134+
135+
Parameters
136+
----------
137+
pos : `~astropy.coordinates.SkyCoord` class or sequence of two floats
138+
the position of the center of the circular search region.
139+
assuming icrs decimal degrees if unit is not specified.
140+
raidus : `~astropy.units.Quantity` class or scalar float
141+
the radius of the circular region around pos in which to search.
142+
assuming icrs decimal degrees if unit is not specified.
143+
band : `~astropy.units.Quantity` class or sequence of two floats
144+
the bandwidth range the observations belong to.
145+
assuming meters if unit is not specified.
146+
time : `~astropy.time.Time` class or sequence of two strings
147+
the datetime range the observations were made in.
148+
assuming iso 8601 if format is not specified.
149+
collection : str
150+
Name of the collection that the data belongs to.
151+
152+
Returns
153+
-------
154+
Results in `pyvo.dal.SSAResults` format.
155+
result.to_table() in Astropy table format
156+
"""
157+
158+
if radius is None:
159+
diameter = None
160+
else:
161+
diameter = 2 * radius
162+
163+
return self.ssa.search(pos=pos, diameter=diameter, band=band, time=time,
164+
format='all', collection=collection)
165+
125166
def list_collections(self, servicetype=None):
126167
"""
127168
Return information of available IRSA SIAv2 collections to be used in ``query_sia`` queries.

astroquery/ipac/irsa/tests/test_irsa_remote.py

+7
Original file line numberDiff line numberDiff line change
@@ -87,3 +87,10 @@ def test_tap(self):
8787
result = Irsa.query_tap(query=query)
8888
assert len(result) == 5
8989
assert result.to_table().colnames == ['ra', 'dec']
90+
91+
def test_ssa(self):
92+
coord = SkyCoord.from_name("Eta Carina")
93+
result = Irsa.query_ssa(pos=coord)
94+
assert len(result) > 260
95+
collections = set(result.to_table()['dataid_collection'])
96+
assert {'champ', 'iso_sws', 'sofia_forcast', 'sofia_great', 'spitzer_sha'}.issubset(collections)

docs/ipac/irsa/irsa.rst

+83-2
Original file line numberDiff line numberDiff line change
@@ -194,8 +194,8 @@ the query is send in asynchronous mode.
194194
>>> from astroquery.ipac.irsa import Irsa
195195
>>> table = Irsa.query_region("HIP 12", catalog="allwise_p3as_psd", spatial="Cone", async_job=True)
196196
>>> print(table)
197-
designation ra dec sigra ... y z spt_ind htm20
198-
deg deg arcsec ...
197+
designation ra dec sigra ... y z spt_ind htm20
198+
deg deg arcsec ...
199199
------------------- --------- ----------- ------ ... ------------------ ------------------- --------- -------------
200200
J000009.78-355736.9 0.0407905 -35.9602605 0.0454 ... 0.0005762523295116 -0.5872239888098030 100102010 8873706189183
201201

@@ -316,6 +316,87 @@ Now plot the cutout.
316316
plt.show()
317317

318318

319+
320+
Simple spectral access queries
321+
------------------------------
322+
323+
`~astroquery.ipac.irsa.IrsaClass.query_ssa` provides a way to access IRSA's Simple
324+
Spectral Access VO service. In the following example we are looking for Spitzer
325+
Enhanced Imaging products in the centre of the COSMOS field as an `~astropy.table.Table`.
326+
327+
.. doctest-remote-data::
328+
329+
>>> from astroquery.ipac.irsa import Irsa
330+
>>> from astropy.coordinates import SkyCoord
331+
>>> from astropy import units as u
332+
>>>
333+
>>> coord = pos = SkyCoord.from_name('Arp 220')
334+
>>> arp220_spectra = Irsa.query_ssa(pos=coord).to_table()
335+
336+
Without specifying the collection, the query returns results from multiple
337+
collections. For example this target has spectra from SOFIA as well as from
338+
Spitzer.
339+
340+
.. doctest-remote-data::
341+
342+
>>> from astropy.table import unique
343+
>>> unique(arp220_spectra, keys='dataid_collection')['dataid_collection']
344+
<MaskedColumn name='dataid_collection' dtype='object' description='IVOA Identifier of collection' length=4>
345+
goals
346+
sofia_fifils
347+
spitzer_irsenh
348+
spitzer_sha
349+
350+
351+
To list available collections for SSA queries, the
352+
`~astroquery.ipac.irsa.IrsaClass.list_collections` method is provided, and
353+
will return a `~astropy.table.Table`.
354+
355+
.. doctest-remote-data::
356+
357+
>>> from astroquery.ipac.irsa import Irsa
358+
>>> Irsa.list_collections(servicetype='SSA')
359+
<Table length=35>
360+
collection
361+
object
362+
------------------------
363+
champ
364+
goals
365+
herschel_gotcplus
366+
herschel_hexos
367+
herschel_hifistars
368+
herschel_hop
369+
herschel_hops
370+
herschel_magcloudscplus
371+
herschel_ppdisks
372+
herschel_prismas
373+
herschel_sag4
374+
herschel_shpdp
375+
herschel_v838mon
376+
herschel_vngs
377+
irtf_mearth
378+
irts
379+
iso_sws
380+
sofia_exes
381+
sofia_fifils
382+
sofia_flitecam
383+
sofia_forcast
384+
sofia_great
385+
spitzer_5muses
386+
spitzer_c2d
387+
spitzer_disks_sh_spectra
388+
spitzer_feps
389+
spitzer_irs_std
390+
spitzer_irsenh
391+
spitzer_m83m33
392+
spitzer_s5
393+
spitzer_sage
394+
spitzer_sass
395+
spitzer_sha
396+
spitzer_sings
397+
spitzer_ssgss
398+
399+
319400
Other Configurations
320401
--------------------
321402

0 commit comments

Comments
 (0)