79def main(repo, collections, skymapName=None, tracts=None, visits=None, physicalFilters=None, bands=None,
80 ccds=None, ccdKey="detector", showPatch=False, saveFile=None, showCcds=False, visitVetoFile=None,
81 minOverlapFraction=None, trimToTracts=False, doUnscaledLimitRatio=False,
82 forceScaledLimitRatio=False):
83 if minOverlapFraction
is not None and tracts
is None:
84 raise RuntimeError(
"Must specify --tracts if --minOverlapFraction is set")
85 logger.info(
"Making butler for collections = %s in repo %s", collections, repo)
86 butler = dafButler.Butler(repo, collections=collections)
87 instrument = butler.find_dataset(
"camera").dataId[
"instrument"]
90 if skymapName
is None:
91 if instrument ==
"HSC":
92 skymapName =
"hsc_rings_v1"
93 detectorSkipList = [9]
94 elif instrument ==
"LSSTCam-imSim":
96 elif instrument ==
"LSSTComCamSim":
97 skymapName =
"ops_rehersal_prep_2k_v1"
98 elif instrument ==
"LATISS":
99 skymapName =
"latiss_v1"
100 elif instrument ==
"DECam":
101 skymapName =
"decam_rings_v1"
103 logger.error(
"Unknown skymapName for instrument: %s. Must specify --skymapName on command line.",
105 logger.info(
"instrument = %s skymapName = %s", instrument, skymapName)
106 camera = butler.get(
"camera", instrument=instrument)
107 skymap = butler.get(
"skyMap", instrument=instrument, skymap=skymapName)
110 if tracts
is not None:
114 if physicalFilters
is not None:
115 physicalFilterStr =
makeWhereInStr(
"physical_filter", physicalFilters, str)
116 whereStr +=
" AND " + physicalFilterStr
if len(whereStr)
else " " + physicalFilterStr
118 if bands
is not None:
120 whereStr +=
" AND " + bandStr
if len(whereStr)
else " " + bandStr
122 if len(whereStr) > 1:
123 whereStr =
"instrument=\'" + instrument +
"\' AND skymap=\'" + skymapName +
"\' AND " + whereStr
125 whereStr =
"instrument=\'" + instrument +
"\' AND skymap=\'" + skymapName +
"\'"
126 logger.info(
"Applying the following where clause in dataId search: %s", whereStr)
129 if visitVetoFile
is not None:
130 with open(visitVetoFile)
as f:
131 content = f.readlines()
132 visitVetoList = [
int(visit.strip())
for visit
in content]
135 dataRefs = list(butler.registry.queryDatasets(
"calexp", where=whereStr).expanded())
137 for dataRef
in dataRefs:
138 visit = dataRef.dataId.visit.id
139 if visit
not in visits
and visit
not in visitVetoList:
143 "List of visits (N=%d) satisfying where and veto clauses: %s", len(visits), visits
146 if len(visitVetoList) > 1:
147 visitListTemp = visits.copy()
148 for visit
in visitListTemp:
149 if visit
in visitVetoList:
151 logger.info(
"List of visits (N=%d) excluding veto list: %s}", len(visits), visits)
152 logger.info(
"List of visits (N=%d): %s", len(visits), visits)
157 if ((ccds
is None or ccdId
in ccds)
and ccd.getType() == cameraGeom.DetectorType.SCIENCE
158 and ccdId
not in detectorSkipList):
159 ccdIdList.append(ccdId)
161 nDetTot = len(ccdIdList)
163 visitIncludeList = []
167 if minOverlapFraction
is not None:
168 for i_v, visit
in enumerate(visits):
171 visitSummary = butler.get(
"visitSummary", visit=visit)
172 except LookupError
as e:
173 logger.warn(
"%s Will try to get wcs from calexp.", e)
175 if tracts
is not None:
177 tractInfo = skymap[tract]
178 sphCorners = tractInfo.wcs.pixelToSky(
geom.Box2D(tractInfo.bbox).getCorners())
179 tractConvexHull = sphgeom.ConvexPolygon.convexHull(
180 [coord.getVector()
for coord
in sphCorners])
181 for ccdId
in ccdIdList:
182 if ccdId
not in ccdOverlapList:
184 ccdKey, ccdId, visit, visitSummary=visitSummary, butler=butler,
186 if raCorners
is not None and decCorners
is not None:
188 for ra, dec
in zip(raCorners, decCorners):
191 detSphCorners.append(pt)
192 detConvexHull = sphgeom.ConvexPolygon.convexHull(
193 [coord.getVector()
for coord
in detSphCorners])
194 if tractConvexHull.contains(detConvexHull):
195 ccdOverlapList.append(ccdId)
197 if len(ccdOverlapList)/nDetTot >= minOverlapFraction:
199 if len(ccdOverlapList)/nDetTot < minOverlapFraction:
200 logger.info(
"Fraction of detectors overlaping any tract for visit %d (%.2f) < "
201 "minimum required (%.2f). Skipping visit...",
202 visit, len(ccdOverlapList)/nDetTot, minOverlapFraction)
204 if visit
not in visitIncludeList:
205 visitIncludeList.append(visit)
207 visitIncludeList = visits
212 cmap =
get_cmap(len(visitIncludeList))
214 maxVisitForLegend = 20
217 includedPhysicalFilters = []
218 for i_v, visit
in enumerate(visitIncludeList):
219 print(
"Working on visit %d [%d of %d]" % (visit, i_v + 1, len(visitIncludeList)), end=
"\r")
222 fillKwargs = {
"fill":
False,
"alpha": alphaEdge,
"facecolor":
None,
"edgecolor": color,
"lw": 0.6}
224 visitSummary = butler.get(
"visitSummary", visit=visit)
225 except Exception
as e:
226 logger.warn(
"%s Will try to get wcs from calexp.", e)
229 band, physicalFilter =
getBand(visitSummary=visitSummary, butler=butler, visit=visit)
230 if band
not in includedBands:
231 includedBands.append(band)
232 if physicalFilter
not in includedPhysicalFilters:
233 includedPhysicalFilters.append(physicalFilter)
235 for ccdId
in ccdIdList:
237 ccdKey, ccdId, visit, visitSummary=visitSummary, butler=butler)
238 if raCorners
is not None and decCorners
is not None:
241 if not inLegend
and len(visitIncludeList) <= maxVisitForLegend:
242 plt.fill(raCorners, decCorners, label=str(visit), **fillKwargs)
245 plt.fill(raCorners, decCorners, **fillKwargs)
246 plt.fill(raCorners, decCorners, fill=
True, alpha=alphaEdge/4, color=color,
248 if visit
not in finalVisitList:
249 finalVisitList.append(visit)
253 deltaRa =
max(raCorners) -
min(raCorners)
254 deltaDec =
max(decCorners) -
min(decCorners)
256 min(decCorners) + overlapFrac*deltaDec)
258 max(decCorners) - overlapFrac*deltaDec)
261 overlaps = [
not bboxDouble.overlaps(otherBbox)
for otherBbox
in bboxesPlotted]
264 str(ccdId), fontsize=6, ha=
"center", va=
"center", color=
"darkblue")
265 bboxesPlotted.append(bboxDouble)
267 logger.info(
"Final list of visits (N=%d) satisfying where and minOverlapFraction clauses: %s",
268 len(finalVisitList), finalVisitList)
270 raToDecLimitRatio =
None
272 tractList = list(set(skymap.findTractIdArray(ras, decs, degrees=
True)))
273 minVisitRa, maxVisitRa =
min(ras),
max(ras)
274 minVisitDec, maxVisitDec =
min(decs),
max(decs)
275 raVisitDiff = maxVisitRa - minVisitRa
276 decVisitDiff = maxVisitDec - minVisitDec
277 midVisitRa = minVisitRa + 0.5*raVisitDiff
278 midVisitDec = minVisitDec + 0.5*decVisitDiff
279 midRa = np.atleast_1d((midVisitRa*units.deg).to(units.radian).value).astype(np.float64)
280 midDec = np.atleast_1d((midVisitDec*units.deg).to(units.radian).value).astype(np.float64)
281 midSkyCoord = SkyCoord(midVisitRa*units.deg, midVisitDec*units.deg)
283 if tracts
is not None:
284 logger.info(
"No calexps were found, but --tracts list was provided, so will go ahead and "
285 "plot the empty tracts.")
289 raise RuntimeError(
"No data to plot (if you want to plot empty tracts, include them as "
290 "a blank-space separated list to the --tracts option.")
292 logger.info(
"List of tracts overlapping data: %s", tractList)
295 if forceScaledLimitRatio:
296 doUnscaledLimitRatio =
False
302 radiusMm = camera.computeMaxFocalPlaneRadius()
304 focalPlaneToFieldAngle = camera.getTransformMap().getTransform(
305 cameraGeom.FOCAL_PLANE, cameraGeom.FIELD_ANGLE
307 fpRadiusDeg = np.rad2deg(focalPlaneToFieldAngle.applyForward(fpRadiusPt))[0]
308 detectorRadiusDeg = fpRadiusDeg/np.sqrt(len(camera))
312 xDelta0 = xLimMax - xLimMin
313 yDelta0 = yLimMax - yLimMin
315 xDelta0 = raVisitDiff
316 yDelta0 = decVisitDiff
317 yLimMin = minVisitDec
318 yLimMax = maxVisitDec
319 raDecScaleThresh = 1.5
321 (xDelta0/yDelta0 > raDecScaleThresh
or yDelta0/xDelta0 > raDecScaleThresh)
322 and max(xDelta0, yDelta0) > 70*detectorRadiusDeg
323 and yLimMin < 75.0
and yLimMax > -75.0
326 "Sky coverage is large (and not too close to a pole), so not scaling to detector coords."
328 doUnscaledLimitRatio =
True
330 if not doUnscaledLimitRatio:
333 minDistToMidCoord = 1e12
336 for i_v, visit
in enumerate(visits):
338 visitSummary = butler.get(
"visitSummary", visit=visit)
339 except Exception
as e:
340 logger.warn(
"%s Will try to get wcs from calexp.", e)
342 for ccdId
in ccdIdList:
344 ccdKey, ccdId, visit, visitSummary=visitSummary, butler=butler, doLogWarn=
False)
345 if raCorners
is not None and decCorners
is not None:
347 for ra, dec
in zip(raCorners, decCorners):
350 detSphCorners.append(pt)
351 ptSkyCoord = SkyCoord(ra*units.deg, dec*units.deg)
352 separation = (midSkyCoord.separation(ptSkyCoord)).degree
353 if separation < minDistToMidCoord:
356 minDistToMidCoord = separation
357 detConvexHull = sphgeom.ConvexPolygon(
358 [coord.getVector()
for coord
in detSphCorners])
359 if detConvexHull.contains(midRa, midDec)
and raToDecLimitRatio
is None:
361 "visit/det overlapping plot coord mid point in RA/Dec: %d %d", visit, ccdId
363 raToDecLimitRatio = (
364 (
max(raCorners) -
min(raCorners))/(
max(decCorners) -
min(decCorners))
367 width = det.getBBox().getWidth()
368 height = det.getBBox().getHeight()
369 if raToDecLimitRatio > 1.0:
370 raToDecLimitRatio /=
max(height/width, width/height)
372 if raToDecLimitRatio < 1.0:
373 raToDecLimitRatio *=
max(height/width, width/height)
375 if raToDecLimitRatio
is not None:
378 if raToDecLimitRatio
is None and minSepVisit
is not None:
380 visitSummary = butler.get(
"visitSummary", visit=minSepVisit)
381 except Exception
as e:
382 logger.warn(
"%s Will try to get wcs from calexp.", e)
385 ccdKey, minSepCcdId, minSepVisit, visitSummary=visitSummary, butler=butler, doLogWarn=
False)
386 for ra, dec
in zip(raCorners, decCorners):
389 detSphCorners.append(pt)
390 detConvexHull = sphgeom.ConvexPolygon([coord.getVector()
for coord
in detSphCorners])
392 "visit/det closest to plot coord mid point in RA/Dec (none actually overlap it): %d %d",
393 minSepVisit, minSepCcdId
395 raToDecLimitRatio = (
max(raCorners) -
min(raCorners))/(
max(decCorners) -
min(decCorners))
396 det = camera[minSepCcdId]
397 width = det.getBBox().getWidth()
398 height = det.getBBox().getHeight()
399 if raToDecLimitRatio > 1.0:
400 raToDecLimitRatio /=
max(height/width, width/height)
402 if raToDecLimitRatio < 1.0:
403 raToDecLimitRatio *=
max(height/width, width/height)
405 if trimToTracts
is True:
406 xlim, ylim =
derivePlotLimits(tractLimitsDict, raToDecLimitRatio=raToDecLimitRatio, buffFrac=0.04)
408 visitLimitsDict = {
"allVisits": {
"ras": [minVisitRa, maxVisitRa],
"decs": [minVisitDec, maxVisitDec]}}
409 xlim, ylim =
derivePlotLimits(visitLimitsDict, raToDecLimitRatio=raToDecLimitRatio, buffFrac=0.04)
411 if doUnscaledLimitRatio:
412 boxAspectRatio = abs((ylim[1] - ylim[0])/(xlim[1] - xlim[0]))
420 if tracts
is not None:
421 tractOutlineList = list(set(tracts + tractList))
423 tractOutlineList = tractList
424 tractOutlineList.sort()
425 logger.info(
"List of tract outlines being plotted: %s", tractOutlineList)
426 for i_t, tract
in enumerate(tractOutlineList):
427 alpha =
max(0.1, alpha0 - i_t*1.0/len(tractOutlineList))
428 tractInfo = skymap[tract]
429 tCenter = tractInfo.ctr_coord
430 tCenterRa = tCenter.getRa().asDegrees()
431 tCenterDec = tCenter.getDec().asDegrees()
432 fracDeltaX = 0.02*abs((xlim[1] - xlim[0]))
433 fracDeltaY = 0.02*abs((ylim[1] - ylim[0]))
434 if (xlim[1] + fracDeltaX < tCenterRa < xlim[0] - fracDeltaX
435 and ylim[0] + fracDeltaY < tCenterDec < ylim[1] - fracDeltaY):
436 if len(tractOutlineList) > 1
or not showPatch:
438 plt.text(tCenterRa, tCenterDec, tract, fontsize=7, alpha=alpha, ha=
"center", va=
"center")
440 plt.text(tCenterRa, tCenterDec, tract, fontsize=7, alpha=1, color=
"white",
441 path_effects=[pathEffects.withStroke(linewidth=3, foreground=
"black")],
442 fontweight=500, ha=
"center", va=
"center", zorder=5)
443 ra, dec =
bboxToRaDec(tractInfo.bbox, tractInfo.getWcs())
444 plt.fill(ra, dec, fill=
False, edgecolor=
"k", lw=1, linestyle=
"dashed", alpha=alpha)
445 tractArtist = matplotlib.patches.Patch(fill=
False, edgecolor=
"k", linestyle=
"dashed", alpha=alpha)
446 tractHandleList.append(tractArtist)
447 tractStrList.append(str(tract))
450 for patch
in tractInfo:
451 ra, dec =
bboxToRaDec(patch.getInnerBBox(), tractInfo.getWcs())
452 plt.fill(ra, dec, fill=
False, edgecolor=patchColor, lw=0.5, linestyle=(0, (5, 6)),
457 str(patch.sequential_index), fontsize=5, color=patchColor,
458 ha=
"center", va=
"center", alpha=alpha)
464 ax.set_box_aspect(boxAspectRatio)
465 if abs(xlim[1] > 99.99):
466 ax.tick_params(
"x", labelrotation=45, pad=0, labelsize=8)
468 ax.tick_params(
"x", labelrotation=0, pad=0, labelsize=8)
469 ax.tick_params(
"y", labelsize=8)
470 ax.set_xlabel(
"RA (deg)", fontsize=9)
471 ax.set_ylabel(
"Dec (deg)", fontsize=9)
473 visitScaleOffset =
None
474 if len(visitIncludeList) > maxVisitForLegend:
475 nz = matplotlib.colors.Normalize()
476 colorBarScale = finalVisitList
477 if max(finalVisitList) > 9999999:
478 visitScaleOffset =
min(finalVisitList)
479 colorBarScale = [visit - visitScaleOffset
for visit
in finalVisitList]
480 nz.autoscale(colorBarScale)
481 cax, _ = matplotlib.colorbar.make_axes(plt.gca(), pad=0.03)
482 cax.tick_params(labelsize=7)
483 cb = matplotlib.colorbar.ColorbarBase(cax, cmap=cmap, norm=nz, alpha=alphaEdge,
484 format=
lambda x, _: f
"{x:.0f}")
485 cb.ax.yaxis.get_offset_text().set_fontsize(7)
486 colorBarLabel =
"visit number"
487 if visitScaleOffset
is not None:
488 colorBarLabel +=
" - {:d}".format(visitScaleOffset)
489 cb.set_label(colorBarLabel, rotation=-90, labelpad=13, fontsize=9)
490 tractLegend = Legend(ax, tractHandleList, tractStrList, loc=
"upper right", fancybox=
True,
491 shadow=
True, fontsize=5, title_fontsize=6, title=
"tracts")
492 ax.add_artist(tractLegend)
494 if len(visitIncludeList) > 0:
495 xBboxAnchor =
min(1.25,
max(1.03, boxAspectRatio*1.15))
496 ax.legend(loc=
"center left", bbox_to_anchor=(xBboxAnchor, 0.5), fancybox=
True,
497 shadow=
True, fontsize=6, title_fontsize=6, title=
"visits")
499 tractLegend = Legend(ax, tractHandleList, tractStrList, loc=
"center left", bbox_to_anchor=(1.0, 0.5),
500 fancybox=
True, shadow=
True, fontsize=6, title_fontsize=6, title=
"tracts")
501 ax.add_artist(tractLegend)
503 titleStr = repo +
"\n" + collections[0]
504 if len(collections) > 1:
505 for collection
in collections[1:]:
506 titleStr +=
"\n" + collection
507 titleStr +=
"\nnVisit: {}".format(str(len(finalVisitList)))
508 if minOverlapFraction
is not None:
509 titleStr +=
" (minOvlpFrac = {:.2f})".format(minOverlapFraction)
510 if len(includedBands) > 0:
511 titleStr +=
" bands: {}".format(str(includedBands).translate({ord(i):
None for i
in "[]\'"}))
512 if len(includedPhysicalFilters) > 0:
513 if len(includedPhysicalFilters[0]) > 9:
515 titleStr +=
" physical filters: {}".format(str(includedPhysicalFilters).translate(
516 {ord(i):
None for i
in "[]\'"}))
517 ax.set_title(
"{}".format(titleStr), fontsize=8)
520 if boxAspectRatio > 1.0:
521 minInches =
max(4.0, 0.3*abs(xlim[1] - xlim[0]))
523 yInches =
min(120.0, boxAspectRatio*minInches)
524 fig.set_size_inches(xInches, yInches)
525 if boxAspectRatio < 1.0:
526 minInches =
max(4.0, 0.3*abs(ylim[1] - ylim[0]))
527 xInches =
min(120.0, minInches/boxAspectRatio)
529 fig.set_size_inches(xInches, yInches)
530 if saveFile
is not None:
531 logger.info(
"Saving file in: %s", saveFile)
532 fig.savefig(saveFile, bbox_inches=
"tight", dpi=150)