Skip to content

Commit

Permalink
Fix LaTeX error due to labels that become empty in simplify_labels
Browse files Browse the repository at this point in the history
Fixes #8004.
  • Loading branch information
lpsinger committed Mar 19, 2020
1 parent 98c3ab2 commit 218d62e
Show file tree
Hide file tree
Showing 3 changed files with 44 additions and 1 deletion.
3 changes: 3 additions & 0 deletions CHANGES.rst
Original file line number Diff line number Diff line change
Expand Up @@ -496,6 +496,9 @@ astropy.visualization
maximum when there are less then ``min_npixels`` in the input array.
[#9913]

- Fix a bug in simplifying axis labels that affected non-rectangular frames.
[#8004, #9991]

astropy.wcs
^^^^^^^^^^^

Expand Down
39 changes: 38 additions & 1 deletion astropy/visualization/wcsaxes/tests/test_misc.py
Original file line number Diff line number Diff line change
Expand Up @@ -17,11 +17,14 @@
from astropy.tests.helper import catch_warnings

from astropy.visualization.wcsaxes.core import WCSAxes
from astropy.visualization.wcsaxes.frame import RectangularFrame, RectangularFrame1D
from astropy.visualization.wcsaxes.frame import (
EllipticalFrame, RectangularFrame, RectangularFrame1D)
from astropy.visualization.wcsaxes.utils import get_coord_meta
from astropy.visualization.wcsaxes.transforms import CurvedTransform

MATPLOTLIB_LT_21 = LooseVersion(matplotlib.__version__) < LooseVersion("2.1")
MATPLOTLIB_LT_22 = LooseVersion(matplotlib.__version__) < LooseVersion("2.2")
TEX_UNAVAILABLE = not matplotlib.checkdep_usetex(True)

DATA = os.path.abspath(os.path.join(os.path.dirname(__file__), 'data'))

Expand Down Expand Up @@ -435,3 +438,37 @@ def test_time_wcs(time_spectral_wcs_2d):
# with a time axis.

plt.subplot(projection=time_spectral_wcs_2d)


@pytest.mark.skipif('MATPLOTLIB_LT_22')
@pytest.mark.skipif('TEX_UNAVAILABLE')
def test_simplify_labels_usetex(ignore_matplotlibrc, tmpdir):
"""Regression test for https://github.com/astropy/astropy/issues/8004."""
plt.rc('text', usetex=True)

header = {
'NAXIS': 2,
'NAXIS1': 360,
'NAXIS2': 180,
'CRPIX1': 180.5,
'CRPIX2': 90.5,
'CRVAL1': 180.0,
'CRVAL2': 0.0,
'CDELT1': -2 * np.sqrt(2) / np.pi,
'CDELT2': 2 * np.sqrt(2) / np.pi,
'CTYPE1': 'RA---MOL',
'CTYPE2': 'DEC--MOL',
'RADESYS': 'ICRS'}

wcs = WCS(header)
fig, ax = plt.subplots(
subplot_kw=dict(frame_class=EllipticalFrame, projection=wcs))
ax.set_xlim(-0.5, header['NAXIS1'] - 0.5)
ax.set_ylim(-0.5, header['NAXIS2'] - 0.5)
ax.coords[0].set_ticklabel(exclude_overlapping=True)
ax.coords[1].set_ticklabel(exclude_overlapping=True)
ax.coords[0].set_ticks(spacing=45 * u.deg)
ax.coords[1].set_ticks(spacing=30 * u.deg)
ax.grid()

fig.savefig(tmpdir / 'plot.png')
3 changes: 3 additions & 0 deletions astropy/visualization/wcsaxes/ticklabels.py
Original file line number Diff line number Diff line change
Expand Up @@ -93,6 +93,9 @@ def simplify_labels(self):
self.text[axis][i] = self.text[axis][i][start:]
if starts_dollar:
self.text[axis][i] = '$' + self.text[axis][i]
# Remove any empty LaTeX inline math mode string
if self.text[axis][i] == '$$':
self.text[axis][i] = ''

def set_pad(self, value):
self._pad = value
Expand Down

0 comments on commit 218d62e

Please sign in to comment.