LSST Applications 27.0.0,g0265f82a02+469cd937ee,g02d81e74bb+21ad69e7e1,g1470d8bcf6+cbe83ee85a,g2079a07aa2+e67c6346a6,g212a7c68fe+04a9158687,g2305ad1205+94392ce272,g295015adf3+81dd352a9d,g2bbee38e9b+469cd937ee,g337abbeb29+469cd937ee,g3939d97d7f+72a9f7b576,g487adcacf7+71499e7cba,g50ff169b8f+5929b3527e,g52b1c1532d+a6fc98d2e7,g591dd9f2cf+df404f777f,g5a732f18d5+be83d3ecdb,g64a986408d+21ad69e7e1,g858d7b2824+21ad69e7e1,g8a8a8dda67+a6fc98d2e7,g99cad8db69+f62e5b0af5,g9ddcbc5298+d4bad12328,ga1e77700b3+9c366c4306,ga8c6da7877+71e4819109,gb0e22166c9+25ba2f69a1,gb6a65358fc+469cd937ee,gbb8dafda3b+69d3c0e320,gc07e1c2157+a98bf949bb,gc120e1dc64+615ec43309,gc28159a63d+469cd937ee,gcf0d15dbbd+72a9f7b576,gdaeeff99f8+a38ce5ea23,ge6526c86ff+3a7c1ac5f1,ge79ae78c31+469cd937ee,gee10cc3b42+a6fc98d2e7,gf1cff7945b+21ad69e7e1,gfbcc870c63+9a11dc8c8f
LSST Data Management Base Package
Loading...
Searching...
No Matches
Functions | Variables
showVisitSkyMap Namespace Reference

Functions

 bboxToRaDec (bbox, wcs)
 
 getValueAtPercentile (values, percentile=0.5)
 
 get_cmap (n, name="hsv")
 
 main (repo, collections, skymapName=None, tracts=None, visits=None, physicalFilters=None, bands=None, ccds=None, ccdKey="detector", showPatch=False, saveFile=None, showCcds=False, visitVetoFile=None, minOverlapFraction=None, trimToTracts=False)
 
 makeWhereInStr (parameterName, parameterList, parameterType)
 
 getTractLimitsDict (skymap, tractList)
 
 derivePlotLimits (limitsDict, raToDecLimitRatio=1.0, buffFrac=0.0)
 
 setLimitsToEqualRatio (xMin, xMax, yMin, yMax)
 
 getDetRaDecCorners (ccdKey, ccdId, visit, visitSummary=None, butler=None, doLogWarn=True)
 
 getBand (visitSummary=None, butler=None, visit=None)
 

Variables

 logger = lsstLog.Log.getLogger("showVisitSkyMap.py")
 
 parser = argparse.ArgumentParser()
 
 type
 
 help
 
 str
 
 nargs
 
 metavar
 
 required
 
 default
 
 None
 
 int
 
 action
 
 float
 
 args = parser.parse_args()
 
 collections
 
 skymapName
 
 tracts
 
 visits
 
 physicalFilters
 
 bands
 
 ccds
 
 ccdKey
 
 showPatch
 
 saveFile
 
 showCcds
 
 visitVetoFile
 
 minOverlapFraction
 
 trimToTracts
 

Function Documentation

◆ bboxToRaDec()

showVisitSkyMap.bboxToRaDec ( bbox,
wcs )
Get the corners of a BBox and convert them to lists of RA and Dec.

Definition at line 42 of file showVisitSkyMap.py.

42def bboxToRaDec(bbox, wcs):
43 """Get the corners of a BBox and convert them to lists of RA and Dec."""
44 sphPoints = wcs.pixelToSky(geom.Box2D(bbox).getCorners())
45 ra = [float(sph.getRa().asDegrees()) for sph in sphPoints]
46 dec = [float(sph.getDec().asDegrees()) for sph in sphPoints]
47 return ra, dec
48
49
A floating-point coordinate rectangle geometry.
Definition Box.h:413

◆ derivePlotLimits()

showVisitSkyMap.derivePlotLimits ( limitsDict,
raToDecLimitRatio = 1.0,
buffFrac = 0.0 )
Derive the axis limits to encompass all points in limitsDict.

Parameters
----------
limitsDict : `dict` [`dict`]
    A dictionary keyed on any id. Each entry includes a `dict`
    keyed on "ras" and "decs" including (at least the minimum
    and maximum) RA and Dec values in units of degrees.
raToDecLimitRatio : `float`, optional
    The aspect ratio between RA and Dec to set the plot limits to.  This
    is to namely to set this ratio to that of the focal plane (i.e. such
    that a square detector appears as a square), but any aspect ratio can,
    in principle, be requested.

Returns
-------
xlim, ylim : `tuple` [`float`]
    Two tuples containing the derived min and max values for the x and
    y-axis limits (in degrees), respectively.

Definition at line 521 of file showVisitSkyMap.py.

