lsst.meas.algorithms 22.0.1-37-g74f620de+8bfad66fb6
loadReferenceObjects.py
Go to the documentation of this file.
2# LSST Data Management System
3#
4# Copyright 2008-2017 AURA/LSST.
5#
6# This product includes software developed by the
7# LSST Project (http://www.lsst.org/).
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 LSST License Statement and
20# the GNU General Public License along with this program. If not,
21# see <https://www.lsstcorp.org/LegalNotices/>.
22#
23
24__all__ = ["getRefFluxField", "getRefFluxKeys", "LoadReferenceObjectsTask", "LoadReferenceObjectsConfig",
25 "ReferenceObjectLoader"]
26
27import abc
28import logging
29
30import astropy.time
31import astropy.units
32import numpy
33
34import lsst.geom as geom
35import lsst.afw.table as afwTable
36import lsst.pex.config as pexConfig
37import lsst.pipe.base as pipeBase
38from lsst import sphgeom
39from lsst.daf.base import PropertyList
40from lsst.utils.timer import timeMethod
41
42
43def isOldFluxField(name, units):
44 """Return True if this name/units combination corresponds to an
45 "old-style" reference catalog flux field.
46 """
47 unitsCheck = units != 'nJy' # (units == 'Jy' or units == '' or units == '?')
48 isFlux = name.endswith('_flux')
49 isFluxSigma = name.endswith('_fluxSigma')
50 isFluxErr = name.endswith('_fluxErr')
51 return (isFlux or isFluxSigma or isFluxErr) and unitsCheck
52
53
55 """Return True if the units of all flux and fluxErr are correct (nJy).
56 """
57 for field in schema:
58 if isOldFluxField(field.field.getName(), field.field.getUnits()):
59 return False
60 return True
61
62
64 """"Return the format version stored in a reference catalog header.
65
66 Parameters
67 ----------
69 Reference catalog to inspect.
70
71 Returns
72 -------
73 version : `int`
74 Format verison integer. Returns `0` if the catalog has no metadata
75 or the metadata does not include a "REFCAT_FORMAT_VERSION" key.
76 """
77 md = refCat.getMetadata()
78 if md is None:
79 return 0
80 try:
81 return md.getScalar("REFCAT_FORMAT_VERSION")
82 except KeyError:
83 return 0
84
85
86def convertToNanojansky(catalog, log, doConvert=True):
87 """Convert fluxes in a catalog from jansky to nanojansky.
88
89 Parameters
90 ----------
92 The catalog to convert.
93 log : `lsst.log.Log` or `logging.Logger`
94 Log to send messages to.
95 doConvert : `bool`, optional
96 Return a converted catalog, or just identify the fields that need to be converted?
97 This supports the "write=False" mode of `bin/convert_to_nJy.py`.
98
99 Returns
100 -------
101 catalog : `lsst.afw.table.SimpleCatalog` or None
102 The converted catalog, or None if ``doConvert`` is False.
103
104 Notes
105 -----
106 Support for old units in reference catalogs will be removed after the
107 release of late calendar year 2019.
108 Use `meas_algorithms/bin/convert_to_nJy.py` to update your reference catalog.
109 """
110 # Do not share the AliasMap: for refcats, that gets created when the
111 # catalog is read from disk and should not be propagated.
112 mapper = afwTable.SchemaMapper(catalog.schema, shareAliasMap=False)
113 mapper.addMinimalSchema(afwTable.SimpleTable.makeMinimalSchema())
114 input_fields = []
115 output_fields = []
116 for field in catalog.schema:
117 oldName = field.field.getName()
118 oldUnits = field.field.getUnits()
119 if isOldFluxField(oldName, oldUnits):
120 units = 'nJy'
121 # remap Sigma flux fields to Err, so we can drop the alias
122 if oldName.endswith('_fluxSigma'):
123 name = oldName.replace('_fluxSigma', '_fluxErr')
124 else:
125 name = oldName
126 newField = afwTable.Field[field.dtype](name, field.field.getDoc(), units)
127 mapper.addMapping(field.getKey(), newField)
128 input_fields.append(field.field)
129 output_fields.append(newField)
130 else:
131 mapper.addMapping(field.getKey())
132
133 fluxFieldsStr = '; '.join("(%s, '%s')" % (field.getName(), field.getUnits()) for field in input_fields)
134
135 if doConvert:
136 newSchema = mapper.getOutputSchema()
137 output = afwTable.SimpleCatalog(newSchema)
138 output.extend(catalog, mapper=mapper)
139 for field in output_fields:
140 output[field.getName()] *= 1e9
141 log.info("Converted refcat flux fields to nJy (name, units): %s", fluxFieldsStr)
142 return output
143 else:
144 log.info("Found old-style refcat flux fields (name, units): %s", fluxFieldsStr)
145 return None
146
147
149 """This is a private helper class which filters catalogs by
150 row based on the row being inside the region used to initialize
151 the class.
152
153 Parameters
154 ----------
155 region : `lsst.sphgeom.Region`
156 The spatial region which all objects should lie within
157 """
158 def __init__(self, region):
159 self.regionregion = region
160
161 def __call__(self, refCat, catRegion):
162 """This call method on an instance of this class takes in a reference
163 catalog, and the region from which the catalog was generated.
164
165 If the catalog region is entirely contained within the region used to
166 initialize this class, then all the entries in the catalog must be
167 within the region and so the whole catalog is returned.
168
169 If the catalog region is not entirely contained, then the location for
170 each record is tested against the region used to initialize the class.
171 Records which fall inside this region are added to a new catalog, and
172 this catalog is then returned.
173
174 Parameters
175 ---------
177 SourceCatalog to be filtered.
178 catRegion : `lsst.sphgeom.Region`
179 Region in which the catalog was created
180 """
181 if catRegion.isWithin(self.regionregion):
182 # no filtering needed, region completely contains refcat
183 return refCat
184
185 filteredRefCat = type(refCat)(refCat.table)
186 for record in refCat:
187 if self.regionregion.contains(record.getCoord().getVector()):
188 filteredRefCat.append(record)
189 return filteredRefCat
190
191
192class LoadReferenceObjectsConfig(pexConfig.Config):
193 pixelMargin = pexConfig.RangeField(
194 doc="Padding to add to 4 all edges of the bounding box (pixels)",
195 dtype=int,
196 default=250,
197 min=0,
198 )
199 anyFilterMapsToThis = pexConfig.Field(
200 doc=("Always use this reference catalog filter, no matter whether or what filter name is "
201 "supplied to the loader. Effectively a trivial filterMap: map all filter names to this filter."
202 " This can be set for purely-astrometric catalogs (e.g. Gaia DR2) where there is only one "
203 "reasonable choice for every camera filter->refcat mapping, but not for refcats used for "
204 "photometry, which need a filterMap and/or colorterms/transmission corrections."),
205 dtype=str,
206 default=None,
207 optional=True
208 )
209 filterMap = pexConfig.DictField(
210 doc=("Mapping of camera filter name: reference catalog filter name; "
211 "each reference filter must exist in the refcat."
212 " Note that this does not perform any bandpass corrections: it is just a lookup."),
213 keytype=str,
214 itemtype=str,
215 default={},
216 )
217 requireProperMotion = pexConfig.Field(
218 doc="Require that the fields needed to correct proper motion "
219 "(epoch, pm_ra and pm_dec) are present?",
220 dtype=bool,
221 default=False,
222 )
223
224 def validate(self):
225 super().validate()
226 if self.filterMapfilterMap != {} and self.anyFilterMapsToThisanyFilterMapsToThis is not None:
227 msg = "`filterMap` and `anyFilterMapsToThis` are mutually exclusive"
228 raise pexConfig.FieldValidationError(LoadReferenceObjectsConfig.anyFilterMapsToThis,
229 self, msg)
230
231
233 """Base class for reference object loaders, to facilitate gen2/gen3 code
234 sharing.
235 """
236 ConfigClass = LoadReferenceObjectsConfig
237
238 def applyProperMotions(self, catalog, epoch):
239 """Apply proper motion correction to a reference catalog.
240
241 Adjust position and position error in the ``catalog``
242 for proper motion to the specified ``epoch``,
243 modifying the catalog in place.
244
245 Parameters
246 ----------
248 Catalog of positions, containing at least these fields:
249
250 - Coordinates, retrieved by the table's coordinate key.
251 - ``coord_raErr`` : Error in Right Ascension (rad).
252 - ``coord_decErr`` : Error in Declination (rad).
253 - ``pm_ra`` : Proper motion in Right Ascension (rad/yr,
254 East positive)
255 - ``pm_raErr`` : Error in ``pm_ra`` (rad/yr), optional.
256 - ``pm_dec`` : Proper motion in Declination (rad/yr,
257 North positive)
258 - ``pm_decErr`` : Error in ``pm_dec`` (rad/yr), optional.
259 - ``epoch`` : Mean epoch of object (an astropy.time.Time)
260 epoch : `astropy.time.Time`
261 Epoch to which to correct proper motion.
262 If None, do not apply PM corrections or raise if
263 ``config.requireProperMotion`` is True.
264
265 Raises
266 ------
267 RuntimeError
268 Raised if ``config.requireProperMotion`` is set but we cannot
269 apply the proper motion correction for some reason.
270 """
271 if epoch is None:
272 if self.config.requireProperMotion:
273 raise RuntimeError("requireProperMotion=True but epoch not provided to loader.")
274 else:
275 self.log.debug("No epoch provided: not applying proper motion corrections to refcat.")
276 return
277
278 # Warn/raise for a catalog in an incorrect format, if epoch was specified.
279 if ("pm_ra" in catalog.schema
280 and not isinstance(catalog.schema["pm_ra"].asKey(), afwTable.KeyAngle)):
281 if self.config.requireProperMotion:
282 raise RuntimeError("requireProperMotion=True but refcat pm_ra field is not an Angle.")
283 else:
284 self.log.warning("Reference catalog pm_ra field is not an Angle; cannot apply proper motion.")
285 return
286
287 if ("epoch" not in catalog.schema or "pm_ra" not in catalog.schema):
288 if self.config.requireProperMotion:
289 raise RuntimeError("requireProperMotion=True but PM data not available from catalog.")
290 else:
291 self.log.warning("Proper motion correction not available for this reference catalog.")
292 return
293
294 applyProperMotionsImpl(self.log, catalog, epoch)
295
296
298 """This class facilitates loading reference catalogs with gen 3 middleware
299
300 The middleware preflight solver will create a list of datarefs that may
301 possibly overlap a given region. These datarefs are then used to construct
302 and instance of this class. The class instance should then be passed into
303 a task which needs reference catalogs. These tasks should then determine
304 the exact region of the sky reference catalogs will be loaded for, and
305 call a corresponding method to load the reference objects.
306 """
307 def __init__(self, dataIds, refCats, config, log=None):
308 """ Constructs an instance of ReferenceObjectLoader
309
310 Parameters
311 ----------
312 dataIds : iterable of `lsst.daf.butler.DataIds`
313 An iterable object of DataSetRefs which point to reference catalogs
314 in a gen 3 repository.
315 refCats : iterable of `lsst.daf.butler.DeferedDatasetHandle`
316 Handles to load refCats on demand
318 Configuration for the loader.
319 log : `lsst.log.Log`, `logging.Logger` or `None`, optional
320 Logger object used to write out messages. If `None` a default
321 logger will be used.
322 """
323 self.dataIdsdataIds = dataIds
324 self.refCatsrefCats = refCats
325 self.loglog = log or logging.getLogger(__name__).getChild("ReferenceObjectLoader")
326 self.configconfig = config
327
328 @staticmethod
329 def _makeBoxRegion(BBox, wcs, BBoxPadding):
330 outerLocalBBox = geom.Box2D(BBox)
331 innerLocalBBox = geom.Box2D(BBox)
332
333 # Grow the bounding box to allow for effects not fully captured by the
334 # wcs provided (which represents the current best-guess wcs solution
335 # associated with the dataset for which the calibration is to be
336 # computed using the loaded and trimmed reference catalog being defined
337 # here). These effects could include pointing errors and/or an
338 # insufficient optical distorition model for the instrument. The idea
339 # is to ensure the spherical geometric region created contains the
340 # entire region covered by the bbox.
341 # Also create an inner region that is sure to be inside the bbox.
342 outerLocalBBox.grow(BBoxPadding)
343 innerLocalBBox.grow(-1*BBoxPadding)
344
345 # Handle the case where the inner bounding box shrank to a zero sized
346 # region (which will be the case if the shrunken size of either
347 # dimension is less than or equal to zero). In this case, the inner
348 # bounding box is set to the original input bounding box. This is
349 # probably not the best way to handle an empty inner bounding box, but
350 # it is what the calling code currently expects.
351 if innerLocalBBox.getDimensions() == geom.Extent2D(0, 0):
352 innerLocalBBox = geom.Box2D(BBox)
353
354 # Convert the corners of the bounding boxes to sky coordinates.
355 innerBoxCorners = innerLocalBBox.getCorners()
356 innerSphCorners = [wcs.pixelToSky(corner).getVector() for corner in innerBoxCorners]
357 innerSkyRegion = sphgeom.ConvexPolygon(innerSphCorners)
358
359 outerBoxCorners = outerLocalBBox.getCorners()
360 outerSphCorners = [wcs.pixelToSky(corner).getVector() for corner in outerBoxCorners]
361 outerSkyRegion = sphgeom.ConvexPolygon(outerSphCorners)
362
363 return innerSkyRegion, outerSkyRegion, innerSphCorners, outerSphCorners
364
365 def loadPixelBox(self, bbox, wcs, filterName, epoch=None, photoCalib=None,
366 bboxToSpherePadding=100):
367 """Load reference objects that are within a pixel-based rectangular
368 region.
369
370 This algorithm works by creating a spherical box whose corners
371 correspond to the WCS converted corners of the input bounding box
372 (possibly padded). It then defines a filtering function which looks at
373 the pixel position of the reference objects and accepts only those that
374 lie within the specified bounding box.
375
376 The spherical box region and filtering function are passed to the
377 generic loadRegion method which loads and filters the reference objects
378 from the datastore and returns a single catalog containing the filtered
379 set of reference objects.
380
381 Parameters
382 ----------
383 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D`
384 Box which bounds a region in pixel space.
386 Wcs object defining the pixel to sky (and inverse) transform for
387 the supplied ``bbox``.
388 filterName : `str`
389 Name of camera filter.
390 epoch : `astropy.time.Time` or `None`, optional
391 Epoch to which to correct proper motion and parallax, or `None`
392 to not apply such corrections.
393 photoCalib : `None`
394 Deprecated and ignored, only included for api compatibility.
395 bboxToSpherePadding : `int`, optional
396 Padding to account for translating a set of corners into a
397 spherical (convex) boundary that is certain to encompase the
398 enitre area covered by the box.
399
400 Returns
401 -------
402 referenceCatalog : `lsst.afw.table.SimpleCatalog`
403 Catalog containing reference objects inside the specified bounding
404 box (padded by self.configconfig.pixelMargin).
405
406 Raises
407 ------
408 RuntimeError
409 Raised if no reference catalogs could be found for the specified
410 region.
411 TypeError
412 Raised if the loaded reference catalogs do not have matching
413 schemas.
414 """
415 paddedBbox = geom.Box2D(bbox)
416 paddedBbox.grow(self.configconfig.pixelMargin)
417 innerSkyRegion, outerSkyRegion, _, _ = self._makeBoxRegion_makeBoxRegion(paddedBbox, wcs, bboxToSpherePadding)
418
419 def _filterFunction(refCat, region):
420 # Perform an initial "pre filter" step based on the refCat coords
421 # and the outerSkyRegion created from the self.config.pixelMargin-
422 # paddedBbox plus an "extra" padding of bboxToSpherePadding and the
423 # raw wcs. This should ensure a large enough projected area on the
424 # sky that accounts for any projection/distortion issues, but small
425 # enough to filter out loaded reference objects that lie well
426 # beyond the projected detector of interest. This step is required
427 # due to the very local nature of the wcs available for the
428 # sky <--> pixel conversions.
429 preFiltFunc = _FilterCatalog(outerSkyRegion)
430 refCat = preFiltFunc(refCat, region)
431
432 # Add columns to the pre-filtered reference catalog relating their
433 # coordinates to equivalent pixel positions for the wcs provided
434 # and use to populate those columns.
435 refCat = self.remapReferenceCatalogSchemaremapReferenceCatalogSchema(refCat, position=True)
436 afwTable.updateRefCentroids(wcs, refCat)
437 # No need to filter the catalog if it is entirely contained in the
438 # region defined by the inner sky region.
439 if innerSkyRegion.contains(region):
440 return refCat
441 # Create a new reference catalog, and populate it only with records
442 # that fall inside the padded bbox.
443 filteredRefCat = type(refCat)(refCat.table)
444 centroidKey = afwTable.Point2DKey(refCat.schema['centroid'])
445 for record in refCat:
446 pixCoords = record[centroidKey]
447 if paddedBbox.contains(geom.Point2D(pixCoords)):
448 filteredRefCat.append(record)
449 return filteredRefCat
450 return self.loadRegionloadRegion(outerSkyRegion, filterName, filtFunc=_filterFunction, epoch=epoch)
451
452 def loadRegion(self, region, filterName, filtFunc=None, epoch=None):
453 """Load reference objects within a specified region.
454
455 This function loads the DataIds used to construct an instance of this
456 class which intersect or are contained within the specified region. The
457 reference catalogs which intersect but are not fully contained within
458 the input region are further filtered by the specified filter function.
459 This function returns a single source catalog containing all reference
460 objects inside the specified region.
461
462 Parameters
463 ----------
464 region : `lsst.sphgeom.Region`
465 This can be any type that is derived from `lsst.sphgeom.Region` and
466 should define the spatial region for which reference objects are to
467 be loaded.
468 filtFunc : callable or `None`, optional
469 This optional parameter should be a callable object that takes a
470 reference catalog and its corresponding region as parameters,
471 filters the catalog by some criteria and returns the filtered
472 reference catalog. If `None`, an internal filter function is used
473 which filters according to if a reference object falls within the
474 input region.
475 filterName : `str`
476 Name of camera filter.
477 epoch : `astropy.time.Time` or `None`, optional
478 Epoch to which to correct proper motion and parallax, or `None` to
479 not apply such corrections.
480
481 Returns
482 -------
483 referenceCatalog : `lsst.afw.table.SourceCatalog`
484 Catalog containing reference objects which intersect the input region,
485 filtered by the specified filter function.
486
487 Raises
488 ------
489 RuntimeError
490 Raised if no reference catalogs could be found for the specified
491 region.
492 TypeError
493 Raised if the loaded reference catalogs do not have matching
494 schemas.
495 """
496 regionLat = region.getBoundingBox().getLat()
497 regionLon = region.getBoundingBox().getLon()
498 self.loglog.info("Loading reference objects from %s in region bounded by "
499 "[%.8f, %.8f], [%.8f, %.8f] RA Dec",
500 # Name of refcat we're loading from is the datasetType.
501 self.refCatsrefCats[0].ref.datasetType.name,
502 regionLon.getA().asDegrees(), regionLon.getB().asDegrees(),
503 regionLat.getA().asDegrees(), regionLat.getB().asDegrees())
504 if filtFunc is None:
505 filtFunc = _FilterCatalog(region)
506 # filter out all the regions supplied by the constructor that do not overlap
507 overlapList = []
508 for dataId, refCat in zip(self.dataIdsdataIds, self.refCatsrefCats):
509 # SphGeom supports some objects intersecting others, but is not symmetric,
510 # try the intersect operation in both directions
511 try:
512 intersects = dataId.region.intersects(region)
513 except TypeError:
514 intersects = region.intersects(dataId.region)
515
516 if intersects:
517 overlapList.append((dataId, refCat))
518
519 if len(overlapList) == 0:
520 raise RuntimeError("No reference tables could be found for input region")
521
522 firstCat = overlapList[0][1].get()
523 refCat = filtFunc(firstCat, overlapList[0][0].region)
524 trimmedAmount = len(firstCat) - len(refCat)
525
526 # Load in the remaining catalogs
527 for dataId, inputRefCat in overlapList[1:]:
528 tmpCat = inputRefCat.get()
529
530 if tmpCat.schema != firstCat.schema:
531 raise TypeError("Reference catalogs have mismatching schemas")
532
533 filteredCat = filtFunc(tmpCat, dataId.region)
534 refCat.extend(filteredCat)
535 trimmedAmount += len(tmpCat) - len(filteredCat)
536
537 self.loglog.debug("Trimmed %d refCat objects lying outside padded region, leaving %d",
538 trimmedAmount, len(refCat))
539 self.loglog.info("Loaded %d reference objects", len(refCat))
540
541 # Ensure that the loaded reference catalog is continuous in memory
542 if not refCat.isContiguous():
543 refCat = refCat.copy(deep=True)
544
545 self.applyProperMotionsapplyProperMotions(refCat, epoch)
546
547 # Verify the schema is in the correct units and has the correct version; automatically convert
548 # it with a warning if this is not the case.
549 if not hasNanojanskyFluxUnits(refCat.schema) or not getFormatVersionFromRefCat(refCat) >= 1:
550 self.loglog.warning("Found version 0 reference catalog with old style units in schema.")
551 self.loglog.warning("run `meas_algorithms/bin/convert_refcat_to_nJy.py` to convert fluxes to nJy.")
552 self.loglog.warning("See RFC-575 for more details.")
553 refCat = convertToNanojansky(refCat, self.loglog)
554
555 expandedCat = self.remapReferenceCatalogSchemaremapReferenceCatalogSchema(refCat, position=True)
556
557 # Add flux aliases
558 expandedCat = self.addFluxAliasesaddFluxAliases(expandedCat, self.configconfig.filterMap)
559
560 # Ensure that the returned reference catalog is continuous in memory
561 if not expandedCat.isContiguous():
562 expandedCat = expandedCat.copy(deep=True)
563
564 fluxField = getRefFluxField(expandedCat.schema, filterName)
565 return pipeBase.Struct(refCat=expandedCat, fluxField=fluxField)
566
567 def loadSkyCircle(self, ctrCoord, radius, filterName, epoch=None):
568 """Load reference objects that lie within a circular region on the sky.
569
570 This method constructs a circular region from an input center and
571 angular radius, loads reference catalogs which are contained in or
572 intersect the circle, and filters reference catalogs which intersect
573 down to objects which lie within the defined circle.
574
575 Parameters
576 ----------
577 ctrCoord : `lsst.geom.SpherePoint`
578 Point defining the center of the circular region.
579 radius : `lsst.geom.Angle`
580 Defines the angular radius of the circular region.
581 filterName : `str`
582 Name of camera filter.
583 epoch : `astropy.time.Time` or `None`, optional
584 Epoch to which to correct proper motion and parallax, or `None` to
585 not apply such corrections.
586
587 Returns
588 -------
589 referenceCatalog : `lsst.afw.table.SourceCatalog`
590 Catalog containing reference objects inside the specified search
591 circle.
592 """
593 centerVector = ctrCoord.getVector()
594 sphRadius = sphgeom.Angle(radius.asRadians())
595 circularRegion = sphgeom.Circle(centerVector, sphRadius)
596 return self.loadRegionloadRegion(circularRegion, filterName, epoch=epoch)
597
598 def joinMatchListWithCatalog(self, matchCat, sourceCat):
599 """Relink an unpersisted match list to sources and reference objects.
600
601 A match list is persisted and unpersisted as a catalog of IDs
602 produced by afw.table.packMatches(), with match metadata
603 (as returned by the astrometry tasks) in the catalog's metadata
604 attribute. This method converts such a match catalog into a match
605 list, with links to source records and reference object records.
606
607 Parameters
608 ----------
609 matchCat : `lsst.afw.table.BaseCatalog`
610 Unpersisted packed match list.
611 ``matchCat.table.getMetadata()`` must contain match metadata,
612 as returned by the astrometry tasks.
613 sourceCat : `lsst.afw.table.SourceCatalog`
614 Source catalog. As a side effect, the catalog will be sorted
615 by ID.
616
617 Returns
618 -------
620 Match list.
621 """
622 return joinMatchListWithCatalogImpl(self, matchCat, sourceCat)
623
624 def getMetadataBox(self, bbox, wcs, filterName, photoCalib=None, epoch=None,
625 bboxToSpherePadding=100):
626 """Return metadata about the load
627
628 This metadata is used for reloading the catalog (e.g., for
629 reconstituting a normalised match list.)
630
631 Parameters
632 ----------
633 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D`
634 Bounding box for the pixels.
636 The WCS object associated with ``bbox``.
637 filterName : `str`
638 Name of the camera filter.
639 photoCalib : `None`
640 Deprecated, only included for api compatibility.
641 epoch : `astropy.time.Time` or `None`, optional
642 Epoch to which to correct proper motion and parallax, or `None` to
643 not apply such corrections.
644 bboxToSpherePadding : `int`, optional
645 Padding to account for translating a set of corners into a
646 spherical (convex) boundary that is certain to encompase the
647 enitre area covered by the box.
648
649 Returns
650 -------
652 The metadata detailing the search parameters used for this
653 dataset.
654 """
655 paddedBbox = geom.Box2D(bbox)
656 paddedBbox.grow(self.configconfig.pixelMargin)
657 _, _, innerCorners, outerCorners = self._makeBoxRegion_makeBoxRegion(paddedBbox, wcs, bboxToSpherePadding)
658 md = PropertyList()
659 for box, corners in zip(("INNER", "OUTER"), (innerCorners, outerCorners)):
660 for (name, corner) in zip(("UPPER_LEFT", "UPPER_RIGHT", "LOWER_LEFT", "LOWER_RIGHT"),
661 corners):
662 md.add(f"{box}_{name}_RA", geom.SpherePoint(corner).getRa().asDegrees(), f"{box}_corner")
663 md.add(f"{box}_{name}_DEC", geom.SpherePoint(corner).getDec().asDegrees(), f"{box}_corner")
664 md.add("SMATCHV", 1, 'SourceMatchVector version number')
665 md.add('FILTER', filterName, 'filter name for photometric data')
666 md.add('EPOCH', "NONE" if epoch is None else epoch.mjd, 'Epoch (TAI MJD) for catalog')
667 return md
668
669 @staticmethod
670 def getMetadataCircle(coord, radius, filterName, photoCalib=None, epoch=None):
671 """Return metadata about the load.
672
673 This metadata is used for reloading the catalog (e.g. for
674 reconstituting a normalized match list.)
675
676 Parameters
677 ----------
678 coord : `lsst.geom.SpherePoint`
679 ICRS center of the search region.
680 radius : `lsst.geom.Angle`
681 Radius of the search region.
682 filterName : `str`
683 Name of the camera filter.
684 photoCalib : `None`
685 Deprecated, only included for api compatibility.
686 epoch : `astropy.time.Time` or `None`, optional
687 Epoch to which to correct proper motion and parallax, or `None` to
688 not apply such corrections.
689
690 Returns
691 -------
693 """
694 md = PropertyList()
695 md.add('RA', coord.getRa().asDegrees(), 'field center in degrees')
696 md.add('DEC', coord.getDec().asDegrees(), 'field center in degrees')
697 md.add('RADIUS', radius.asDegrees(), 'field radius in degrees, minimum')
698 md.add('SMATCHV', 1, 'SourceMatchVector version number')
699 md.add('FILTER', filterName, 'filter name for photometric data')
700 md.add('EPOCH', "NONE" if epoch is None else epoch.mjd, 'Epoch (TAI MJD) for catalog')
701 return md
702
703 @staticmethod
704 def addFluxAliases(refCat, filterReferenceMap):
705 """Add flux columns and aliases for camera to reference mapping.
706
707 Creates a new catalog containing the information of the input refCat
708 as well as added flux columns and aliases between camera and reference
709 fluxes.
710
711 Parameters
712 ----------
714 Catalog of reference objects
715 filterReferenceMap : `dict` of `str`
716 Dictionary with keys corresponding to a filter name and values
717 which correspond to the name of the reference filter.
718
719 Returns
720 -------
722 Reference catalog with columns added to track reference filters.
723
724 Raises
725 ------
726 `RuntimeError`
727 If the specified reference filter name is not specifed as a
728 key in the reference filter map.
729 """
730 refCat = ReferenceObjectLoader.remapReferenceCatalogSchema(refCat,
731 filterNameList=filterReferenceMap.keys())
732 aliasMap = refCat.schema.getAliasMap()
733 if filterReferenceMap is None:
734 filterReferenceMap = {}
735 for filterName, refFilterName in filterReferenceMap.items():
736 if refFilterName:
737 camFluxName = filterName + "_camFlux"
738 refFluxName = refFilterName + "_flux"
739 if refFluxName not in refCat.schema:
740 raise RuntimeError("Unknown reference filter %s" % (refFluxName,))
741 aliasMap.set(camFluxName, refFluxName)
742
743 refFluxErrName = refFluxName + "Err"
744 camFluxErrName = camFluxName + "Err"
745 aliasMap.set(camFluxErrName, refFluxErrName)
746
747 return refCat
748
749 @staticmethod
750 def remapReferenceCatalogSchema(refCat, *, filterNameList=None, position=False, photometric=False):
751 """This function takes in a reference catalog and creates a new catalog with additional
752 columns defined the remaining function arguments.
753
754 Parameters
755 ----------
757 Reference catalog to map to new catalog
758
759 Returns
760 -------
761 expandedCat : `lsst.afw.table.SimpleCatalog`
762 Deep copy of input reference catalog with additional columns added
763 """
764 mapper = afwTable.SchemaMapper(refCat.schema, True)
765 mapper.addMinimalSchema(refCat.schema, True)
766 mapper.editOutputSchema().disconnectAliases()
767 if filterNameList:
768 for filterName in filterNameList:
769 mapper.editOutputSchema().addField(f"{filterName}_flux",
770 type=numpy.float64,
771 doc=f"flux in filter {filterName}",
772 units="Jy"
773 )
774 mapper.editOutputSchema().addField(f"{filterName}_fluxErr",
775 type=numpy.float64,
776 doc=f"flux uncertanty in filter {filterName}",
777 units="Jy"
778 )
779
780 if position:
781 mapper.editOutputSchema().addField("centroid_x", type=float, doReplace=True)
782 mapper.editOutputSchema().addField("centroid_y", type=float, doReplace=True)
783 mapper.editOutputSchema().addField("hasCentroid", type="Flag", doReplace=True)
784 mapper.editOutputSchema().getAliasMap().set("slot_Centroid", "centroid")
785
786 if photometric:
787 mapper.editOutputSchema().addField("photometric",
788 type="Flag",
789 doc="set if the object can be used for photometric"
790 "calibration",
791 )
792 mapper.editOutputSchema().addField("resolved",
793 type="Flag",
794 doc="set if the object is spatially resolved"
795 )
796 mapper.editOutputSchema().addField("variable",
797 type="Flag",
798 doc="set if the object has variable brightness"
799 )
800
801 expandedCat = afwTable.SimpleCatalog(mapper.getOutputSchema())
802 expandedCat.setMetadata(refCat.getMetadata())
803 expandedCat.extend(refCat, mapper=mapper)
804
805 return expandedCat
806
807
808def getRefFluxField(schema, filterName):
809 """Get the name of a flux field from a schema.
810
811 return the alias of "anyFilterMapsToThis", if present
812 else:
813 return "*filterName*_camFlux" if present
814 else return "*filterName*_flux" if present (camera filter name
815 matches reference filter name)
816 else throw RuntimeError
817
818 Parameters
819 ----------
820 schema : `lsst.afw.table.Schema`
821 Reference catalog schema.
822 filterName : `str`
823 Name of camera filter.
824
825 Returns
826 -------
827 fluxFieldName : `str`
828 Name of flux field.
829
830 Raises
831 ------
832 RuntimeError
833 If an appropriate field is not found.
834 """
835 if not isinstance(schema, afwTable.Schema):
836 raise RuntimeError("schema=%s is not a schema" % (schema,))
837 try:
838 return schema.getAliasMap().get("anyFilterMapsToThis")
839 except LookupError:
840 pass # try the filterMap next
841
842 fluxFieldList = [filterName + "_camFlux", filterName + "_flux"]
843 for fluxField in fluxFieldList:
844 if fluxField in schema:
845 return fluxField
846
847 raise RuntimeError("Could not find flux field(s) %s" % (", ".join(fluxFieldList)))
848
849
850def getRefFluxKeys(schema, filterName):
851 """Return keys for flux and flux error.
852
853 Parameters
854 ----------
855 schema : `lsst.afw.table.Schema`
856 Reference catalog schema.
857 filterName : `str`
858 Name of camera filter.
859
860 Returns
861 -------
862 keys : `tuple` of (`lsst.afw.table.Key`, `lsst.afw.table.Key`)
863 Two keys:
864
865 - flux key
866 - flux error key, if present, else None
867
868 Raises
869 ------
870 RuntimeError
871 If flux field not found.
872 """
873 fluxField = getRefFluxField(schema, filterName)
874 fluxErrField = fluxField + "Err"
875 fluxKey = schema[fluxField].asKey()
876 try:
877 fluxErrKey = schema[fluxErrField].asKey()
878 except Exception:
879 fluxErrKey = None
880 return (fluxKey, fluxErrKey)
881
882
883class LoadReferenceObjectsTask(pipeBase.Task, ReferenceObjectLoaderBase, metaclass=abc.ABCMeta):
884 """Abstract base class to load objects from reference catalogs.
885 """
886 _DefaultName = "LoadReferenceObjects"
887
888 def __init__(self, butler=None, *args, **kwargs):
889 """Construct a LoadReferenceObjectsTask
890
891 Parameters
892 ----------
893 butler : `lsst.daf.persistence.Butler`
894 Data butler, for access reference catalogs.
895 """
896 pipeBase.Task.__init__(self, *args, **kwargs)
897 self.butlerbutler = butler
898
899 @timeMethod
900 def loadPixelBox(self, bbox, wcs, filterName, photoCalib=None, epoch=None):
901 """Load reference objects that overlap a rectangular pixel region.
902
903 Parameters
904 ----------
905 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D`
906 Bounding box for pixels.
908 WCS; used to convert pixel positions to sky coordinates
909 and vice-versa.
910 filterName : `str`
911 Name of filter. This can be used for flux limit comparisons.
912 photoCalib : `None`
913 Deprecated, only included for api compatibility.
914 epoch : `astropy.time.Time` or `None`, optional
915 Epoch to which to correct proper motion and parallax, or `None` to
916 not apply such corrections.
917
918 Returns
919 -------
920 results : `lsst.pipe.base.Struct`
921 A Struct containing the following fields:
922 refCat : `lsst.afw.catalog.SimpleCatalog`
923 A catalog of reference objects with the standard
924 schema, as documented in the main doc string for
925 `LoadReferenceObjects`.
926 The catalog is guaranteed to be contiguous.
927 fluxField : `str`
928 Name of flux field for specified `filterName`.
929
930 Notes
931 -----
932 The search algorithm works by searching in a region in sky
933 coordinates whose center is the center of the bbox and radius
934 is large enough to just include all 4 corners of the bbox.
935 Stars that lie outside the bbox are then trimmed from the list.
936 """
937 circle = self._calculateCircle_calculateCircle(bbox, wcs)
938
939 # find objects in circle
940 self.log.info("Loading reference objects from %s using center %s and radius %s deg",
941 self.config.ref_dataset_name, circle.coord, circle.radius.asDegrees())
942 loadRes = self.loadSkyCircleloadSkyCircle(circle.coord, circle.radius, filterName, epoch=epoch,
943 centroids=True)
944 refCat = loadRes.refCat
945 numFound = len(refCat)
946
947 # trim objects outside bbox
948 refCat = self._trimToBBox_trimToBBox(refCat=refCat, bbox=circle.bbox, wcs=wcs)
949 numTrimmed = numFound - len(refCat)
950 self.log.debug("trimmed %d out-of-bbox objects, leaving %d", numTrimmed, len(refCat))
951 self.log.info("Loaded %d reference objects", len(refCat))
952
953 # make sure catalog is contiguous
954 if not refCat.isContiguous():
955 loadRes.refCat = refCat.copy(deep=True)
956
957 return loadRes
958
959 @abc.abstractmethod
960 def loadSkyCircle(self, ctrCoord, radius, filterName, epoch=None, centroids=False):
961 """Load reference objects that overlap a circular sky region.
962
963 Parameters
964 ----------
965 ctrCoord : `lsst.geom.SpherePoint`
966 ICRS center of search region.
967 radius : `lsst.geom.Angle`
968 Radius of search region.
969 filterName : `str`
970 Name of filter. This can be used for flux limit comparisons.
971 epoch : `astropy.time.Time` or `None`, optional
972 Epoch to which to correct proper motion and parallax, or `None` to
973 not apply such corrections.
974 centroids : `bool`, optional
975 Add centroid fields to the loaded Schema. ``loadPixelBox`` expects
976 these fields to exist.
977
978 Returns
979 -------
980 results : `lsst.pipe.base.Struct`
981 A Struct containing the following fields:
982 refCat : `lsst.afw.catalog.SimpleCatalog`
983 A catalog of reference objects with the standard
984 schema, as documented in the main doc string for
985 `LoadReferenceObjects`.
986 The catalog is guaranteed to be contiguous.
987 fluxField : `str`
988 Name of flux field for specified `filterName`.
989
990 Notes
991 -----
992 Note that subclasses are responsible for performing the proper motion
993 correction, since this is the lowest-level interface for retrieving
994 the catalog.
995 """
996 return
997
998 @staticmethod
999 def _trimToBBox(refCat, bbox, wcs):
1000 """Remove objects outside a given pixel bounding box and set
1001 centroid and hasCentroid fields.
1002
1003 Parameters
1004 ----------
1006 A catalog of objects. The schema must include fields
1007 "coord", "centroid" and "hasCentroid".
1008 The "coord" field is read.
1009 The "centroid" and "hasCentroid" fields are set.
1010 bbox : `lsst.geom.Box2D`
1011 Pixel region
1012 wcs : `lsst.afw.geom.SkyWcs`
1013 WCS; used to convert sky coordinates to pixel positions.
1014
1015 Returns
1016 -------
1018 Reference objects in the bbox, with centroid and
1019 hasCentroid fields set.
1020 """
1021 afwTable.updateRefCentroids(wcs, refCat)
1022 centroidKey = afwTable.Point2DKey(refCat.schema["centroid"])
1023 retStarCat = type(refCat)(refCat.table)
1024 for star in refCat:
1025 point = star.get(centroidKey)
1026 if bbox.contains(point):
1027 retStarCat.append(star)
1028 return retStarCat
1029
1030 def _addFluxAliases(self, schema):
1031 """Add aliases for camera filter fluxes to the schema.
1032
1033 For each camFilter: refFilter in self.config.filterMap adds these aliases:
1034 <camFilter>_camFlux: <refFilter>_flux
1035 <camFilter>_camFluxErr: <refFilter>_fluxErr, if the latter exists
1036
1037 Parameters
1038 ----------
1039 schema : `lsst.afw.table.Schema`
1040 Schema for reference catalog.
1041
1042 Raises
1043 ------
1044 RuntimeError
1045 If any reference flux field is missing from the schema.
1046 """
1047 aliasMap = schema.getAliasMap()
1048
1049 if self.config.anyFilterMapsToThis is not None:
1050 refFluxName = self.config.anyFilterMapsToThis + "_flux"
1051 if refFluxName not in schema:
1052 msg = f"Unknown reference filter for anyFilterMapsToThis='{refFluxName}'"
1053 raise RuntimeError(msg)
1054 aliasMap.set("anyFilterMapsToThis", refFluxName)
1055 return # this is mutually exclusive with filterMap
1056
1057 def addAliasesForOneFilter(filterName, refFilterName):
1058 """Add aliases for a single filter
1059
1060 Parameters
1061 ----------
1062 filterName : `str` (optional)
1063 Camera filter name. The resulting alias name is
1064 <filterName>_camFlux
1065 refFilterName : `str`
1066 Reference catalog filter name; the field
1067 <refFilterName>_flux must exist.
1068 """
1069 camFluxName = filterName + "_camFlux"
1070 refFluxName = refFilterName + "_flux"
1071 if refFluxName not in schema:
1072 raise RuntimeError("Unknown reference filter %s" % (refFluxName,))
1073 aliasMap.set(camFluxName, refFluxName)
1074 refFluxErrName = refFluxName + "Err"
1075 if refFluxErrName in schema:
1076 camFluxErrName = camFluxName + "Err"
1077 aliasMap.set(camFluxErrName, refFluxErrName)
1078
1079 for filterName, refFilterName in self.config.filterMap.items():
1080 addAliasesForOneFilter(filterName, refFilterName)
1081
1082 @staticmethod
1083 def makeMinimalSchema(filterNameList, *, addCentroid=False,
1084 addIsPhotometric=False, addIsResolved=False,
1085 addIsVariable=False, coordErrDim=2,
1086 addProperMotion=False, properMotionErrDim=2,
1087 addParallax=False):
1088 """Make a standard schema for reference object catalogs.
1089
1090 Parameters
1091 ----------
1092 filterNameList : `list` of `str`
1093 List of filter names. Used to create <filterName>_flux fields.
1094 addIsPhotometric : `bool`
1095 If True then add field "photometric".
1096 addIsResolved : `bool`
1097 If True then add field "resolved".
1098 addIsVariable : `bool`
1099 If True then add field "variable".
1100 coordErrDim : `int`
1101 Number of coord error fields; must be one of 0, 2, 3:
1102
1103 - If 2 or 3: add fields "coord_raErr" and "coord_decErr".
1104 - If 3: also add field "coord_radecErr".
1105 addProperMotion : `bool`
1106 If True add fields "epoch", "pm_ra", "pm_dec" and "pm_flag".
1107 properMotionErrDim : `int`
1108 Number of proper motion error fields; must be one of 0, 2, 3;
1109 ignored if addProperMotion false:
1110 - If 2 or 3: add fields "pm_raErr" and "pm_decErr".
1111 - If 3: also add field "pm_radecErr".
1112 addParallax : `bool`
1113 If True add fields "epoch", "parallax", "parallaxErr"
1114 and "parallax_flag".
1115
1116 Returns
1117 -------
1118 schema : `lsst.afw.table.Schema`
1119 Schema for reference catalog, an
1121
1122 Notes
1123 -----
1124 Reference catalogs support additional covariances, such as
1125 covariance between RA and proper motion in declination,
1126 that are not supported by this method, but can be added after
1127 calling this method.
1128 """
1129 schema = afwTable.SimpleTable.makeMinimalSchema()
1130 if addCentroid:
1131 afwTable.Point2DKey.addFields(
1132 schema,
1133 "centroid",
1134 "centroid on an exposure, if relevant",
1135 "pixel",
1136 )
1137 schema.addField(
1138 field="hasCentroid",
1139 type="Flag",
1140 doc="is position known?",
1141 )
1142 for filterName in filterNameList:
1143 schema.addField(
1144 field="%s_flux" % (filterName,),
1145 type=numpy.float64,
1146 doc="flux in filter %s" % (filterName,),
1147 units="nJy",
1148 )
1149 for filterName in filterNameList:
1150 schema.addField(
1151 field="%s_fluxErr" % (filterName,),
1152 type=numpy.float64,
1153 doc="flux uncertainty in filter %s" % (filterName,),
1154 units="nJy",
1155 )
1156 if addIsPhotometric:
1157 schema.addField(
1158 field="photometric",
1159 type="Flag",
1160 doc="set if the object can be used for photometric calibration",
1161 )
1162 if addIsResolved:
1163 schema.addField(
1164 field="resolved",
1165 type="Flag",
1166 doc="set if the object is spatially resolved",
1167 )
1168 if addIsVariable:
1169 schema.addField(
1170 field="variable",
1171 type="Flag",
1172 doc="set if the object has variable brightness",
1173 )
1174 if coordErrDim not in (0, 2, 3):
1175 raise ValueError("coordErrDim={}; must be (0, 2, 3)".format(coordErrDim))
1176 if coordErrDim > 0:
1177 afwTable.CovarianceMatrix2fKey.addFields(
1178 schema=schema,
1179 prefix="coord",
1180 names=["ra", "dec"],
1181 units=["rad", "rad"],
1182 diagonalOnly=(coordErrDim == 2),
1183 )
1184
1185 if addProperMotion or addParallax:
1186 schema.addField(
1187 field="epoch",
1188 type=numpy.float64,
1189 doc="date of observation (TAI, MJD)",
1190 units="day",
1191 )
1192
1193 if addProperMotion:
1194 schema.addField(
1195 field="pm_ra",
1196 type="Angle",
1197 doc="proper motion in the right ascension direction = dra/dt * cos(dec)",
1198 units="rad/year",
1199 )
1200 schema.addField(
1201 field="pm_dec",
1202 type="Angle",
1203 doc="proper motion in the declination direction",
1204 units="rad/year",
1205 )
1206 if properMotionErrDim not in (0, 2, 3):
1207 raise ValueError("properMotionErrDim={}; must be (0, 2, 3)".format(properMotionErrDim))
1208 if properMotionErrDim > 0:
1209 afwTable.CovarianceMatrix2fKey.addFields(
1210 schema=schema,
1211 prefix="pm",
1212 names=["ra", "dec"],
1213 units=["rad/year", "rad/year"],
1214 diagonalOnly=(properMotionErrDim == 2),
1215 )
1216 schema.addField(
1217 field="pm_flag",
1218 type="Flag",
1219 doc="Set if proper motion or proper motion error is bad",
1220 )
1221
1222 if addParallax:
1223 schema.addField(
1224 field="parallax",
1225 type="Angle",
1226 doc="parallax",
1227 units="rad",
1228 )
1229 schema.addField(
1230 field="parallaxErr",
1231 type="Angle",
1232 doc="uncertainty in parallax",
1233 units="rad",
1234 )
1235 schema.addField(
1236 field="parallax_flag",
1237 type="Flag",
1238 doc="Set if parallax or parallax error is bad",
1239 )
1240 return schema
1241
1242 def _calculateCircle(self, bbox, wcs):
1243 """Compute on-sky center and radius of search region.
1244
1245 Parameters
1246 ----------
1247 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D`
1248 Pixel bounding box.
1249 wcs : `lsst.afw.geom.SkyWcs`
1250 WCS; used to convert pixel positions to sky coordinates.
1251
1252 Returns
1253 -------
1254 results : `lsst.pipe.base.Struct`
1255 A Struct containing:
1256
1257 - coord : `lsst.geom.SpherePoint`
1258 ICRS center of the search region.
1259 - radius : `lsst.geom.Angle`
1260 Radius of the search region.
1261 - bbox : `lsst.geom.Box2D`
1262 Bounding box used to compute the circle.
1263 """
1264 bbox = geom.Box2D(bbox) # make sure bbox is double and that we have a copy
1265 bbox.grow(self.config.pixelMargin)
1266 coord = wcs.pixelToSky(bbox.getCenter())
1267 radius = max(coord.separation(wcs.pixelToSky(pp)) for pp in bbox.getCorners())
1268 return pipeBase.Struct(coord=coord, radius=radius, bbox=bbox)
1269
1270 def getMetadataBox(self, bbox, wcs, filterName, photoCalib=None, epoch=None):
1271 """Return metadata about the load.
1272
1273 This metadata is used for reloading the catalog (e.g., for
1274 reconstituting a normalised match list.
1275
1276 Parameters
1277 ----------
1278 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D`
1279 Pixel bounding box.
1280 wcs : `lsst.afw.geom.SkyWcs`
1281 WCS; used to convert pixel positions to sky coordinates.
1282 filterName : `str`
1283 Name of camera filter.
1284 photoCalib : `None`
1285 Deprecated, only included for api compatibility.
1286 epoch : `astropy.time.Time` or `None`, optional
1287 Epoch to which to correct proper motion and parallax,
1288 or None to not apply such corrections.
1289
1290 Returns
1291 -------
1292 metadata : `lsst.daf.base.PropertyList`
1293 Metadata about the load.
1294 """
1295 circle = self._calculateCircle_calculateCircle(bbox, wcs)
1296 return self.getMetadataCirclegetMetadataCircle(circle.coord, circle.radius, filterName, epoch=epoch)
1297
1298 def getMetadataCircle(self, coord, radius, filterName, photoCalib=None, epoch=None):
1299 """Return metadata about the load.
1300
1301 This metadata is used for reloading the catalog (e.g., for
1302 reconstituting a normalised match list.
1303
1304 Parameters
1305 ----------
1306 coord : `lsst.geom.SpherePoint`
1307 ICRS center of the search region.
1308 radius : `lsst.geom.Angle`
1309 Radius of the search region.
1310 filterName : `str`
1311 Name of camera filter.
1312 photoCalib : `None`
1313 Deprecated, only included for api compatibility.
1314 epoch : `astropy.time.Time` (optional)
1315 Epoch to which to correct proper motion and parallax, or `None` to
1316 not apply such corrections.
1317
1318 Returns
1319 -------
1321 Metadata about the load
1322 """
1323 md = PropertyList()
1324 md.add('RA', coord.getRa().asDegrees(), 'field center in degrees')
1325 md.add('DEC', coord.getDec().asDegrees(), 'field center in degrees')
1326 md.add('RADIUS', radius.asDegrees(), 'field radius in degrees, minimum')
1327 md.add('SMATCHV', 1, 'SourceMatchVector version number')
1328 md.add('FILTER', filterName, 'filter name for photometric data')
1329 md.add('EPOCH', "NONE" if epoch is None else epoch.mjd, 'Epoch (TAI MJD) for catalog')
1330 return md
1331
1332 def joinMatchListWithCatalog(self, matchCat, sourceCat):
1333 """Relink an unpersisted match list to sources and reference
1334 objects.
1335
1336 A match list is persisted and unpersisted as a catalog of IDs
1337 produced by afw.table.packMatches(), with match metadata
1338 (as returned by the astrometry tasks) in the catalog's metadata
1339 attribute. This method converts such a match catalog into a match
1340 list, with links to source records and reference object records.
1341
1342 Parameters
1343 ----------
1344 matchCat : `lsst.afw.table.BaseCatalog`
1345 Unperisted packed match list.
1346 ``matchCat.table.getMetadata()`` must contain match metadata,
1347 as returned by the astrometry tasks.
1348 sourceCat : `lsst.afw.table.SourceCatalog`
1349 Source catalog. As a side effect, the catalog will be sorted
1350 by ID.
1351
1352 Returns
1353 -------
1355 Match list.
1356 """
1357 return joinMatchListWithCatalogImpl(self, matchCat, sourceCat)
1358
1359
1360def joinMatchListWithCatalogImpl(refObjLoader, matchCat, sourceCat):
1361 """Relink an unpersisted match list to sources and reference
1362 objects.
1363
1364 A match list is persisted and unpersisted as a catalog of IDs
1365 produced by afw.table.packMatches(), with match metadata
1366 (as returned by the astrometry tasks) in the catalog's metadata
1367 attribute. This method converts such a match catalog into a match
1368 list, with links to source records and reference object records.
1369
1370 Parameters
1371 ----------
1372 refObjLoader
1373 Reference object loader to use in getting reference objects
1374 matchCat : `lsst.afw.table.BaseCatalog`
1375 Unperisted packed match list.
1376 ``matchCat.table.getMetadata()`` must contain match metadata,
1377 as returned by the astrometry tasks.
1378 sourceCat : `lsst.afw.table.SourceCatalog`
1379 Source catalog. As a side effect, the catalog will be sorted
1380 by ID.
1381
1382 Returns
1383 -------
1385 Match list.
1386 """
1387 matchmeta = matchCat.table.getMetadata()
1388 version = matchmeta.getInt('SMATCHV')
1389 if version != 1:
1390 raise ValueError('SourceMatchVector version number is %i, not 1.' % version)
1391 filterName = matchmeta.getString('FILTER').strip()
1392 try:
1393 epoch = matchmeta.getDouble('EPOCH')
1394 except (LookupError, TypeError):
1395 epoch = None # Not present, or not correct type means it's not set
1396 if 'RADIUS' in matchmeta:
1397 # This is a circle style metadata, call loadSkyCircle
1398 ctrCoord = geom.SpherePoint(matchmeta.getDouble('RA'),
1399 matchmeta.getDouble('DEC'), geom.degrees)
1400 rad = matchmeta.getDouble('RADIUS')*geom.degrees
1401 refCat = refObjLoader.loadSkyCircle(ctrCoord, rad, filterName, epoch=epoch).refCat
1402 elif "INNER_UPPER_LEFT_RA" in matchmeta:
1403 # This is the sky box type (only triggers in the LoadReferenceObject class, not task)
1404 # Only the outer box is required to be loaded to get the maximum region, all filtering
1405 # will be done by the unpackMatches function, and no spatial filtering needs to be done
1406 # by the refObjLoader
1407 box = []
1408 for place in ("UPPER_LEFT", "UPPER_RIGHT", "LOWER_LEFT", "LOWER_RIGHT"):
1409 coord = geom.SpherePoint(matchmeta.getDouble(f"OUTER_{place}_RA"),
1410 matchmeta.getDouble(f"OUTER_{place}_DEC"),
1411 geom.degrees).getVector()
1412 box.append(coord)
1413 outerBox = sphgeom.ConvexPolygon(box)
1414 refCat = refObjLoader.loadRegion(outerBox, filterName, epoch=epoch).refCat
1415
1416 refCat.sort()
1417 sourceCat.sort()
1418 return afwTable.unpackMatches(matchCat, refCat, sourceCat)
1419
1420
1421def applyProperMotionsImpl(log, catalog, epoch):
1422 """Apply proper motion correction to a reference catalog.
1423
1424 Adjust position and position error in the ``catalog``
1425 for proper motion to the specified ``epoch``,
1426 modifying the catalog in place.
1427
1428 Parameters
1429 ----------
1430 log : `lsst.log.Log` or `logging.getLogger`
1431 Log object to write to.
1433 Catalog of positions, containing:
1434
1435 - Coordinates, retrieved by the table's coordinate key.
1436 - ``coord_raErr`` : Error in Right Ascension (rad).
1437 - ``coord_decErr`` : Error in Declination (rad).
1438 - ``pm_ra`` : Proper motion in Right Ascension (rad/yr,
1439 East positive)
1440 - ``pm_raErr`` : Error in ``pm_ra`` (rad/yr), optional.
1441 - ``pm_dec`` : Proper motion in Declination (rad/yr,
1442 North positive)
1443 - ``pm_decErr`` : Error in ``pm_dec`` (rad/yr), optional.
1444 - ``epoch`` : Mean epoch of object (an astropy.time.Time)
1445 epoch : `astropy.time.Time`
1446 Epoch to which to correct proper motion.
1447 """
1448 if "epoch" not in catalog.schema or "pm_ra" not in catalog.schema or "pm_dec" not in catalog.schema:
1449 log.warning("Proper motion correction not available from catalog")
1450 return
1451 if not catalog.isContiguous():
1452 raise RuntimeError("Catalog must be contiguous")
1453 catEpoch = astropy.time.Time(catalog["epoch"], scale="tai", format="mjd")
1454 log.info("Correcting reference catalog for proper motion to %r", epoch)
1455 # Use `epoch.tai` to make sure the time difference is in TAI
1456 timeDiffsYears = (epoch.tai - catEpoch).to(astropy.units.yr).value
1457 coordKey = catalog.table.getCoordKey()
1458 # Compute the offset of each object due to proper motion
1459 # as components of the arc of a great circle along RA and Dec
1460 pmRaRad = catalog["pm_ra"]
1461 pmDecRad = catalog["pm_dec"]
1462 offsetsRaRad = pmRaRad*timeDiffsYears
1463 offsetsDecRad = pmDecRad*timeDiffsYears
1464 # Compute the corresponding bearing and arc length of each offset
1465 # due to proper motion, and apply the offset
1466 # The factor of 1e6 for computing bearing is intended as
1467 # a reasonable scale for typical values of proper motion
1468 # in order to avoid large errors for small values of proper motion;
1469 # using the offsets is another option, but it can give
1470 # needlessly large errors for short duration
1471 offsetBearingsRad = numpy.arctan2(pmDecRad*1e6, pmRaRad*1e6)
1472 offsetAmountsRad = numpy.hypot(offsetsRaRad, offsetsDecRad)
1473 for record, bearingRad, amountRad in zip(catalog, offsetBearingsRad, offsetAmountsRad):
1474 record.set(coordKey,
1475 record.get(coordKey).offset(bearing=bearingRad*geom.radians,
1476 amount=amountRad*geom.radians))
1477 # Increase error in RA and Dec based on error in proper motion
1478 if "coord_raErr" in catalog.schema:
1479 catalog["coord_raErr"] = numpy.hypot(catalog["coord_raErr"],
1480 catalog["pm_raErr"]*timeDiffsYears)
1481 if "coord_decErr" in catalog.schema:
1482 catalog["coord_decErr"] = numpy.hypot(catalog["coord_decErr"],
1483 catalog["pm_decErr"]*timeDiffsYears)
def loadSkyCircle(self, ctrCoord, radius, filterName, epoch=None, centroids=False)
def getMetadataBox(self, bbox, wcs, filterName, photoCalib=None, epoch=None)
def loadPixelBox(self, bbox, wcs, filterName, photoCalib=None, epoch=None)
def makeMinimalSchema(filterNameList, *addCentroid=False, addIsPhotometric=False, addIsResolved=False, addIsVariable=False, coordErrDim=2, addProperMotion=False, properMotionErrDim=2, addParallax=False)
def getMetadataCircle(self, coord, radius, filterName, photoCalib=None, epoch=None)
def getMetadataCircle(coord, radius, filterName, photoCalib=None, epoch=None)
def remapReferenceCatalogSchema(refCat, *filterNameList=None, position=False, photometric=False)
def loadSkyCircle(self, ctrCoord, radius, filterName, epoch=None)
def loadRegion(self, region, filterName, filtFunc=None, epoch=None)
def loadPixelBox(self, bbox, wcs, filterName, epoch=None, photoCalib=None, bboxToSpherePadding=100)
def getMetadataBox(self, bbox, wcs, filterName, photoCalib=None, epoch=None, bboxToSpherePadding=100)
def convertToNanojansky(catalog, log, doConvert=True)
def joinMatchListWithCatalogImpl(refObjLoader, matchCat, sourceCat)