LSST Applications 29.0.1,g0fba68d861+132dd21e0a,g107a963962+1bb9f809a9,g1fd858c14a+005be21cae,g21d47ad084+8a07b29876,g325378336f+5d73323c8f,g330003fc43+40b4eaffc6,g35bb328faa+fcb1d3bbc8,g36ff55ed5b+9c28a42a87,g4e0f332c67+5fbd1e3e73,g53246c7159+fcb1d3bbc8,g60b5630c4e+9c28a42a87,g67b6fd64d1+a38b34ea13,g78460c75b0+2f9a1b4bcd,g786e29fd12+cf7ec2a62a,g7b71ed6315+fcb1d3bbc8,g86c591e316+6b2b2d0295,g8852436030+bf14db0e33,g89139ef638+a38b34ea13,g8b8da53e10+e3777245af,g9125e01d80+fcb1d3bbc8,g989de1cb63+a38b34ea13,g9f1445be69+9c28a42a87,g9f33ca652e+52c8f07962,ga9baa6287d+9c28a42a87,ga9e4eb89a6+9f84bd6575,gabe3b4be73+1e0a283bba,gb037a4e798+f3cbcd26c0,gb1101e3267+e7be8da0f8,gb58c049af0+f03b321e39,gb89ab40317+a38b34ea13,gcf25f946ba+bf14db0e33,gd6cbbdb0b4+bce7f7457e,gd9a9a58781+fcb1d3bbc8,gde0f65d7ad+53d424b1ae,ge278dab8ac+222406d50a,ge410e46f29+a38b34ea13,ge80e9994a3+664d6357dc,gf67bdafdda+a38b34ea13
LSST Data Management Base Package
Loading...
Searching...
No Matches
getTemplate.py
Go to the documentation of this file.
1# This file is part of ip_diffim.
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/>.
21import collections
22
23import numpy as np
24from deprecated.sphinx import deprecated
25
26import lsst.afw.image as afwImage
27import lsst.geom as geom
28import lsst.afw.geom as afwGeom
29import lsst.afw.table as afwTable
30import lsst.afw.math as afwMath
31import lsst.pex.config as pexConfig
32import lsst.pipe.base as pipeBase
33from lsst.skymap import BaseSkyMap
34from lsst.ip.diffim.dcrModel import DcrModel
35from lsst.meas.algorithms import CoaddPsf, CoaddPsfConfig
36from lsst.utils.timer import timeMethod
37
38__all__ = ["GetTemplateTask", "GetTemplateConfig",
39 "GetDcrTemplateTask", "GetDcrTemplateConfig"]
40
41
42class GetTemplateConnections(pipeBase.PipelineTaskConnections,
43 dimensions=("instrument", "visit", "detector", "skymap"),
44 defaultTemplates={"coaddName": "goodSeeing",
45 "warpTypeSuffix": "",
46 "fakesType": ""}):
47 bbox = pipeBase.connectionTypes.Input(
48 doc="Bounding box of exposure to determine the geometry of the output template.",
49 name="{fakesType}calexp.bbox",
50 storageClass="Box2I",
51 dimensions=("instrument", "visit", "detector"),
52 )
53 wcs = pipeBase.connectionTypes.Input(
54 doc="WCS of the exposure that we will construct the template for.",
55 name="{fakesType}calexp.wcs",
56 storageClass="Wcs",
57 dimensions=("instrument", "visit", "detector"),
58 )
59 skyMap = pipeBase.connectionTypes.Input(
60 doc="Geometry of the tracts and patches that the coadds are defined on.",
61 name=BaseSkyMap.SKYMAP_DATASET_TYPE_NAME,
62 dimensions=("skymap", ),
63 storageClass="SkyMap",
64 )
65 coaddExposures = pipeBase.connectionTypes.Input(
66 doc="Coadds that may overlap the desired region, as possible inputs to the template."
67 " Will be restricted to those that directly overlap the projected bounding box.",
68 dimensions=("tract", "patch", "skymap", "band"),
69 storageClass="ExposureF",
70 name="{fakesType}{coaddName}Coadd{warpTypeSuffix}",
71 multiple=True,
72 deferLoad=True,
73 deferGraphConstraint=True
74 )
75
76 template = pipeBase.connectionTypes.Output(
77 doc="Warped template, pixel matched to the bounding box and WCS.",
78 dimensions=("instrument", "visit", "detector"),
79 storageClass="ExposureF",
80 name="{fakesType}{coaddName}Diff_templateExp{warpTypeSuffix}",
81 )
82
83
84class GetTemplateConfig(pipeBase.PipelineTaskConfig,
85 pipelineConnections=GetTemplateConnections):
86 templateBorderSize = pexConfig.Field(
87 dtype=int,
88 default=20,
89 doc="Number of pixels to grow the requested template image to account for warping"
90 )
91 warp = pexConfig.ConfigField(
92 dtype=afwMath.Warper.ConfigClass,
93 doc="warper configuration",
94 )
95 coaddPsf = pexConfig.ConfigField(
96 doc="Configuration for CoaddPsf",
97 dtype=CoaddPsfConfig,
98 )
99
100 def setDefaults(self):
101 # Use a smaller cache: per SeparableKernel.computeCache, this should
102 # give a warping error of a fraction of a count (these must match).
103 self.warp.cacheSize = 100000
104 self.coaddPsf.cacheSize = self.warp.cacheSize
105 # The WCS for LSST should be smoothly varying, so we can use a longer
106 # interpolation length for WCS evaluations.
107 self.warp.interpLength = 100
108 self.warp.warpingKernelName = 'lanczos3'
109 self.coaddPsf.warpingKernelName = self.warp.warpingKernelName
110
111
112class GetTemplateTask(pipeBase.PipelineTask):
113 ConfigClass = GetTemplateConfig
114 _DefaultName = "getTemplate"
115
116 def __init__(self, *args, **kwargs):
117 super().__init__(*args, **kwargs)
118 self.warper = afwMath.Warper.fromConfig(self.config.warp)
119 self.schema = afwTable.ExposureTable.makeMinimalSchema()
120 self.schema.addField('tract', type=np.int32, doc='Which tract this exposure came from.')
121 self.schema.addField('patch', type=np.int32, doc='Which patch in the tract this exposure came from.')
122 self.schema.addField('weight', type=float,
123 doc='Weight for each exposure, used to make the CoaddPsf; should always be 1.')
124
125 def runQuantum(self, butlerQC, inputRefs, outputRefs):
126 inputs = butlerQC.get(inputRefs)
127 bbox = inputs.pop("bbox")
128 wcs = inputs.pop("wcs")
129 coaddExposures = inputs.pop('coaddExposures')
130 skymap = inputs.pop("skyMap")
131
132 # This should not happen with a properly configured execution context.
133 assert not inputs, "runQuantum got more inputs than expected"
134
135 results = self.getExposures(coaddExposures, bbox, skymap, wcs)
136 physical_filter = butlerQC.quantum.dataId["physical_filter"]
137 outputs = self.run(coaddExposures=results.coaddExposures,
138 bbox=bbox,
139 wcs=wcs,
140 dataIds=results.dataIds,
141 physical_filter=physical_filter)
142 butlerQC.put(outputs, outputRefs)
143
144 @deprecated(reason="Replaced by getExposures, which uses explicit arguments instead of a kwargs dict. "
145 "This method will be removed after v29.",
146 version="v29.0", category=FutureWarning)
147 def getOverlappingExposures(self, inputs):
148 return self.getExposures(inputs["coaddExposures"],
149 inputs["bbox"],
150 inputs["skyMap"],
151 inputs["wcs"])
152
153 def getExposures(self, coaddExposureHandles, bbox, skymap, wcs):
154 """Return a data structure containing the coadds that overlap the
155 specified bbox projected onto the sky, and a corresponding data
156 structure of their dataIds.
157 These are the appropriate inputs to this task's `run` method.
158
159 The spatial index in the butler registry has generous padding and often
160 supplies patches near, but not directly overlapping the desired region.
161 This method filters the inputs so that `run` does not have to read in
162 all possibly-matching coadd exposures.
163
164 Parameters
165 ----------
166 coaddExposureHandles : `iterable` \
167 [`lsst.daf.butler.DeferredDatasetHandle` of \
168 `lsst.afw.image.Exposure`]
169 Dataset handles to exposures that might overlap the desired
170 region.
171 bbox : `lsst.geom.Box2I`
172 Template bounding box of the pixel geometry onto which the
173 coaddExposures will be resampled.
174 skymap : `lsst.skymap.SkyMap`
175 Geometry of the tracts and patches the coadds are defined on.
176 wcs : `lsst.afw.geom.SkyWcs`
177 Template WCS onto which the coadds will be resampled.
178
179 Returns
180 -------
181 result : `lsst.pipe.base.Struct`
182 A struct with attributes:
183
184 ``coaddExposures``
185 Dict of coadd exposures that overlap the projected bbox,
186 indexed on tract id
187 (`dict` [`int`, `list` [`lsst.afw.image.Exposure`] ]).
188 ``dataIds``
189 Dict of data IDs of the coadd exposures that overlap the
190 projected bbox, indexed on tract id
191 (`dict` [`int`, `list [`lsst.daf.butler.DataCoordinate`] ]).
192
193 Raises
194 ------
195 NoWorkFound
196 Raised if no patches overlap the input detector bbox, or the input
197 WCS is None.
198 """
199 if wcs is None:
200 raise pipeBase.NoWorkFound("WCS is None; cannot find overlapping exposures.")
201
202 # Exposure's validPolygon would be more accurate
203 detectorPolygon = geom.Box2D(bbox)
204 overlappingArea = 0
205 coaddExposures = collections.defaultdict(list)
206 dataIds = collections.defaultdict(list)
207
208 for coaddRef in coaddExposureHandles:
209 dataId = coaddRef.dataId
210 patchWcs = skymap[dataId['tract']].getWcs()
211 patchBBox = skymap[dataId['tract']][dataId['patch']].getOuterBBox()
212 patchCorners = patchWcs.pixelToSky(geom.Box2D(patchBBox).getCorners())
213 patchPolygon = afwGeom.Polygon(wcs.skyToPixel(patchCorners))
214 if patchPolygon.intersection(detectorPolygon):
215 overlappingArea += patchPolygon.intersectionSingle(detectorPolygon).calculateArea()
216 self.log.info("Using template input tract=%s, patch=%s", dataId['tract'], dataId['patch'])
217 coaddExposures[dataId['tract']].append(coaddRef.get())
218 dataIds[dataId['tract']].append(dataId)
219
220 if not overlappingArea:
221 raise pipeBase.NoWorkFound('No patches overlap detector')
222
223 return pipeBase.Struct(coaddExposures=coaddExposures,
224 dataIds=dataIds)
225
226 @timeMethod
227 def run(self, *, coaddExposures, bbox, wcs, dataIds, physical_filter):
228 """Warp coadds from multiple tracts and patches to form a template to
229 subtract from a science image.
230
231 Tract and patch overlap regions are combined by a variance-weighted
232 average, and the variance planes are combined with the same weights,
233 not added in quadrature; the overlap regions are not statistically
234 independent, because they're derived from the same original data.
235 The PSF on the template is created by combining the CoaddPsf on each
236 template image into a meta-CoaddPsf.
237
238 Parameters
239 ----------
240 coaddExposures : `dict` [`int`, `list` [`lsst.afw.image.Exposure`]]
241 Coadds to be mosaicked, indexed on tract id.
242 bbox : `lsst.geom.Box2I`
243 Template Bounding box of the detector geometry onto which to
244 resample the ``coaddExposures``. Modified in-place to include the
245 template border.
246 wcs : `lsst.afw.geom.SkyWcs`
247 Template WCS onto which to resample the ``coaddExposures``.
248 dataIds : `dict` [`int`, `list` [`lsst.daf.butler.DataCoordinate`]]
249 Record of the tract and patch of each coaddExposure, indexed on
250 tract id.
251 physical_filter : `str`
252 Physical filter of the science image.
253
254 Returns
255 -------
256 result : `lsst.pipe.base.Struct`
257 A struct with attributes:
258
259 ``template``
260 A template coadd exposure assembled out of patches
261 (`lsst.afw.image.ExposureF`).
262
263 Raises
264 ------
265 NoWorkFound
266 If no coadds are found with sufficient un-masked pixels.
267 """
268 band, photoCalib = self._checkInputs(dataIds, coaddExposures)
269
270 bbox.grow(self.config.templateBorderSize)
271
272 warped = {}
273 catalogs = []
274 for tract in coaddExposures:
275 maskedImages, catalog, totalBox = self._makeExposureCatalog(coaddExposures[tract],
276 dataIds[tract])
277 # Combine images from individual patches together.
278 unwarped = self._merge(maskedImages, totalBox, catalog[0].wcs)
279 potentialInput = self.warper.warpExposure(wcs, unwarped, destBBox=bbox)
280 if not np.any(np.isfinite(potentialInput.image.array)):
281 self.log.info("No overlap from coadds in tract %s; not including in output.", tract)
282 continue
283 catalogs.append(catalog)
284 warped[tract] = potentialInput
285 warped[tract].setWcs(wcs)
286
287 if len(warped) == 0:
288 raise pipeBase.NoWorkFound("No patches found to overlap science exposure.")
289 template = self._merge([x.maskedImage for x in warped.values()], bbox, wcs)
290
291 # Make a single catalog containing all the inputs that were accepted.
292 catalog = afwTable.ExposureCatalog(self.schema)
293 catalog.reserve(sum([len(c) for c in catalogs]))
294 for c in catalogs:
295 catalog.extend(c)
296
297 template.setPsf(self._makePsf(template, catalog, wcs))
298 template.setFilter(afwImage.FilterLabel(band, physical_filter))
299 template.setPhotoCalib(photoCalib)
300 return pipeBase.Struct(template=template)
301
302 @staticmethod
303 def _checkInputs(dataIds, coaddExposures):
304 """Check that the all the dataIds are from the same band and that
305 the exposures all have the same photometric calibration.
306
307 Parameters
308 ----------
309 dataIds : `dict` [`int`, `list` [`lsst.daf.butler.DataCoordinate`]]
310 Record of the tract and patch of each coaddExposure.
311 coaddExposures : `dict` [`int`, `list` [`lsst.afw.image.Exposure`]]
312 Coadds to be mosaicked.
313
314 Returns
315 -------
316 band : `str`
317 Filter band of all the input exposures.
318 photoCalib : `lsst.afw.image.PhotoCalib`
319 Photometric calibration of all of the input exposures.
320
321 Raises
322 ------
323 RuntimeError
324 Raised if the bands or calibrations of the input exposures are not
325 all the same.
326 """
327 bands = set(dataId["band"] for tract in dataIds for dataId in dataIds[tract])
328 if len(bands) > 1:
329 raise RuntimeError(f"GetTemplateTask called with multiple bands: {bands}")
330 band = bands.pop()
331 photoCalibs = [exposure.photoCalib for exposures in coaddExposures.values() for exposure in exposures]
332 if not all([photoCalibs[0] == x for x in photoCalibs]):
333 msg = f"GetTemplateTask called with exposures with different photoCalibs: {photoCalibs}"
334 raise RuntimeError(msg)
335 photoCalib = photoCalibs[0]
336 return band, photoCalib
337
338 def _makeExposureCatalog(self, exposures, dataIds):
339 """Make an exposure catalog for one tract.
340
341 Parameters
342 ----------
343 exposures : `list` [`lsst.afw.image.Exposuref`]
344 Exposures to include in the catalog.
345 dataIds : `list` [`lsst.daf.butler.DataCoordinate`]
346 Data ids of each of the included exposures; must have "tract" and
347 "patch" entries.
348
349 Returns
350 -------
351 images : `list` [`lsst.afw.image.MaskedImage`]
352 MaskedImages of each of the input exposures, for warping.
353 catalog : `lsst.afw.table.ExposureCatalog`
354 Catalog of metadata for each exposure
355 totalBox : `lsst.geom.Box2I`
356 The union of the bounding boxes of all the input exposures.
357 """
358 catalog = afwTable.ExposureCatalog(self.schema)
359 catalog.reserve(len(exposures))
360 images = [exposure.maskedImage for exposure in exposures]
361 totalBox = geom.Box2I()
362 for coadd, dataId in zip(exposures, dataIds):
363 totalBox = totalBox.expandedTo(coadd.getBBox())
364 record = catalog.addNew()
365 record.setPsf(coadd.psf)
366 record.setWcs(coadd.wcs)
367 record.setPhotoCalib(coadd.photoCalib)
368 record.setBBox(coadd.getBBox())
369 record.setValidPolygon(afwGeom.Polygon(geom.Box2D(coadd.getBBox()).getCorners()))
370 record.set("tract", dataId["tract"])
371 record.set("patch", dataId["patch"])
372 # Weight is used by CoaddPsf, but the PSFs from overlapping patches
373 # should be very similar, so this value mostly shouldn't matter.
374 record.set("weight", 1)
375
376 return images, catalog, totalBox
377
378 @staticmethod
379 def _merge(maskedImages, bbox, wcs):
380 """Merge the images that came from one tract into one larger image,
381 ignoring NaN pixels and non-finite variance pixels from individual
382 exposures.
383
384 Parameters
385 ----------
386 maskedImages : `list` [`lsst.afw.image.MaskedImage`]
387 Images to be merged into one larger bounding box.
388 bbox : `lsst.geom.Box2I`
389 Bounding box defining the image to merge into.
390 wcs : `lsst.afw.geom.SkyWcs`
391 WCS of all of the input images to set on the output image.
392
393 Returns
394 -------
395 merged : `lsst.afw.image.MaskedImage`
396 Merged image with all of the inputs at their respective bbox
397 positions.
398 """
399 merged = afwImage.ExposureF(bbox, wcs)
400 weights = afwImage.ImageF(bbox)
401 for maskedImage in maskedImages:
402 # Catch both zero-value and NaN variance plane pixels
403 good = maskedImage.variance.array > 0
404 weight = afwImage.ImageF(maskedImage.getBBox())
405 weight.array[good] = maskedImage.variance.array[good]**(-0.5)
406 bad = np.isnan(maskedImage.image.array) | ~good
407 # Note that modifying the patch MaskedImage in place is fine;
408 # we're throwing it away at the end anyway.
409 maskedImage.image.array[bad] = 0.0
410 maskedImage.variance.array[bad] = 0.0
411 # Reset mask, too, since these pixels don't contribute to sum.
412 maskedImage.mask.array[bad] = 0
413 # Cannot use `merged.maskedImage *= weight` because that operator
414 # multiplies the variance by the weight twice; in this case
415 # `weight` are the exact values we want to scale by.
416 maskedImage.image *= weight
417 maskedImage.variance *= weight
418 merged.maskedImage[maskedImage.getBBox()] += maskedImage
419 weights[maskedImage.getBBox()] += weight
420
421 inverseWeights = np.zeros_like(weights.array)
422 good = weights.array > 0
423 inverseWeights[good] = 1/weights.array[good]
424
425 # Cannot use `merged.maskedImage *= inverseWeights` because that
426 # operator divides the variance by the weight twice; in this case
427 # `inverseWeights` are the exact values we want to scale by.
428 merged.image.array *= inverseWeights
429 merged.variance.array *= inverseWeights
430 merged.mask.array |= merged.mask.getPlaneBitMask("NO_DATA") * (inverseWeights == 0)
431
432 return merged
433
434 def _makePsf(self, template, catalog, wcs):
435 """Return a PSF containing the PSF at each of the input regions.
436
437 Note that although this includes all the exposures from the catalog,
438 the PSF knows which part of the template the inputs came from, so when
439 evaluated at a given position it will not include inputs that never
440 went in to those pixels.
441
442 Parameters
443 ----------
444 template : `lsst.afw.image.Exposure`
445 Generated template the PSF is for.
446 catalog : `lsst.afw.table.ExposureCatalog`
447 Catalog of exposures that went into the template that contains all
448 of the input PSFs.
449 wcs : `lsst.afw.geom.SkyWcs`
450 WCS of the template, to warp the PSFs to.
451
452 Returns
453 -------
454 coaddPsf : `lsst.meas.algorithms.CoaddPsf`
455 The meta-psf constructed from all of the input catalogs.
456 """
457 # CoaddPsf centroid not only must overlap image, but must overlap the
458 # part of image with data. Use centroid of region with data.
459 boolmask = template.mask.array & template.mask.getPlaneBitMask('NO_DATA') == 0
460 maskx = afwImage.makeMaskFromArray(boolmask.astype(afwImage.MaskPixel))
461 centerCoord = afwGeom.SpanSet.fromMask(maskx, 1).computeCentroid()
462
463 ctrl = self.config.coaddPsf.makeControl()
464 coaddPsf = CoaddPsf(catalog, wcs, centerCoord, ctrl.warpingKernelName, ctrl.cacheSize)
465 return coaddPsf
466
467
469 dimensions=("instrument", "visit", "detector", "skymap"),
470 defaultTemplates={"coaddName": "dcr",
471 "warpTypeSuffix": "",
472 "fakesType": ""}):
473 visitInfo = pipeBase.connectionTypes.Input(
474 doc="VisitInfo of calexp used to determine observing conditions.",
475 name="{fakesType}calexp.visitInfo",
476 storageClass="VisitInfo",
477 dimensions=("instrument", "visit", "detector"),
478 )
479 dcrCoadds = pipeBase.connectionTypes.Input(
480 doc="Input DCR template to match and subtract from the exposure",
481 name="{fakesType}dcrCoadd{warpTypeSuffix}",
482 storageClass="ExposureF",
483 dimensions=("tract", "patch", "skymap", "band", "subfilter"),
484 multiple=True,
485 deferLoad=True
486 )
487
488 def __init__(self, *, config=None):
489 super().__init__(config=config)
490 self.inputs.remove("coaddExposures")
491
492
493class GetDcrTemplateConfig(GetTemplateConfig,
494 pipelineConnections=GetDcrTemplateConnections):
495 numSubfilters = pexConfig.Field(
496 doc="Number of subfilters in the DcrCoadd.",
497 dtype=int,
498 default=3,
499 )
500 effectiveWavelength = pexConfig.Field(
501 doc="Effective wavelength of the filter in nm.",
502 optional=False,
503 dtype=float,
504 )
505 bandwidth = pexConfig.Field(
506 doc="Bandwidth of the physical filter.",
507 optional=False,
508 dtype=float,
509 )
510
511 def validate(self):
512 if self.effectiveWavelength is None or self.bandwidth is None:
513 raise ValueError("The effective wavelength and bandwidth of the physical filter "
514 "must be set in the getTemplate config for DCR coadds. "
515 "Required until transmission curves are used in DM-13668.")
516
517
518class GetDcrTemplateTask(GetTemplateTask):
519 ConfigClass = GetDcrTemplateConfig
520 _DefaultName = "getDcrTemplate"
521
522 def runQuantum(self, butlerQC, inputRefs, outputRefs):
523 inputs = butlerQC.get(inputRefs)
524 bbox = inputs.pop("bbox")
525 wcs = inputs.pop("wcs")
526 dcrCoaddExposureHandles = inputs.pop('dcrCoadds')
527 skymap = inputs.pop("skyMap")
528 visitInfo = inputs.pop("visitInfo")
529
530 # This should not happen with a properly configured execution context.
531 assert not inputs, "runQuantum got more inputs than expected"
532
533 results = self.getExposures(dcrCoaddExposureHandles, bbox, skymap, wcs, visitInfo)
534 physical_filter = butlerQC.quantum.dataId["physical_filter"]
535 outputs = self.run(coaddExposures=results.coaddExposures,
536 bbox=bbox,
537 wcs=wcs,
538 dataIds=results.dataIds,
539 physical_filter=physical_filter)
540 butlerQC.put(outputs, outputRefs)
541
542 @deprecated(reason="Replaced by getExposures, which uses explicit arguments instead of a kwargs dict. "
543 "This method will be removed after v29.",
544 version="v29.0", category=FutureWarning)
545 def getOverlappingExposures(self, inputs):
546 return self.getExposures(inputs["dcrCoadds"],
547 inputs["bbox"],
548 inputs["skyMap"],
549 inputs["wcs"],
550 inputs["visitInfo"])
551
552 def getExposures(self, dcrCoaddExposureHandles, bbox, skymap, wcs, visitInfo):
553 """Return lists of coadds and their corresponding dataIds that overlap
554 the detector.
555
556 The spatial index in the registry has generous padding and often
557 supplies patches near, but not directly overlapping the detector.
558 Filters inputs so that we don't have to read in all input coadds.
559
560 Parameters
561 ----------
562 dcrCoaddExposureHandles : `list` \
563 [`lsst.daf.butler.DeferredDatasetHandle` of \
564 `lsst.afw.image.Exposure`]
565 Data references to exposures that might overlap the detector.
566 bbox : `lsst.geom.Box2I`
567 Template Bounding box of the detector geometry onto which to
568 resample the coaddExposures.
569 skymap : `lsst.skymap.SkyMap`
570 Input definition of geometry/bbox and projection/wcs for
571 template exposures.
572 wcs : `lsst.afw.geom.SkyWcs`
573 Template WCS onto which to resample the coaddExposures.
574 visitInfo : `lsst.afw.image.VisitInfo`
575 Metadata for the science image.
576
577 Returns
578 -------
579 result : `lsst.pipe.base.Struct`
580 A struct with attibutes:
581
582 ``coaddExposures``
583 Dict of coadd exposures that overlap the projected bbox,
584 indexed on tract id
585 (`dict` [`int`, `list` [`lsst.afw.image.Exposure`] ]).
586 ``dataIds``
587 Dict of data IDs of the coadd exposures that overlap the
588 projected bbox, indexed on tract id
589 (`dict` [`int`, `list [`lsst.daf.butler.DataCoordinate`] ]).
590
591 Raises
592 ------
593 pipeBase.NoWorkFound
594 Raised if no patches overlatp the input detector bbox.
595 """
596 # Check that the patches actually overlap the detector
597 # Exposure's validPolygon would be more accurate
598 if wcs is None:
599 raise pipeBase.NoWorkFound("Exposure has no WCS; cannot create a template.")
600
601 detectorPolygon = geom.Box2D(bbox)
602 overlappingArea = 0
603 dataIds = collections.defaultdict(list)
604 patchList = dict()
605 for coaddRef in dcrCoaddExposureHandles:
606 dataId = coaddRef.dataId
607 patchWcs = skymap[dataId['tract']].getWcs()
608 patchBBox = skymap[dataId['tract']][dataId['patch']].getOuterBBox()
609 patchCorners = patchWcs.pixelToSky(geom.Box2D(patchBBox).getCorners())
610 patchPolygon = afwGeom.Polygon(wcs.skyToPixel(patchCorners))
611 if patchPolygon.intersection(detectorPolygon):
612 overlappingArea += patchPolygon.intersectionSingle(detectorPolygon).calculateArea()
613 self.log.info("Using template input tract=%s, patch=%s, subfilter=%s" %
614 (dataId['tract'], dataId['patch'], dataId["subfilter"]))
615 if dataId['tract'] in patchList:
616 patchList[dataId['tract']].append(dataId['patch'])
617 else:
618 patchList[dataId['tract']] = [dataId['patch'], ]
619 dataIds[dataId['tract']].append(dataId)
620
621 if not overlappingArea:
622 raise pipeBase.NoWorkFound('No patches overlap detector')
623
624 self.checkPatchList(patchList)
625
626 coaddExposures = self.getDcrModel(patchList, dcrCoaddExposureHandles, visitInfo)
627 return pipeBase.Struct(coaddExposures=coaddExposures,
628 dataIds=dataIds)
629
630 def checkPatchList(self, patchList):
631 """Check that all of the DcrModel subfilters are present for each
632 patch.
633
634 Parameters
635 ----------
636 patchList : `dict`
637 Dict of the patches containing valid data for each tract.
638
639 Raises
640 ------
641 RuntimeError
642 If the number of exposures found for a patch does not match the
643 number of subfilters.
644 """
645 for tract in patchList:
646 for patch in set(patchList[tract]):
647 if patchList[tract].count(patch) != self.config.numSubfilters:
648 raise RuntimeError("Invalid number of DcrModel subfilters found: %d vs %d expected",
649 patchList[tract].count(patch), self.config.numSubfilters)
650
651 def getDcrModel(self, patchList, coaddRefs, visitInfo):
652 """Build DCR-matched coadds from a list of exposure references.
653
654 Parameters
655 ----------
656 patchList : `dict`
657 Dict of the patches containing valid data for each tract.
658 coaddRefs : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
659 Data references to `~lsst.afw.image.Exposure` representing
660 DcrModels that overlap the detector.
661 visitInfo : `lsst.afw.image.VisitInfo`
662 Metadata for the science image.
663
664 Returns
665 -------
666 coaddExposures : `list` [`lsst.afw.image.Exposure`]
667 Coadd exposures that overlap the detector.
668 """
669 coaddExposures = collections.defaultdict(list)
670 for tract in patchList:
671 for patch in set(patchList[tract]):
672 coaddRefList = [coaddRef for coaddRef in coaddRefs
673 if _selectDataRef(coaddRef, tract, patch)]
674
675 dcrModel = DcrModel.fromQuantum(coaddRefList,
676 self.config.effectiveWavelength,
677 self.config.bandwidth,
678 self.config.numSubfilters)
679 coaddExposures[tract].append(dcrModel.buildMatchedExposure(visitInfo=visitInfo))
680 return coaddExposures
681
682
683def _selectDataRef(coaddRef, tract, patch):
684 condition = (coaddRef.dataId['tract'] == tract) & (coaddRef.dataId['patch'] == patch)
685 return condition
Cartesian polygons.
Definition Polygon.h:59
A group of labels for a filter in an exposure or coadd.
Definition FilterLabel.h:58
A floating-point coordinate rectangle geometry.
Definition Box.h:413
An integer coordinate rectangle.
Definition Box.h:55
CoaddPsf is the Psf derived to be used for non-PSF-matched Coadd images.
Definition CoaddPsf.h:58
_selectDataRef(coaddRef, tract, patch)
checkPatchList(self, patchList)
run(self, *, coaddExposures, bbox, wcs, dataIds, physical_filter)
_makeExposureCatalog(self, exposures, dataIds)
getDcrModel(self, patchList, coaddRefs, visitInfo)
_makePsf(self, template, catalog, wcs)