521def derivePlotLimits(limitsDict, raToDecLimitRatio=1.0, buffFrac=0.0):
522 """Derive the axis limits to encompass all points in limitsDict.
523
524 Parameters
525 ----------
526 limitsDict : `dict` [`dict`]
527 A dictionary keyed on any id. Each entry includes a `dict`
528 keyed on "ras" and "decs" including (at least the minimum
529 and maximum) RA and Dec values in units of degrees.
530 raToDecLimitRatio : `float`, optional
531 The aspect ratio between RA and Dec to set the plot limits to. This
532 is to namely to set this ratio to that of the focal plane (i.e. such
533 that a square detector appears as a square), but any aspect ratio can,
534 in principle, be requested.
535
536 Returns
537 -------
538 xlim, ylim : `tuple` [`float`]
539 Two tuples containing the derived min and max values for the x and
540 y-axis limits (in degrees), respectively.
541 """
542 xLimMin, yLimMin = 1e12, 1e12
543 xLimMax, yLimMax = -1e12, -1e12
544 for limitId, limits in limitsDict.items():
545 xLimMin = min(xLimMin, min(limits["ras"]))
546 xLimMax = max(xLimMax, max(limits["ras"]))
547 yLimMin = min(yLimMin, min(limits["decs"]))
548 yLimMax = max(yLimMax, max(limits["decs"]))
549 xDelta0 = xLimMax - xLimMin
550 yDelta0 = yLimMax - yLimMin
551
552 if raToDecLimitRatio == 1.0:
553 if xDelta0 > yDelta0:
554 xLimMin -= buffFrac*yDelta0
555 xLimMax += buffFrac*yDelta0
556 else:
557 yLimMin -= buffFrac*yDelta0
558 yLimMax += buffFrac*yDelta0
559 xLimMin, xLimMax, yLimMin, yLimMax = setLimitsToEqualRatio(xLimMin, xLimMax, yLimMin, yLimMax)
560 else:
561 xLimMin -= buffFrac*xDelta0
562 xLimMax += buffFrac*xDelta0
563 yLimMin -= buffFrac*yDelta0
564 yLimMax += buffFrac*yDelta0
565 xLimMin, xLimMax, yLimMin, yLimMax = setLimitsToEqualRatio(xLimMin, xLimMax, yLimMin, yLimMax)
566 xDelta = xLimMax - xLimMin
567 yDelta = yLimMax - yLimMin
568 if raToDecLimitRatio > 1.0:
569 if yDelta0 > xDelta:
570 xMid = xLimMin + 0.5*(xDelta)
571 xLimMin = xMid - 0.5*yDelta*raToDecLimitRatio
572 xLimMax = xMid + 0.5*yDelta*raToDecLimitRatio
573 else:
574 yMid = yLimMin + 0.5*(yDelta)
575 yLimMin = yMid - 0.5*xDelta/raToDecLimitRatio
576 yLimMax = yMid + 0.5*xDelta/raToDecLimitRatio
577 else:
578 if xDelta0 > yDelta0:
579 yMid = yLimMin + 0.5*(yDelta)
580 yLimMin = yMid - 0.5*xDelta/raToDecLimitRatio
581 yLimMax = yMid + 0.5*xDelta/raToDecLimitRatio
582 else:
583 xMid = xLimMin + 0.5*(xDelta)
584 xLimMin = xMid - 0.5*yDelta*raToDecLimitRatio
585 xLimMax = xMid + 0.5*yDelta*raToDecLimitRatio
586 xlim = xLimMax, xLimMin
587 ylim = yLimMin, yLimMax
588 return xlim, ylim
589
590
int min
int max

◆ get_cmap()

showVisitSkyMap.get_cmap ( n,
name = "hsv" )
Returns a function that maps each index in 0, 1, ..., n-1 to a distinct
RGB color; the keyword argument name must be a standard mpl colormap name.

Definition at line 72 of file showVisitSkyMap.py.

72def get_cmap(n, name="hsv"):
73 """Returns a function that maps each index in 0, 1, ..., n-1 to a distinct
74 RGB color; the keyword argument name must be a standard mpl colormap name.
75 """
76 return matplotlib.colormaps[name].resampled(n)
77
78

◆ getBand()

showVisitSkyMap.getBand ( visitSummary = None,
butler = None,
visit = None )
Determine band and physical filter for given visit.

Parameters
----------
visitSummary : `lsst.afw.table.ExposureCatalog` or `None`, optional
    The visitSummary table for the visit for which to determine the band.
butler : `lsst.daf.butler.Butler` or `None`, optional
    The butler from which to look up the Dimension Records. Only needed
    if ``visitSummary`` is `None`.
