Skip to content

Commit 128628d

Browse files
committed
Merge pull request #446 from loicseguin/simbad-radius
Raise exception when radius=0 passed to Simbad query_region
2 parents 356438a + ac3dc3e commit 128628d

File tree

3 files changed

+44
-5
lines changed

3 files changed

+44
-5
lines changed

astroquery/simbad/core.py

Lines changed: 6 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -682,7 +682,7 @@ def _args_to_payload(self, *args, **kwargs):
682682
ra, dec, frame = _parse_coordinates(coordinates)
683683
args = [ra, dec]
684684
kwargs['frame'] = frame
685-
if kwargs.get('radius'):
685+
if kwargs.get('radius') is not None:
686686
kwargs['radius'] = _parse_radius(kwargs['radius'])
687687
# rename equinox to equi as required by SIMBAD script
688688
if kwargs.get('equinox'):
@@ -785,7 +785,11 @@ def _parse_radius(radius):
785785
try:
786786
angle = commons.parse_radius(radius)
787787
# find the most appropriate unit - d, m or s
788-
index = min([i for (i, val) in enumerate(angle.dms) if int(val) > 0])
788+
nonzero_indices = [i for (i, val) in enumerate(angle.dms) if int(val) > 0]
789+
if len(nonzero_indices) > 0:
790+
index = min(nonzero_indices)
791+
else:
792+
index = 2 # use arcseconds when radius smaller than 1 arcsecond
789793
unit = ('d', 'm', 's')[index]
790794
if unit == 'd':
791795
return str(int(angle.degree)) + unit

astroquery/simbad/tests/test_simbad.py

Lines changed: 24 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -251,6 +251,30 @@ def test_query_region(patch_post, coordinates, radius, equinox, epoch):
251251
assert isinstance(result2, Table)
252252

253253

254+
@pytest.mark.parametrize(('coordinates', 'radius', 'equinox', 'epoch'),
255+
[(ICRS_COORDS, 0, None, None)])
256+
def test_query_region_radius_error(patch_post, coordinates, radius, equinox, epoch):
257+
with pytest.raises(u.UnitsError):
258+
result1 = simbad.core.Simbad.query_region(coordinates, radius=radius,
259+
equinox=equinox, epoch=epoch)
260+
with pytest.raises(u.UnitsError):
261+
result2 = simbad.core.Simbad().query_region(coordinates, radius=radius,
262+
equinox=equinox, epoch=epoch)
263+
264+
265+
@pytest.mark.parametrize(('coordinates', 'radius', 'equinox', 'epoch'),
266+
[(ICRS_COORDS, "0d", None, None),
267+
(GALACTIC_COORDS, 1.0 * u.marcsec, 2000.0, 'J2000')
268+
])
269+
def test_query_region_small_radius(patch_post, coordinates, radius, equinox, epoch):
270+
result1 = simbad.core.Simbad.query_region(coordinates, radius=radius,
271+
equinox=equinox, epoch=epoch)
272+
result2 = simbad.core.Simbad().query_region(coordinates, radius=radius,
273+
equinox=equinox, epoch=epoch)
274+
assert isinstance(result1, Table)
275+
assert isinstance(result2, Table)
276+
277+
254278
@pytest.mark.parametrize(('object_name', 'wildcard'),
255279
[("m1", None),
256280
("m [0-9]", True)

astroquery/simbad/tests/test_simbad_remote.py

Lines changed: 14 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -122,8 +122,19 @@ def test_query_objects_null(self):
122122

123123
# Special case of null test: zero-sized region
124124
def test_query_region_null(self):
125-
# This test will fail if some object is discovered within 1 arcsec of the origin
126-
# (Due to other bugs I could not set a smaller radius, or radius=0)
127-
result = simbad.core.Simbad.query_region(coord.ICRS("00h00m0.0s 00h00m0.0s"), radius=1 * u.arcsec,
125+
result = simbad.core.Simbad.query_region(coord.ICRS("00h00m0.0s 00h00m0.0s"), radius="0d",
128126
equinox=2000.0, epoch='J2000')
129127
assert result is None
128+
129+
# Special case of null test: very small region
130+
def test_query_small_region_null(self):
131+
result = simbad.core.Simbad.query_region(coord.ICRS("00h00m0.0s 00h00m0.0s"), radius=1.0 * u.marcsec,
132+
equinox=2000.0, epoch='J2000')
133+
assert result is None
134+
135+
# Special case : zero-sized region with one object
136+
def test_query_zero_sized_region(self):
137+
result = simbad.core.Simbad.query_region(coord.ICRS("20h54m05.6889s 37d01m17.380s"), radius="0d",
138+
equinox=2000.0, epoch='J2000')
139+
# This should find a single star, BD+36 4308
140+
assert len(result) == 1

0 commit comments

Comments
 (0)