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