visit : `int` or `None, optional
    The visit number for which to determine the band. Only needed
    if ``visitSummary`` is `None`.

Returns
-------
band, physicalFilter : `str`
    The band and physical filter for the given visit.

Definition at line 645 of file showVisitSkyMap.py.

645def getBand(visitSummary=None, butler=None, visit=None):
646 """Determine band and physical filter for given visit.
647
648 Parameters
649 ----------
650 visitSummary : `lsst.afw.table.ExposureCatalog` or `None`, optional
651 The visitSummary table for the visit for which to determine the band.
652 butler : `lsst.daf.butler.Butler` or `None`, optional
653 The butler from which to look up the Dimension Records. Only needed
654 if ``visitSummary`` is `None`.
655 visit : `int` or `None, optional
656 The visit number for which to determine the band. Only needed
657 if ``visitSummary`` is `None`.
658
659 Returns
660 -------
661 band, physicalFilter : `str`
662 The band and physical filter for the given visit.
663 """
664 if visitSummary is not None:
665 band = visitSummary[0]["band"]
666 physicalFilter = visitSummary[0]["physical_filter"]
667 else:
668 record = list(butler.registry.queryDimensionRecords("band", visit=visit))[0]
669 band = record.name
670 record = list(butler.registry.queryDimensionRecords("physical_filter", visit=visit))[0]
671 physicalFilter = record.name
672 return band, physicalFilter
673
674

◆ getDetRaDecCorners()

showVisitSkyMap.getDetRaDecCorners ( ccdKey,
ccdId,
visit,
visitSummary = None,
butler = None,
doLogWarn = True )
Compute the RA/Dec corners lists for a given detector in a visit.

Definition at line 620 of file showVisitSkyMap.py.

620def getDetRaDecCorners(ccdKey, ccdId, visit, visitSummary=None, butler=None, doLogWarn=True):
621 """Compute the RA/Dec corners lists for a given detector in a visit.
622 """
623 raCorners, decCorners = None, None
624 if visitSummary is not None:
625 row = visitSummary.find(ccdId)
626 if row is None:
627 if doLogWarn:
628 logger.warn("No row found for %d in visitSummary of visit %d. "
629 "Skipping and continuing...", ccdId, visit)
630 else:
631 raCorners = list(row["raCorners"])
632 decCorners = list(row["decCorners"])
633 else:
634 try:
635 dataId = {"visit": visit, ccdKey: ccdId}
636 wcs = butler.get("calexp.wcs", dataId)
637 bbox = butler.get("calexp.bbox", dataId)
638 raCorners, decCorners = bboxToRaDec(bbox, wcs)
639 except LookupError as e:
640 logger.warn("%s Skipping and continuing...", e)
641
642 return raCorners, decCorners
643
644

◆ getTractLimitsDict()

showVisitSkyMap.getTractLimitsDict ( skymap,
tractList )
Return a dict containing tract limits needed for outline plotting.

Parameters
----------
skymap : `lsst.skymap.BaseSkyMap`
    The sky map used for this dataset. Used to obtain tract
    parameters.
tractList : `list` [`int`]
    The list of tract ids (as integers) for which to determine the
    limits.

Returns
-------
tractLimitsDict : `dict` [`dict`]
    A dictionary keyed on tract id. Each entry includes a `dict`
    including the tract RA corners, Dec corners, and the tract center,
    all in units of degrees. These are used for plotting the tract
    outlines.

Definition at line 483 of file showVisitSkyMap.py.

483def getTractLimitsDict(skymap, tractList):
484 """Return a dict containing tract limits needed for outline plotting.
485
486 Parameters
487 ----------
488 skymap : `lsst.skymap.BaseSkyMap`
489 The sky map used for this dataset. Used to obtain tract
490 parameters.
491 tractList : `list` [`int`]
492 The list of tract ids (as integers) for which to determine the
493 limits.
494
495 Returns
496 -------
497 tractLimitsDict : `dict` [`dict`]
498 A dictionary keyed on tract id. Each entry includes a `dict`
499 including the tract RA corners, Dec corners, and the tract center,
500 all in units of degrees. These are used for plotting the tract
501 outlines.
502 """
503 tractLimitsDict = {}
504 for tract in tractList:
505 tractInfo = skymap[tract]
506 tractBbox = tractInfo.outer_sky_polygon.getBoundingBox()
507 tractCenter = tractBbox.getCenter()
508 tractRa0 = (tractCenter[0] - tractBbox.getWidth() / 2).asDegrees()
509 tractRa1 = (tractCenter[0] + tractBbox.getWidth() / 2).asDegrees()
510 tractDec0 = (tractCenter[1] - tractBbox.getHeight() / 2).asDegrees()
511 tractDec1 = (tractCenter[1] + tractBbox.getHeight() / 2).asDegrees()
512 tractLimitsDict[tract] = {
513 "ras": [tractRa0, tractRa1, tractRa1, tractRa0, tractRa0],
514 "decs": [tractDec0, tractDec0, tractDec1, tractDec1, tractDec0],
515 "center": [tractCenter[0].asDegrees(), tractCenter[1].asDegrees()],
516 }
517
518 return tractLimitsDict
519
520

◆ getValueAtPercentile()

showVisitSkyMap.getValueAtPercentile ( values,
percentile = 0.5 )
Return a value a fraction of the way between the min and max values in a
list.

Parameters
----------
values : `list` [`float`]
    The list of values under consideration.
percentile : `float`, optional
    The percentile (expressed as a number between 0.0 and 1.0) at which
    to determine the value in `values`.

Returns
-------
result : `float`
    The value at the given percentile of ``values``.

Definition at line 50 of file showVisitSkyMap.py.

50def getValueAtPercentile(values, percentile=0.5):
51 """Return a value a fraction of the way between the min and max values in a
52 list.
53
54 Parameters
55 ----------
56 values : `list` [`float`]
57 The list of values under consideration.
58 percentile : `float`, optional
59 The percentile (expressed as a number between 0.0 and 1.0) at which
60 to determine the value in `values`.
61
62 Returns
63 -------
64 result : `float`
65 The value at the given percentile of ``values``.
66 """
67 m = min(values)
68 interval = max(values) - m
69 return m + percentile*interval
70
71

◆ main()

showVisitSkyMap.main ( repo,
collections,
skymapName = None,
tracts = None,
visits = None,
physicalFilters = None,
bands = None,
ccds = None,
ccdKey = "detector",
showPatch = False,
saveFile = None,
showCcds = False,
visitVetoFile = None,
minOverlapFraction = None,
trimToTracts = False )

Definition at line 79 of file showVisitSkyMap.py.

81 minOverlapFraction=None, trimToTracts=False):
82 if minOverlapFraction is not None and tracts is None:
83 raise RuntimeError("Must specify --tracts if --minOverlapFraction is set")
84 logger.info("Making butler for collections = {} in repo {}".format(collections, repo))
85 butler = dafButler.Butler(repo, collections=collections)
86 instrument = butler.find_dataset("camera").dataId["instrument"]
87 detectorSkipList = []
88 # Make a guess at the skymapName if not provided
89 if skymapName is None:
90 if instrument == "HSC":
91 skymapName = "hsc_rings_v1"
92 detectorSkipList = [9] # detector 9 has long been dead for HSC
93 elif instrument == "LSSTCam-imSim":
94 skymapName = "DC2"
95 elif instrument == "LATISS":
96 skymapName = "latiss_v1"
97 elif instrument == "DECam":
98 skymapName = "decam_rings_v1"
99 else:
100 logger.error("Unknown skymapName for instrument: %s. Must specify --skymapName on command line.",
101 instrument)
102 logger.info("instrument = {} skymapName = {}".format(instrument, skymapName))
103 camera = butler.get("camera", instrument=instrument)
104 skymap = butler.get("skyMap", instrument=instrument, skymap=skymapName)
105
106 whereStr = ""
107 if tracts is not None:
108 tractStr = makeWhereInStr("tract", tracts, int)
109 whereStr += tractStr
110
111 if physicalFilters is not None:
112 physicalFilterStr = makeWhereInStr("physical_filter", physicalFilters, str)
113 whereStr += " AND " + physicalFilterStr if len(whereStr) else " " + physicalFilterStr
114
115 if bands is not None:
116 bandStr = makeWhereInStr("band", bands, str)
117 whereStr += " AND " + bandStr if len(whereStr) else " " + bandStr
118
119 if len(whereStr) > 1:
120 whereStr = "instrument=\'" + instrument + "\' AND skymap=\'" + skymapName + "\' AND " + whereStr
121 else:
122 whereStr = "instrument=\'" + instrument + "\' AND skymap=\'" + skymapName + "\'"
123 logger.info("Applying the following where clause in dataId search: {} ".format(whereStr))
124
125 visitVetoList = []
126 if visitVetoFile is not None:
127 with open(visitVetoFile) as f:
128 content = f.readlines()
129 visitVetoList = [int(visit.strip()) for visit in content]
130
131 if visits is None:
132 dataRefs = list(butler.registry.queryDatasets("calexp", where=whereStr).expanded())
133 visits = []
134 for dataRef in dataRefs:
135 visit = dataRef.dataId.visit.id
136 if visit not in visits and visit not in visitVetoList:
137 visits.append(visit)
138 visits.sort()
139 logger.info("List of visits (N={}) satisfying where and veto clauses: {}".format(len(visits),
140 visits))
141 else:
142 if len(visitVetoList) > 1:
143 visitListTemp = visits.copy()
144 for visit in visitListTemp:
145 if visit in visitVetoList:
146 visits.remove(visit)
147 logger.info("List of visits (N={}) excluding veto list: {}".format(len(visits), visits))
148
149 ccdIdList = []
150 for ccd in camera:
151 ccdId = ccd.getId()
152 if ((ccds is None or ccdId in ccds) and ccd.getType() == cameraGeom.DetectorType.SCIENCE
153 and ccdId not in detectorSkipList):
154 ccdIdList.append(ccdId)
155 ccdIdList.sort()
156 nDetTot = len(ccdIdList)
157
158 visitIncludeList = []
159 # Determine the fraction of detectors that overlap any tract under
160 # consideration. If this fraction does not exceed minOverlapFraction,
161 # skip this visit.
162 if minOverlapFraction is not None:
163 for i_v, visit in enumerate(visits):
164 ccdOverlapList = []
165 try:
166 visitSummary = butler.get("visitSummary", visit=visit)
167 except LookupError as e:
168 logger.warn("%s Will try to get wcs from calexp.", e)
169 visitSummary = None
170 if tracts is not None:
171 for tract in tracts:
172 tractInfo = skymap[tract]
173 sphCorners = tractInfo.wcs.pixelToSky(geom.Box2D(tractInfo.bbox).getCorners())
174 tractConvexHull = sphgeom.ConvexPolygon.convexHull(
175 [coord.getVector() for coord in sphCorners])
176 for ccdId in ccdIdList:
177 if ccdId not in ccdOverlapList:
178 raCorners, decCorners = getDetRaDecCorners(
179 ccdKey, ccdId, visit, visitSummary=visitSummary, butler=butler,
180 doLogWarn=False)
181 if raCorners is not None and decCorners is not None:
182 detSphCorners = []
183 for ra, dec in zip(raCorners, decCorners):
184 pt = geom.SpherePoint(geom.Angle(ra, geom.degrees),
185 geom.Angle(dec, geom.degrees))
186 detSphCorners.append(pt)
187 detConvexHull = sphgeom.ConvexPolygon.convexHull(
188 [coord.getVector() for coord in detSphCorners])
189 if tractConvexHull.contains(detConvexHull):
190 ccdOverlapList.append(ccdId)
191
192 if len(ccdOverlapList)/nDetTot >= minOverlapFraction:
193 break
194 if len(ccdOverlapList)/nDetTot < minOverlapFraction:
195 logger.info("Fraction of detectors overlaping any tract for visit %d (%.2f) < "
196 "minimum required (%.2f). Skipping visit...",
197 visit, len(ccdOverlapList)/nDetTot, minOverlapFraction)
198 else:
199 if visit not in visitIncludeList:
200 visitIncludeList.append(visit)
201 else:
202 visitIncludeList = visits
203
204 # draw the CCDs
205 ras, decs = [], []
206 bboxesPlotted = []
207 cmap = get_cmap(len(visitIncludeList))
208 alphaEdge = 0.7
209 maxVisitForLegend = 20
210 finalVisitList = []
211 includedBands = []
212 includedPhysicalFilters = []
213 for i_v, visit in enumerate(visitIncludeList):
214 print("Working on visit %d [%d of %d]" % (visit, i_v + 1, len(visitIncludeList)), end="\r")
215 inLegend = False
216 color = cmap(i_v)
217 fillKwargs = {"fill": False, "alpha": alphaEdge, "facecolor": None, "edgecolor": color, "lw": 0.6}
218 try:
219 visitSummary = butler.get("visitSummary", visit=visit)
220 except LookupError as e:
221 logger.warn("%s Will try to get wcs from calexp.", e)
222 visitSummary = None
223
224 band, physicalFilter = getBand(visitSummary=visitSummary, butler=butler, visit=visit)
225 if band not in includedBands:
226 includedBands.append(band)
227 if physicalFilter not in includedPhysicalFilters:
228 includedPhysicalFilters.append(physicalFilter)
229
230 for ccdId in ccdIdList:
231 raCorners, decCorners = getDetRaDecCorners(
232 ccdKey, ccdId, visit, visitSummary=visitSummary, butler=butler)
233 if raCorners is not None and decCorners is not None:
234 ras += raCorners
235 decs += decCorners
236 if not inLegend and len(visitIncludeList) <= maxVisitForLegend:
237 plt.fill(raCorners, decCorners, label=str(visit), **fillKwargs)
238 inLegend = True
239 else:
240 plt.fill(raCorners, decCorners, **fillKwargs)
241 plt.fill(raCorners, decCorners, fill=True, alpha=alphaEdge/4, color=color,
242 edgecolor=color)
243 if visit not in finalVisitList:
244 finalVisitList.append(visit)
245 # add CCD serial numbers
246 if showCcds:
247 overlapFrac = 0.2
248 deltaRa = max(raCorners) - min(raCorners)
249 deltaDec = max(decCorners) - min(decCorners)
250 minPoint = geom.Point2D(min(raCorners) + overlapFrac*deltaRa,
251 min(decCorners) + overlapFrac*deltaDec)
252 maxPoint = geom.Point2D(max(raCorners) - overlapFrac*deltaRa,
253 max(decCorners) - overlapFrac*deltaDec)
254 # Use doubles in Box2D to check overlap
255 bboxDouble = geom.Box2D(minPoint, maxPoint)
256 overlaps = [not bboxDouble.overlaps(otherBbox) for otherBbox in bboxesPlotted]
257 if all(overlaps):
258 plt.text(getValueAtPercentile(raCorners), getValueAtPercentile(decCorners),
259 str(ccdId), fontsize=6, ha="center", va="center", color="darkblue")
260 bboxesPlotted.append(bboxDouble)
261
262 logger.info("Final list of visits (N={}) satisfying where and minOverlapFraction clauses: {}"
263 .format(len(finalVisitList), finalVisitList))
264
265 raToDecLimitRatio = None
266 if len(ras) > 0:
267 tractList = list(set(skymap.findTractIdArray(ras, decs, degrees=True)))
268 minVisitRa, maxVisitRa = min(ras), max(ras)
269 minVisitDec, maxVisitDec = min(decs), max(decs)
270 raVisitDiff = maxVisitRa - minVisitRa
271 decVisitDiff = maxVisitDec - minVisitDec
272 midVisitRa = minVisitRa + 0.5*raVisitDiff
273 midVisitDec = minVisitDec + 0.5*decVisitDiff
274 midRa = np.atleast_1d((midVisitRa*units.deg).to(units.radian).value).astype(np.float64)
275 midDec = np.atleast_1d((midVisitDec*units.deg).to(units.radian).value).astype(np.float64)
276 midSkyCoord = SkyCoord(midVisitRa*units.deg, midVisitDec*units.deg)
277 else:
278 if tracts is not None:
279 logger.info("No calexps were found, but --tracts list was provided, so will go ahead and "
280 "plot the empty tracts.")
281 tractList = tracts
282 trimToTracts = True
283 raToDecLimitRatio = 1.0
284 else:
285 raise RuntimeError("No data to plot (if you want to plot empty tracts, include them as "
286 "a blank-space separated list to the --tracts option.")
287 tractList.sort()
288 logger.info("List of tracts overlapping data: {}".format(tractList))
289 tractLimitsDict = getTractLimitsDict(skymap, tractList)
290
291 # Find a detector that contains the mid point in RA/Dec (or the closest
292 # one) to set the plot aspect ratio.
293 minDistToMidCood = 1e12
294 minSepVisit = None
295 minSepCcdId = None
296 for i_v, visit in enumerate(visits):
297 try:
298 visitSummary = butler.get("visitSummary", visit=visit)
299 except LookupError as e:
300 logger.warn("%s Will try to get wcs from calexp.", e)
301 visitSummary = None
302 for ccdId in ccdIdList:
303 raCorners, decCorners = getDetRaDecCorners(
304 ccdKey, ccdId, visit, visitSummary=visitSummary, butler=butler, doLogWarn=False)
305 if raCorners is not None and decCorners is not None:
306 detSphCorners = []
307 for ra, dec in zip(raCorners, decCorners):
308 pt = geom.SpherePoint(geom.Angle(ra, geom.degrees),
309 geom.Angle(dec, geom.degrees))
310 detSphCorners.append(pt)
311 ptSkyCoord = SkyCoord(ra*units.deg, dec*units.deg)
312 separation = (midSkyCoord.separation(ptSkyCoord)).degree
313 if separation < minDistToMidCood:
314 minSepVisit = visit
315 minSepCcdId = ccdId
316 minDistToMidCood = separation
317 detConvexHull = sphgeom.ConvexPolygon(
318 [coord.getVector() for coord in detSphCorners])
319 if detConvexHull.contains(midRa, midDec):
320 logger.info("visit/det overlapping plot coord mid point in RA/Dec: %d %d", visit, ccdId)
321 raToDecLimitRatio = (max(raCorners) - min(raCorners))/(max(decCorners) - min(decCorners))
322 det = camera[ccdId]
323 width = det.getBBox().getWidth()
324 height = det.getBBox().getHeight()
325 if raToDecLimitRatio > 1.0:
326 raToDecLimitRatio /= max(height/width, width/height)
327 else:
328 if raToDecLimitRatio < 1.0:
329 raToDecLimitRatio *= max(height/width, width/height)
330 break
331 if raToDecLimitRatio is not None:
332 break
333
334 if raToDecLimitRatio is None:
335 try:
336 visitSummary = butler.get("visitSummary", visit=minSepVisit)
337 except LookupError as e:
338 logger.warn("%s Will try to get wcs from calexp.", e)
339 visitSummary = None
340 raCorners, decCorners = getDetRaDecCorners(
341 ccdKey, minSepCcdId, minSepVisit, visitSummary=visitSummary, butler=butler, doLogWarn=False)
342 for ra, dec in zip(raCorners, decCorners):
343 pt = geom.SpherePoint(geom.Angle(ra, geom.degrees),
344 geom.Angle(dec, geom.degrees))
345 detSphCorners.append(pt)
346 detConvexHull = sphgeom.ConvexPolygon([coord.getVector() for coord in detSphCorners])
347 logger.info("visit/det closest to plot coord mid point in RA/Dec (none actually overlat it): %d %d",
348 minSepVisit, minSepCcdId)
349 raToDecLimitRatio = (max(raCorners) - min(raCorners))/(max(decCorners) - min(decCorners))
350 det = camera[minSepCcdId]
351 width = det.getBBox().getWidth()
352 height = det.getBBox().getHeight()
353 if raToDecLimitRatio > 1.0:
354 raToDecLimitRatio /= max(height/width, width/height)
355 else:
356 if raToDecLimitRatio < 1.0:
357 raToDecLimitRatio *= max(height/width, width/height)
358
359 if trimToTracts is True:
360 xlim, ylim = derivePlotLimits(tractLimitsDict, raToDecLimitRatio=raToDecLimitRatio, buffFrac=0.04)
361 else:
362 visitLimitsDict = {"allVisits": {"ras": [minVisitRa, maxVisitRa], "decs": [minVisitDec, maxVisitDec]}}
363 xlim, ylim = derivePlotLimits(visitLimitsDict, raToDecLimitRatio=raToDecLimitRatio, buffFrac=0.04)
364
365 # draw the skymap
366 alpha0 = 1.0
367 tractHandleList = []
368 tractStrList = []
369 if tracts is not None:
370 tractOutlineList = list(set(tracts + tractList))
371 else:
372 tractOutlineList = tractList
373 tractOutlineList.sort()
374 logger.info("List of tract outlines being plotted: {}".format(tractOutlineList))
375 for i_t, tract in enumerate(tractOutlineList):
376 alpha = max(0.1, alpha0 - i_t*1.0/len(tractOutlineList))
377 tractInfo = skymap[tract]
378 tCenter = tractInfo.ctr_coord
379 tCenterRa = tCenter.getRa().asDegrees()
380 tCenterDec = tCenter.getDec().asDegrees()
381 fracDeltaX = 0.02*abs((xlim[1] - xlim[0]))
382 fracDeltaY = 0.02*abs((ylim[1] - ylim[0]))
383 if (xlim[1] + fracDeltaX < tCenterRa < xlim[0] - fracDeltaX
384 and ylim[0] + fracDeltaY < tCenterDec < ylim[1] - fracDeltaY):
385 if len(tractOutlineList) > 1 or not showPatch:
386 plt.text(tCenterRa, tCenterDec, tract, fontsize=7, alpha=alpha, ha="center", va="center")
387 ra, dec = bboxToRaDec(tractInfo.bbox, tractInfo.getWcs())
388 plt.fill(ra, dec, fill=False, edgecolor="k", lw=1, linestyle="dashed", alpha=alpha)
389 tractArtist = mpatches.Patch(fill=False, edgecolor="k", linestyle="dashed", alpha=alpha)
390 tractHandleList.append(tractArtist)
391 tractStrList.append(str(tract))
392 if showPatch:
393 patchColor = "k"
394 for patch in tractInfo:
395 ra, dec = bboxToRaDec(patch.getInnerBBox(), tractInfo.getWcs())
396 plt.fill(ra, dec, fill=False, edgecolor=patchColor, lw=0.5, linestyle=(0, (5, 6)),
397 alpha=alpha)
398 if (xlim[1] + fracDeltaX < getValueAtPercentile(ra) < xlim[0] - fracDeltaX
399 and ylim[0] + fracDeltaY < getValueAtPercentile(dec) < ylim[1] - fracDeltaY):
400 plt.text(getValueAtPercentile(ra), getValueAtPercentile(dec),
401 str(patch.sequential_index), fontsize=5, color=patchColor,
402 ha="center", va="center", alpha=alpha)
403
404 # add labels and save
405 ax = plt.gca()
406 ax.set_xlim(xlim)
407 ax.set_ylim(ylim)
408 ax.set_box_aspect(1)
409 if abs(xlim[1] > 99.99):
410 ax.tick_params("x", labelrotation=45, pad=0, labelsize=8)
411 else:
412 ax.tick_params("x", labelrotation=0, pad=0, labelsize=8)
413 ax.tick_params("y", labelsize=8)
414 ax.set_xlabel("RA (deg)", fontsize=9)
415 ax.set_ylabel("Dec (deg)", fontsize=9)
416
417 visitScaleOffset = None
418 if len(visitIncludeList) > maxVisitForLegend:
419 nz = matplotlib.colors.Normalize()
420 colorBarScale = finalVisitList
421 if max(finalVisitList) > 9999999:
422 visitScaleOffset = min(finalVisitList)
423 colorBarScale = [visit - visitScaleOffset for visit in finalVisitList]
424 nz.autoscale(colorBarScale)
425 cax, _ = matplotlib.colorbar.make_axes(plt.gca(), pad=0.03)
426 cax.tick_params(labelsize=7)
427 cb = matplotlib.colorbar.ColorbarBase(cax, cmap=cmap, norm=nz, alpha=alphaEdge,
428 format=lambda x, _: f"{x:.0f}")
429 cb.ax.yaxis.get_offset_text().set_fontsize(7)
430 colorBarLabel = "visit number"
431 if visitScaleOffset is not None:
432 colorBarLabel += " - {:d}".format(visitScaleOffset)
433 cb.set_label(colorBarLabel, rotation=-90, labelpad=13, fontsize=9)
434 tractLegend = Legend(ax, tractHandleList, tractStrList, loc="upper right",
435 fancybox=True, shadow=True, fontsize=5, title_fontsize=6, title="tracts")
436 ax.add_artist(tractLegend)
437 else:
438 if len(visitIncludeList) > 0:
439 ax.legend(loc="center left", bbox_to_anchor=(1.15, 0.5), fancybox=True, shadow=True,
440 fontsize=6, title_fontsize=6, title="visits")
441 # Create the second legend and add the artist manually.
442 tractLegend = Legend(ax, tractHandleList, tractStrList, loc="center left", bbox_to_anchor=(1.0, 0.5),
443 fancybox=True, shadow=True, fontsize=6, title_fontsize=6, title="tracts")
444 ax.add_artist(tractLegend)
445
446 titleStr = repo + "\n" + collections[0]
447 if len(collections) > 1:
448 for collection in collections[1:]:
449 titleStr += "\n" + collection
450 titleStr += "\nnVisit: {}".format(str(len(finalVisitList)))
451 if minOverlapFraction is not None:
452 titleStr += " (minOvlpFrac = {:.2f})".format(minOverlapFraction)
453 if len(includedBands) > 0:
454 titleStr += " bands: {}".format(str(includedBands).translate({ord(i): None for i in "[]\'"}))
455 if len(includedPhysicalFilters) > 0:
456 if len(includedPhysicalFilters[0]) > 9:
457 titleStr += "\n"
458 titleStr += " physical filters: {}".format(str(includedPhysicalFilters).translate(
459 {ord(i): None for i in "[]\'"}))
460 ax.set_title("{}".format(titleStr), fontsize=8)
461
462 fig = plt.gcf()
463 if saveFile is not None:
464 logger.info("Saving file in: %s", saveFile)
465 fig.savefig(saveFile, bbox_inches="tight", dpi=150)
466 else:
467 fig.show()
468
469
A class representing an angle.
Definition Angle.h:128
Point in an unspecified spherical coordinate system.
Definition SpherePoint.h:57
daf::base::PropertySet * set
Definition fits.cc:931

◆ makeWhereInStr()

showVisitSkyMap.makeWhereInStr ( parameterName,
parameterList,
parameterType )
Create the string to be used in the where clause for registry lookup.

Definition at line 470 of file showVisitSkyMap.py.

470def makeWhereInStr(parameterName, parameterList, parameterType):
471 """Create the string to be used in the where clause for registry lookup.
472 """
473 typeStr = "\'" if parameterType is str else ""
474 whereInStr = parameterName + " IN (" + typeStr + str(parameterList[0]) + typeStr
475 if len(parameterList) > 1:
476 for param in parameterList[1:]:
477 whereInStr += ", " + typeStr + str(param) + typeStr
478 whereInStr += ")"
479
480 return whereInStr
481
482

◆ setLimitsToEqualRatio()

showVisitSkyMap.setLimitsToEqualRatio ( xMin,
xMax,
yMin,
yMax )
For a given set of x/y min/max, redefine to have equal aspect ratio.

The limits are extended on both ends such that the central value is
preserved.

Parameters
----------
xMin, xMax, yMin, yMax : `float`
    The min/max values of the x/y ranges for which to match in dynamic
    range while perserving the central values.

Returns
-------
xMin, xMax, yMin, yMax : `float`
    The adjusted min/max values of the x/y ranges with equal aspect ratios.

Definition at line 591 of file showVisitSkyMap.py.

591def setLimitsToEqualRatio(xMin, xMax, yMin, yMax):
592 """For a given set of x/y min/max, redefine to have equal aspect ratio.
593
594 The limits are extended on both ends such that the central value is
595 preserved.
596
597 Parameters
598 ----------
599 xMin, xMax, yMin, yMax : `float`
600 The min/max values of the x/y ranges for which to match in dynamic
601 range while perserving the central values.
602
603 Returns
604 -------
605 xMin, xMax, yMin, yMax : `float`
606 The adjusted min/max values of the x/y ranges with equal aspect ratios.
607 """
608 xDelta = xMax - xMin
609 yDelta = yMax - yMin
610 deltaDiff = yDelta - xDelta
611 if deltaDiff > 0:
612 xMin -= 0.5 * deltaDiff
613 xMax += 0.5 * deltaDiff
614 elif deltaDiff < 0:
615 yMin -= 0.5 * abs(deltaDiff)
616 yMax += 0.5 * abs(deltaDiff)
617 return xMin, xMax, yMin, yMax
618
619

Variable Documentation

◆ action

showVisitSkyMap.action

Definition at line 697 of file showVisitSkyMap.py.

◆ args

showVisitSkyMap.args = parser.parse_args()

Definition at line 710 of file showVisitSkyMap.py.

◆ bands

showVisitSkyMap.bands

Definition at line 712 of file showVisitSkyMap.py.

◆ ccdKey

showVisitSkyMap.ccdKey

Definition at line 712 of file showVisitSkyMap.py.

◆ ccds

showVisitSkyMap.ccds

Definition at line 712 of file showVisitSkyMap.py.

◆ collections

showVisitSkyMap.collections

Definition at line 711 of file showVisitSkyMap.py.

◆ default

showVisitSkyMap.default

Definition at line 682 of file showVisitSkyMap.py.

◆ float

showVisitSkyMap.float

Definition at line 706 of file showVisitSkyMap.py.

◆ help

showVisitSkyMap.help

Definition at line 678 of file showVisitSkyMap.py.

◆ int

showVisitSkyMap.int

Definition at line 683 of file showVisitSkyMap.py.

◆ logger

showVisitSkyMap.logger = lsstLog.Log.getLogger("showVisitSkyMap.py")

Definition at line 39 of file showVisitSkyMap.py.

◆ metavar

showVisitSkyMap.metavar

Definition at line 681 of file showVisitSkyMap.py.

◆ minOverlapFraction

showVisitSkyMap.minOverlapFraction

Definition at line 714 of file showVisitSkyMap.py.

◆ nargs

showVisitSkyMap.nargs

Definition at line 679 of file showVisitSkyMap.py.

◆ None

showVisitSkyMap.None

Definition at line 682 of file showVisitSkyMap.py.

◆ parser

showVisitSkyMap.parser = argparse.ArgumentParser()

Definition at line 676 of file showVisitSkyMap.py.

◆ physicalFilters

showVisitSkyMap.physicalFilters

Definition at line 712 of file showVisitSkyMap.py.

◆ required

showVisitSkyMap.required

Definition at line 681 of file showVisitSkyMap.py.

◆ saveFile

showVisitSkyMap.saveFile

Definition at line 713 of file showVisitSkyMap.py.

◆ showCcds

showVisitSkyMap.showCcds

Definition at line 713 of file showVisitSkyMap.py.

◆ showPatch

showVisitSkyMap.showPatch

Definition at line 713 of file showVisitSkyMap.py.

◆ skymapName

showVisitSkyMap.skymapName

Definition at line 711 of file showVisitSkyMap.py.

◆ str

showVisitSkyMap.str

Definition at line 679 of file showVisitSkyMap.py.

◆ tracts

showVisitSkyMap.tracts

Definition at line 711 of file showVisitSkyMap.py.

◆ trimToTracts

showVisitSkyMap.trimToTracts

Definition at line 715 of file showVisitSkyMap.py.

◆ type

showVisitSkyMap.type

Definition at line 677 of file showVisitSkyMap.py.

◆ visits

showVisitSkyMap.visits

Definition at line 711 of file showVisitSkyMap.py.

◆ visitVetoFile

showVisitSkyMap.visitVetoFile

Definition at line 714 of file showVisitSkyMap.py.