LSST Applications 26.0.0,g0265f82a02+6660c170cc,g07994bdeae+30b05a742e,g0a0026dc87+17526d298f,g0a60f58ba1+17526d298f,g0e4bf8285c+96dd2c2ea9,g0ecae5effc+c266a536c8,g1e7d6db67d+6f7cb1f4bb,g26482f50c6+6346c0633c,g2bbee38e9b+6660c170cc,g2cc88a2952+0a4e78cd49,g3273194fdb+f6908454ef,g337abbeb29+6660c170cc,g337c41fc51+9a8f8f0815,g37c6e7c3d5+7bbafe9d37,g44018dc512+6660c170cc,g4a941329ef+4f7594a38e,g4c90b7bd52+5145c320d2,g58be5f913a+bea990ba40,g635b316a6c+8d6b3a3e56,g67924a670a+bfead8c487,g6ae5381d9b+81bc2a20b4,g93c4d6e787+26b17396bd,g98cecbdb62+ed2cb6d659,g98ffbb4407+81bc2a20b4,g9ddcbc5298+7f7571301f,ga1e77700b3+99e9273977,gae46bcf261+6660c170cc,gb2715bf1a1+17526d298f,gc86a011abf+17526d298f,gcf0d15dbbd+96dd2c2ea9,gdaeeff99f8+0d8dbea60f,gdb4ec4c597+6660c170cc,ge23793e450+96dd2c2ea9,gf041782ebf+171108ac67
LSST Data Management Base Package
Loading...
Searching...
No Matches
astrowidgets.py
Go to the documentation of this file.
1# This file is part of display_astrowidgets.
2#
3# Developed for the LSST Data Management System.
4# This product includes software developed by the LSST Project
5# (https://www.lsst.org).
6# See the COPYRIGHT file at the top-level directory of this distribution
7# for details of code ownership.
8#
9# This program is free software: you can redistribute it and/or modify
10# it under the terms of the GNU General Public License as published by
11# the Free Software Foundation, either version 3 of the License, or
12# (at your option) any later version.
13#
14# This program is distributed in the hope that it will be useful,
15# but WITHOUT ANY WARRANTY; without even the implied warranty of
16# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17# GNU General Public License for more details.
18#
19# You should have received a copy of the GNU General Public License
20# along with this program. If not, see <https://www.gnu.org/licenses/>.
21
22__all__ = ["AstroWidgetsVersion", "DisplayImpl"]
23
24import sys
25from astropy.table import Table
26
27import lsst.afw.display.interface as interface
28import lsst.afw.display.virtualDevice as virtualDevice
29import lsst.afw.display.ds9Regions as ds9Regions
30import lsst.afw.geom as afwGeom
31
32try:
33 from ginga.misc.log import get_logger
34 from ginga.AstroImage import AstroImage
35 from ginga.util.wcsmod.wcs_astropy import AstropyWCS
36 haveGinga = True
37except ImportError:
38 import logging
39 logging.getLogger("lsst.afw.display.astrowidgets").warning("Cannot import ginga libraries.")
40
42 def skyToPixel(*args, **kwargs):
43 pass
44
45 def pixelToSky(*args, **kwargs):
46 pass
47
48 haveGinga = False
49
50
51try:
52 import astrowidgets
53 haveAstrowidgets = True
54except ImportError:
55 haveAstrowidgets = False
56
57try:
58 _maskTransparency
59except NameError:
60 _maskTransparency = None
61
62
64 """Get the version of DS9 in use.
65
66 Returns
67 -------
68 version : `str`
69 Version of DS9 in use.
70 """
71 return astrowidgets.__version__
72
73
74class AstroWidgetsEvent(interface.Event):
75 """An event generated by a mouse or key click on the display"""
76
77 def __int__(self, k, x, y):
78 interface.Event.__init__(self, k, x, y)
79
80
81class DisplayImpl(virtualDevice.DisplayImpl):
82 """Virtual device display implementation.
83
84 Parameters
85 ----------
87 Display object to connect to.
88 dims : `tuple` [`int`, `int`], optional
89 Dimensions of the viewer window.
90 use_opencv : `bool`, optional
91 Should openCV be used to speed drawing?
92 verbose : `bool`, optional
93 Increase log verbosity?
94 """
95 markerDict = {'+': 'plus', 'x': 'cross', '.': 'circle', '*': 'circle', 'o': 'circle'}
96
97 def __init__(self, display, dims=None, use_opencv=False, verbose=False, *args, **kwargs):
98 virtualDevice.DisplayImpl.__init__(self, display, verbose)
99 if dims is None:
100 width, height = 1024, 768
101 else:
102 width, height = dims
103 if haveGinga:
104 self.logger = get_logger("ginga", log_stderr=True, level=40)
105 else:
106 self.logger = None
107 self._viewer = astrowidgets.ImageWidget(image_width=width, image_height=height,
108 use_opencv=use_opencv, logger=self.logger)
110 self._callbackDict = dict()
111
112 # We want to display the IW, but ginga has all the handles
113 self._gingaViewer = self._viewer._viewer
114
115 bd = self._gingaViewer.get_bindings()
116 bd.enable_all(True)
117 self._canvas = self._viewer.canvas
118 self._canvas.enable_draw(False)
120 self._redraw = True
121
122 def embed(self):
123 """Attach this display to the output of the current cell."""
124 return self._viewer
125
126 def get_viewer(self):
127 """Return the ginga viewer"""
128 return self._viewer
129
130 def show_color_bar(self, show=True):
131 """Show (or hide) the colour bar.
132
133 Parameters
134 ----------
135 show : `bool`, optional
136 Should the color bar be shown?
137 """
139
140 def show_pan_mark(self, show=True, color='red'):
141 """Show (or hide) the pan mark.
142
143 Parameters
144 ----------
145 show : `bool`, optional
146 Should the pan marker be shown?
147 color : `str`, optional
148 What color should the pan mark be?
149 """
150 self._gingaViewer.show_pan_mark(show, color)
151
152 def _setMaskTransparency(self, transparency, maskplane=None):
153 """Specify mask transparency (percent); or None to not set it when loading masks.
154
155 Parameters
156 ----------
157 transparency : `float`
158 Transparency of the masks in percent (0-100).
159 maskplane : `str`, optional
160 Unsupported option to only change the transparency of
161 certain masks.
162 """
163 if maskplane is not None:
164 print("display_astrowidgets is not yet able to set transparency for individual maskplanes" % maskplane, # noqa E501
165 file=sys.stderr)
166 return
167
168 self._maskTransparency = 0.01*transparency
169
170 def _getMaskTransparency(self, maskplane=None):
171 """Return the current mask transparency."""
172 return self._maskTransparency
173
174 def _mtv(self, image, mask=None, wcs=None, title=""):
175 """Display an Image and/or Mask on a ginga display
176
177 Parameters
178 ----------
180 Image to display.
181 mask : `lsst.afw.image.Mask`, optional
182 Mask to use, if the input does not contain one.
183 wcs : `ginga.util.wcsmod.wcs_astropy`
184 WCS to use, if the input does not contain one.
185 title : `str`, optional
186 Unsupported display title.
187 """
188 self._erase()
189 self._canvas.delete_all_objects()
190 self._buffer()
191 if haveGinga:
192 Aimage = AstroImage(inherit_primary_header=True)
193 Aimage.set_data(image.getArray())
194
195 self._gingaViewer.set_image(Aimage)
196
197 if wcs is not None:
198 if haveGinga:
199 _wcs = AstropyWCS(self.logger)
200 Aimage.lsst_wcs = WcsAdaptorForGinga(wcs)
201 _wcs.pixtoradec = Aimage.lsst_wcs.pixtoradec
202 _wcs.pixtosystem = Aimage.lsst_wcs.pixtosystem
203 _wcs.radectopix = Aimage.lsst_wcs.radectopix
204
205 Aimage.set_wcs(_wcs)
206 Aimage.wcs.wcs = Aimage.lsst_wcs
207
208 if mask:
209 maskColorFromName = {'BAD': 'red',
210 'SAT': 'green',
211 'INTRP': 'green',
212 'CR': 'magenta',
213 'EDGE': 'yellow',
214 'DETECTED': 'blue',
215 'DETECTED_NEGATIVE': 'cyan',
216 'SUSPECT': 'yellow',
217 'NO_DATA': 'orange',
218 'CROSSTALK': None,
219 'UNMASKEDNAN': None}
220 maskDict = dict()
221 for plane, bit in mask.getMaskPlaneDict().items():
222 color = maskColorFromName.get(plane, None)
223 if color:
224 maskDict[1 << bit] = color
225 # This value of 0.9 is pretty thick for the alpha.
226 self.overlay_mask(mask, maskDict,
228 self._buffer(enable=False)
229 self._flush()
230
231 def overlay_mask(self, maskImage, maskDict, maskAlpha):
232 """Draw mask onto the image display.
233
234 Parameters
235 ----------
236 maskImage : `lsst.afw.image.Mask`
237 Mask to display.
238 maskDict : `dict` [`str`, `str`]
239 Dictionary of mask plane names to colors.
240 maskAlpha : `float`
241 Transparency to display the mask.
242 """
243 import numpy as np
244 from ginga.RGBImage import RGBImage
245 from ginga import colors
246
247 maskArray = maskImage.getArray()
248 height, width = maskArray.shape
249 maskRGBA = np.zeros((height, width, 4), dtype=np.uint8)
250 nSet = np.zeros_like(maskArray, dtype=np.uint8)
251
252 for maskValue, maskColor in maskDict.items():
253 r, g, b = colors.lookup_color(maskColor)
254 isSet = (maskArray & maskValue) != 0
255 if (isSet == 0).all():
256 continue
257
258 maskRGBA[:, :, 0][isSet] = 255 * r
259 maskRGBA[:, :, 1][isSet] = 255 * g
260 maskRGBA[:, :, 2][isSet] = 255 * b
261
262 nSet[isSet] += 1
263
264 maskRGBA[:, :, 3][nSet == 0] = 0
265 maskRGBA[:, :, 3][nSet != 0] = 255 * maskAlpha
266
267 nSet[nSet == 0] = 1
268 for C in (0, 1, 2):
269 maskRGBA[:, :, C] //= nSet
270
271 rgb_img = RGBImage(data_np=maskRGBA)
272 Image = self._viewer.canvas.get_draw_class('image')
273 maskImageRGBA = Image(0, 0, rgb_img)
274
275 if "mask_overlay" in self._gingaViewer.canvas.get_tags():
276 self._gingaViewer.canvas.delete_object_by_tag("mask_overlay")
277 self._gingaViewer.canvas.add(maskImageRGBA, tag="mask_overlay")
278
279 def _buffer(self, enable=True):
280 self._redraw = not enable
281
282 def _flush(self):
283 self._gingaViewer.redraw(whence=3)
284
285 def _erase(self):
286 """Erase the display"""
287 self._canvas.delete_all_objects()
288
289 def _dot(self, symb, c, r, size, ctype, fontFamily="helvetica", textAngle=None, label='_dot'):
290 """Draw a symbol at (col,row) = (c,r) [0-based coordinates].
291
292 Parameters
293 ----------
294 symb : `str`
295 Symbol to draw. Should be one of '+', 'x', '*', 'o', '.'.
296 c : `int`
297 Image column for dot center (0-based coordinates).
298 r : `int`
299 Image row for dot center (0-based coordinate).
300 size : `int`
301 Size of dot.
302 fontFamily : `str`, optional
303 Font to use for text symbols.
304 textAngle : `float`, optional
305 Text rotation angle.
306 label : `str`, optional
307 Label to store this dot in the internal list.
308 """
309 dataTable = Table([{'x': c, 'y': r}])
310 if symb in '+x*.o':
311 self._viewer.marker = {'type': self.markerDict[symb], 'color': ctype, 'radius': size}
312 self._viewer.add_markers(dataTable, marker_name=label)
313 self._flush()
314 else:
315 Line = self._canvas.get_draw_class('line')
316 Text = self._canvas.get_draw_class('text')
317
318 for ds9Cmd in ds9Regions.dot(symb, c, r, size, fontFamily="helvetica", textAngle=None):
319 tmp = ds9Cmd.split('#')
320 cmd = tmp.pop(0).split()
321 comment = tmp.pop(0) if tmp else ""
322
323 cmd, args = cmd[0], cmd[1:]
324 if cmd == "line":
325 self._gingaViewer.canvas.add(Line(*[float(p) - 1 for p in args], color=ctype),
326 redraw=self._redraw)
327 elif cmd == "text":
328 x, y = [float(p) - 1 for p in args[0:2]]
329 self._gingaViewer.canvas.add(Text(x, y, symb, color=ctype), redraw=self._redraw)
330 else:
331 raise RuntimeError(ds9Cmd)
332 if comment:
333 print(comment)
334
335 def _drawLines(self, points, ctype):
336 """Connect the points, a list of (col,row).
337
338 Parameters
339 ----------
340 points : `list` [`tuple` [`int`, `int`]]
341 Points to connect with lines.
342 ctype : `str`
343 Color to use.
344 """
345 Line = self._gingaViewer.canvas.get_draw_class('line')
346 p0 = points[0]
347 for p in points[1:]:
348 self._gingaViewer.canvas.add(Line(p0[0], p0[1], p[0], p[1], color=ctype), redraw=self._redraw)
349 p0 = p
350
351 def beginMarking(self, symb='+', ctype='cyan', size=10, label='interactive'):
352 """Begin interactive mark adding.
353
354 Parameters
355 ----------
356 symb : `str`, optional
357 Symbol to use. Should be one of '+', 'x', '*', 'o', '.'.
358 ctype : `str`, optional
359 Color of markers.
360 size : `float`, optional
361 Size of marker.
362 label : `str`
363 Label to store this marker in the internal list.
364 """
365 self._viewer.start_marking(marker_name=label,
366 marker={'type': self.markerDict[symb], 'color': ctype, 'radius': size})
367
368 def endMarking(self):
369 """End interactive mark adding."""
370 self._viewer.stop_marking()
371
372 def getMarkers(self, label='interactive'):
373 """Get list of markers.
374
375 Parameters
376 ----------
377 label : `str`, optional
378 Marker label to return.
379
380 Returns
381 -------
382 table : `astropy.table.Table`
383 Table of markers with the given label.
384 """
385 return self._viewer.get_markers(marker_name=label)
386
387 def clearMarkers(self, label=None):
388 """Clear markers.
389
390 Parameters
391 ----------
392 label : `str`, optional
393 Marker label to clear. If None, all markers are cleared.
394 """
395 if label:
396 self._viewer.remove_markers(label)
397 else:
398 self._viewer.reset_markers()
399
400 def linkMarkers(self, ctype='brown', label='interactive'):
401 """Connect markers with lines.
402
403 Parameters
404 ----------
405 ctype : `str`, optional
406 Color to draw the lines.
407 label : `str`, optional
408 Marker label to connect. Lines are drawn in the order
409 found in the table.
410 """
411 Line = self._gingaViewer.canvas.get_draw_class('line')
412 table = self._viewer.get_markers(marker_name=label)
413
414 x0, y0 = (0, 0)
415 for rowCount, (x, y) in enumerate(table.iterrows('x', 'y')):
416 if rowCount != 0:
417 self._gingaViewer.canvas.add(Line(x0, y0, x, y, color=ctype), redraw=self._redraw)
418 x0 = x
419 y0 = y
420
421 def clearLines(self):
422 """Remove all lines from the display."""
423 self._gingaViewer.canvas.deleteObjects(list(self._gingaViewer.canvas.get_objects_by_kind('line')))
424
425 def _scale(self, algorithm, min, max, unit, *args, **kwargs):
426 """Set greyscale values.
427
428 Parameters
429 ----------
430 algorithm : `str`
431 Image scaling algorithm to use.
432 min : `float` or `str`
433 Minimum value to set to black. If a string, should be one of 'zscale' or 'minmax'.
434 max : `float`
435 Maximum value to set to white.
436 unit : `str`
437 Scaling units. This is ignored.
438 """
439 self._gingaViewer.set_color_map('gray')
440 self._gingaViewer.set_color_algorithm(algorithm)
441
442 if min == "zscale":
443 self._gingaViewer.set_autocut_params('zscale', contrast=0.25)
444 self._gingaViewer.auto_levels()
445 elif min == "minmax":
446 self._gingaViewer.set_autocut_params('minmax')
447 self._gingaViewer.auto_levels()
448 else:
449 if unit:
450 print("ginga: ignoring scale unit %s" % unit, file=sys.stderr)
451
452 self._gingaViewer.cut_levels(min, max)
453
454 def _show(self):
455 """Show the requested display.
456
457 In this case, embed it in the notebook (equivalent to
458 Display.get_viewer().show(); see also
459 Display.get_viewer().embed() N.b. These command *must* be the
460 last entry in their cell
461 """
462 return self._gingaViewer.show()
463
464 #
465 # Zoom and Pan
466 #
467 def _zoom(self, zoomfac):
468 """Zoom by specified amount
469
470 Parameters
471 ----------
472 zoomfac : `float`
473 Zoom factor to use.
474 """
475 self._gingaViewer.scale_to(zoomfac, zoomfac)
476
477 def _pan(self, colc, rowc):
478 """Pan to (colc, rowc)
479
480 Parameters
481 ----------
482 colc : `int`
483 Column to center in viewer (0-based coordinate).
484 rowc : `int`
485 Row to center in viewer (0-based coordinate).
486 """
487 self._gingaViewer.set_pan(colc, rowc)
488
489 def _getEvent(self):
490 """Listen for a key press on a frame in DS9 and return an event.
491
492 Returns
493 -------
494 event : `Ds9Event`
495 Event with (key, x, y).
496 """
497 pass
498
499
500# Copy ginga's WCS implementation
502 """A class to adapt the LSST Wcs class for Ginga.
503
504 This was taken largely from the afw.display.ginga package.
505
506 Parameters
507 ----------
508 wcs : `ginga.util.wcsmod.wcs_astropy`
509 WCS to adapt for Ginga.
510 """
511 def __init__(self, wcs):
512 self._wcs = wcs
513
514 def pixtoradec(self, idxs, coords='data'):
515 """Return (ra, dec) in degrees given a position in pixels.
516
517 Parameters
518 ----------
519 idxs : `list` [`tuple` [`float`, `float`]]
520 Pixel locations to convert.
521 coords : `str`, optional
522 This parameter is ignored.
523 Returns
524 -------
525 ra : `list`
526 RA position in degrees.
527 dec : `list`
528 DEC position in degrees.
529 """
530 ra, dec = self._wcs.pixelToSky(*idxs)
531
532 return ra.asDegrees(), dec.asDegrees()
533
534 def pixtosystem(self, idxs, system=None, coords='data'):
535 """Return (ra, dec) in degrees given a position in pixels.
536
537 Parameters
538 ----------
539 idxs : `list` [`tuple` [`float`, `float`]]
540 Pixel locations to convert.
541 system : `str`, optional
542 This parameter is ignored.
543 coords : `str`, optional
544 This parameter is ignored.
545
546 Returns
547 -------
548 ra : `list`
549 RA position in degrees.
550 dec : `list`
551 DEC position in degrees.
552 """
553 return self.pixtoradec(idxs, coords=coords)
554
555 def radectopix(self, ra_deg, dec_deg, coords='data', naxispath=None):
556 """Return (x, y) in pixels given (ra, dec) in degrees
557
558 Parameters
559 ----------
560 ra_deg : `list` [`float`]
561 RA position in degrees.
562 dec_deg : `list` [`float`]
563 DEC position in degrees.
564 coords : `str`, optional
565 This parameter is ignored.
566 naxispath : `str`, optional
567 This parameter is ignored.
568
569 Returns
570 -------
571 out : `tuple` [`list` [`float, `float`]]
572 Image coordates for input positions.
573 """
574 return self._wcs.skyToPixel(ra_deg*afwGeom.degrees, dec_deg*afwGeom.degrees)
575
576 def all_pix2world(self, *args, **kwargs):
577 out = []
578 print(f"{args}")
579 for pos in args[0]:
580 r, d = self.pixtoradec(pos)
581 out.append([r, d])
582 return tuple(out)
583
584 def datapt_to_wcspt(self, *args):
585 return (0.0, 0.0)
586
587 def wcspt_to_datapt(self, *args):
588 return (0.0, 0.0)
std::vector< SchemaItem< Flag > > * items
A class to contain the data, WCS, and other information needed to describe an image of the sky.
Definition Exposure.h:72
A class to represent a 2-dimensional array of pixels.
Definition Image.h:51
Represent a 2-dimensional array of bitmask pixels.
Definition Mask.h:77
beginMarking(self, symb='+', ctype='cyan', size=10, label='interactive')
overlay_mask(self, maskImage, maskDict, maskAlpha)
_setMaskTransparency(self, transparency, maskplane=None)
__init__(self, display, dims=None, use_opencv=False, verbose=False, *args, **kwargs)
_scale(self, algorithm, min, max, unit, *args, **kwargs)
_mtv(self, image, mask=None, wcs=None, title="")
linkMarkers(self, ctype='brown', label='interactive')
_dot(self, symb, c, r, size, ctype, fontFamily="helvetica", textAngle=None, label='_dot')
radectopix(self, ra_deg, dec_deg, coords='data', naxispath=None)
pixtosystem(self, idxs, system=None, coords='data')
daf::base::PropertyList * list
Definition fits.cc:928