File: custom_frames.rst

package info (click to toggle)
astropy 7.0.1-3
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 35,328 kB
  • sloc: python: 233,437; ansic: 55,264; javascript: 17,680; lex: 8,621; sh: 3,317; xml: 2,287; makefile: 191
file content (136 lines) | stat: -rw-r--r-- 4,242 bytes parent folder | download | duplicates (4)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
********************
Using a custom frame
********************

By default, `~astropy.visualization.wcsaxes.WCSAxes` will make use of a rectangular
frame for a plot, but this can be changed to provide any custom frame. The
following example shows how to use the built-in
:class:`~astropy.visualization.wcsaxes.frame.EllipticalFrame` class, which is an ellipse which extends to the same limits as the built-in rectangular frame:

.. plot::
   :context: reset
   :include-source:
   :align: center

    from astropy.wcs import WCS
    from astropy.io import fits
    from astropy.utils.data import get_pkg_data_filename
    from astropy.visualization.wcsaxes.frame import EllipticalFrame
    import matplotlib.pyplot as plt

    filename = get_pkg_data_filename('galactic_center/gc_msx_e.fits')
    hdu = fits.open(filename)[0]
    wcs = WCS(hdu.header)

    ax = plt.subplot(projection=wcs, frame_class=EllipticalFrame)

    ax.coords.grid(color='white')

    im = ax.imshow(hdu.data, vmin=-2.e-5, vmax=2.e-4, origin='lower')

    # Clip the image to the frame
    im.set_clip_path(ax.coords.frame.patch)

The :class:`~astropy.visualization.wcsaxes.frame.EllipticalFrame` class is especially useful for
all-sky plots such as Aitoff projections:

.. plot::
   :context: reset
   :include-source:
   :align: center

    from astropy.wcs import WCS
    from astropy.io import fits
    from astropy.utils.data import get_pkg_data_filename
    from astropy.visualization.wcsaxes.frame import EllipticalFrame
    from matplotlib import patheffects
    import matplotlib.pyplot as plt

    filename = get_pkg_data_filename('allsky/allsky_rosat.fits')
    hdu = fits.open(filename)[0]
    wcs = WCS(hdu.header)

    ax = plt.subplot(projection=wcs, frame_class=EllipticalFrame)

    path_effects=[patheffects.withStroke(linewidth=3, foreground='black')]
    ax.coords.grid(color='white')
    ax.coords['glon'].set_ticklabel(color='white', path_effects=path_effects)

    im = ax.imshow(hdu.data, vmin=0., vmax=300., origin='lower')

    # Clip the image to the frame
    im.set_clip_path(ax.coords.frame.patch)

However, you can also write your own frame class. The idea is to set up any
number of connecting spines that define the frame. You can define a frame as a
spine, but if you define it as multiple spines you will be able to control on
which spine the tick labels and ticks should appear.

The following example shows how you could for example define a hexagonal frame:

.. plot::
   :context: reset
   :include-source:
   :nofigs:

    import numpy as np
    from astropy.visualization.wcsaxes.frame import BaseFrame

    class HexagonalFrame(BaseFrame):

        spine_names = 'abcdef'

        def update_spines(self):

            xmin, xmax = self.parent_axes.get_xlim()
            ymin, ymax = self.parent_axes.get_ylim()

            ymid = 0.5 * (ymin + ymax)
            xmid1 = (xmin + xmax) / 4.
            xmid2 = (xmin + xmax) * 3. / 4.

            self['a'].data = np.array(([xmid1, ymin], [xmid2, ymin]))
            self['b'].data = np.array(([xmid2, ymin], [xmax, ymid]))
            self['c'].data = np.array(([xmax, ymid], [xmid2, ymax]))
            self['d'].data = np.array(([xmid2, ymax], [xmid1, ymax]))
            self['e'].data = np.array(([xmid1, ymax], [xmin, ymid]))
            self['f'].data = np.array(([xmin, ymid], [xmid1, ymin]))

which we can then use:

.. plot::
    :context:
    :include-source:
    :align: center

     from astropy.wcs import WCS
     from astropy.io import fits
     from astropy.utils.data import get_pkg_data_filename
     import matplotlib.pyplot as plt

     filename = get_pkg_data_filename('galactic_center/gc_msx_e.fits')
     hdu = fits.open(filename)[0]
     wcs = WCS(hdu.header)

     ax = plt.subplot(projection=wcs, frame_class=HexagonalFrame)

     ax.coords.grid(color='white')

     im = ax.imshow(hdu.data, vmin=-2.e-5, vmax=2.e-4, origin='lower')

     # Clip the image to the frame
     im.set_clip_path(ax.coords.frame.patch)


Frame properties
****************

The color and linewidth of the frame can also be set by

.. plot::
    :context:
    :include-source:
    :align: center

    ax.coords.frame.set_color('red')
    ax.coords.frame.set_linewidth(2)