LSST Applications g0b6bd0c080+a72a5dd7e6,g1182afd7b4+2a019aa3bb,g17e5ecfddb+2b8207f7de,g1d67935e3f+06cf436103,g38293774b4+ac198e9f13,g396055baef+6a2097e274,g3b44f30a73+6611e0205b,g480783c3b1+98f8679e14,g48ccf36440+89c08d0516,g4b93dc025c+98f8679e14,g5c4744a4d9+a302e8c7f0,g613e996a0d+e1c447f2e0,g6c8d09e9e7+25247a063c,g7271f0639c+98f8679e14,g7a9cd813b8+124095ede6,g9d27549199+a302e8c7f0,ga1cf026fa3+ac198e9f13,ga32aa97882+7403ac30ac,ga786bb30fb+7a139211af,gaa63f70f4e+9994eb9896,gabf319e997+ade567573c,gba47b54d5d+94dc90c3ea,gbec6a3398f+06cf436103,gc6308e37c7+07dd123edb,gc655b1545f+ade567573c,gcc9029db3c+ab229f5caf,gd01420fc67+06cf436103,gd877ba84e5+06cf436103,gdb4cecd868+6f279b5b48,ge2d134c3d5+cc4dbb2e3f,ge448b5faa6+86d1ceac1d,gecc7e12556+98f8679e14,gf3ee170dca+25247a063c,gf4ac96e456+ade567573c,gf9f5ea5b4d+ac198e9f13,gff490e6085+8c2580be5c,w.2022.27
LSST Data Management Base Package
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", "ReferenceObjectLoaderBase"]
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.reserve(len(catalog))
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
147
148
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.
153
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.regionregion = region
161
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.
165
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.
169
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.
174
175 Parameters
176 ---------
178 SourceCatalog to be filtered.
179 catRegion : `lsst.sphgeom.Region`
180 Region in which the catalog was created
181 """
182 if catRegion.isWithin(self.regionregion):
183 # no filtering needed, region completely contains refcat
184 return refCat
185
186 filteredRefCat = type(refCat)(refCat.table)
187 for record in refCat:
188 if self.regionregion.contains(record.getCoord().getVector()):
189 filteredRefCat.append(record)
190 return filteredRefCat
191
192
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 anyFilterMapsToThis = pexConfig.Field(
201 doc=("Always use this reference catalog filter, no matter whether or what filter name is "
202 "supplied to the loader. Effectively a trivial filterMap: map all filter names to this filter."
203 " This can be set for purely-astrometric catalogs (e.g. Gaia DR2) where there is only one "
204 "reasonable choice for every camera filter->refcat mapping, but not for refcats used for "
205 "photometry, which need a filterMap and/or colorterms/transmission corrections."),
206 dtype=str,
207 default=None,
208 optional=True
209 )
210 filterMap = pexConfig.DictField(
211 doc=("Mapping of camera filter name: reference catalog filter name; "
212 "each reference filter must exist in the refcat."
213 " Note that this does not perform any bandpass corrections: it is just a lookup."),
214 keytype=str,
215 itemtype=str,
216 default={},
217 )
218 requireProperMotion = pexConfig.Field(
219 doc="Require that the fields needed to correct proper motion "
220 "(epoch, pm_ra and pm_dec) are present?",
221 dtype=bool,
222 default=False,
223 )
224
225 def validate(self):
226 super().validate()
227 if self.filterMapfilterMap != {} and self.anyFilterMapsToThisanyFilterMapsToThis is not None:
228 msg = "`filterMap` and `anyFilterMapsToThis` are mutually exclusive"
229 raise pexConfig.FieldValidationError(LoadReferenceObjectsConfig.anyFilterMapsToThis,
230 self, msg)
231
232
234 """Base class for reference object loaders, to facilitate gen2/gen3 code
235 sharing.
236
237 Parameters
238 ----------
239 config : `lsst.pex.config.Config`
240 Configuration for the loader.
241 """
242 ConfigClass = LoadReferenceObjectsConfig
243
244 def __init__(self, config=None, *args, **kwargs):
245 self.configconfig = config
246
247 def applyProperMotions(self, catalog, epoch):
248 """Apply proper motion correction to a reference catalog.
249
250 Adjust position and position error in the ``catalog``
251 for proper motion to the specified ``epoch``,
252 modifying the catalog in place.
253
254 Parameters
255 ----------
257 Catalog of positions, containing at least these fields:
258
259 - Coordinates, retrieved by the table's coordinate key.
260 - ``coord_raErr`` : Error in Right Ascension (rad).
261 - ``coord_decErr`` : Error in Declination (rad).
262 - ``pm_ra`` : Proper motion in Right Ascension (rad/yr,
263 East positive)
264 - ``pm_raErr`` : Error in ``pm_ra`` (rad/yr), optional.
265 - ``pm_dec`` : Proper motion in Declination (rad/yr,
266 North positive)
267 - ``pm_decErr`` : Error in ``pm_dec`` (rad/yr), optional.
268 - ``epoch`` : Mean epoch of object (an astropy.time.Time)
269 epoch : `astropy.time.Time`
270 Epoch to which to correct proper motion.
271 If None, do not apply PM corrections or raise if
272 ``config.requireProperMotion`` is True.
273
274 Raises
275 ------
276 RuntimeError
277 Raised if ``config.requireProperMotion`` is set but we cannot
278 apply the proper motion correction for some reason.
279 """
280 if epoch is None:
281 if self.configconfig.requireProperMotion:
282 raise RuntimeError("requireProperMotion=True but epoch not provided to loader.")
283 else:
284 self.log.debug("No epoch provided: not applying proper motion corrections to refcat.")
285 return
286
287 # Warn/raise for a catalog in an incorrect format, if epoch was specified.
288 if ("pm_ra" in catalog.schema
289 and not isinstance(catalog.schema["pm_ra"].asKey(), afwTable.KeyAngle)):
290 if self.configconfig.requireProperMotion:
291 raise RuntimeError("requireProperMotion=True but refcat pm_ra field is not an Angle.")
292 else:
293 self.log.warning("Reference catalog pm_ra field is not an Angle; cannot apply proper motion.")
294 return
295
296 if ("epoch" not in catalog.schema or "pm_ra" not in catalog.schema):
297 if self.configconfig.requireProperMotion:
298 raise RuntimeError("requireProperMotion=True but PM data not available from catalog.")
299 else:
300 self.log.warning("Proper motion correction not available for this reference catalog.")
301 return
302
303 applyProperMotionsImpl(self.log, catalog, epoch)
304
305 @staticmethod
306 def makeMinimalSchema(filterNameList, *, addCentroid=False,
307 addIsPhotometric=False, addIsResolved=False,
308 addIsVariable=False, coordErrDim=2,
309 addProperMotion=False, properMotionErrDim=2,
310 addParallax=False):
311 """Make a standard schema for reference object catalogs.
312
313 Parameters
314 ----------
315 filterNameList : `list` of `str`
316 List of filter names. Used to create <filterName>_flux fields.
317 addIsPhotometric : `bool`
318 If True then add field "photometric".
319 addIsResolved : `bool`
320 If True then add field "resolved".
321 addIsVariable : `bool`
322 If True then add field "variable".
323 coordErrDim : `int`
324 Number of coord error fields; must be one of 0, 2, 3:
325
326 - If 2 or 3: add fields "coord_raErr" and "coord_decErr".
327 - If 3: also add field "coord_radecErr".
328 addProperMotion : `bool`
329 If True add fields "epoch", "pm_ra", "pm_dec" and "pm_flag".
330 properMotionErrDim : `int`
331 Number of proper motion error fields; must be one of 0, 2, 3;
332 ignored if addProperMotion false:
333 - If 2 or 3: add fields "pm_raErr" and "pm_decErr".
334 - If 3: also add field "pm_radecErr".
335 addParallax : `bool`
336 If True add fields "epoch", "parallax", "parallaxErr"
337 and "parallax_flag".
338
339 Returns
340 -------
341 schema : `lsst.afw.table.Schema`
342 Schema for reference catalog, an
344
345 Notes
346 -----
347 Reference catalogs support additional covariances, such as
348 covariance between RA and proper motion in declination,
349 that are not supported by this method, but can be added after
350 calling this method.
351 """
352 schema = afwTable.SimpleTable.makeMinimalSchema()
353 if addCentroid:
354 afwTable.Point2DKey.addFields(
355 schema,
356 "centroid",
357 "centroid on an exposure, if relevant",
358 "pixel",
359 )
360 schema.addField(
361 field="hasCentroid",
362 type="Flag",
363 doc="is position known?",
364 )
365 for filterName in filterNameList:
366 schema.addField(
367 field="%s_flux" % (filterName,),
368 type=numpy.float64,
369 doc="flux in filter %s" % (filterName,),
370 units="nJy",
371 )
372 for filterName in filterNameList:
373 schema.addField(
374 field="%s_fluxErr" % (filterName,),
375 type=numpy.float64,
376 doc="flux uncertainty in filter %s" % (filterName,),
377 units="nJy",
378 )
379 if addIsPhotometric:
380 schema.addField(
381 field="photometric",
382 type="Flag",
383 doc="set if the object can be used for photometric calibration",
384 )
385 if addIsResolved:
386 schema.addField(
387 field="resolved",
388 type="Flag",
389 doc="set if the object is spatially resolved",
390 )
391 if addIsVariable:
392 schema.addField(
393 field="variable",
394 type="Flag",
395 doc="set if the object has variable brightness",
396 )
397 if coordErrDim not in (0, 2, 3):
398 raise ValueError("coordErrDim={}; must be (0, 2, 3)".format(coordErrDim))
399 if coordErrDim > 0:
400 afwTable.CovarianceMatrix2fKey.addFields(
401 schema=schema,
402 prefix="coord",
403 names=["ra", "dec"],
404 units=["rad", "rad"],
405 diagonalOnly=(coordErrDim == 2),
406 )
407
408 if addProperMotion or addParallax:
409 schema.addField(
410 field="epoch",
411 type=numpy.float64,
412 doc="date of observation (TAI, MJD)",
413 units="day",
414 )
415
416 if addProperMotion:
417 schema.addField(
418 field="pm_ra",
419 type="Angle",
420 doc="proper motion in the right ascension direction = dra/dt * cos(dec)",
421 units="rad/year",
422 )
423 schema.addField(
424 field="pm_dec",
425 type="Angle",
426 doc="proper motion in the declination direction",
427 units="rad/year",
428 )
429 if properMotionErrDim not in (0, 2, 3):
430 raise ValueError("properMotionErrDim={}; must be (0, 2, 3)".format(properMotionErrDim))
431 if properMotionErrDim > 0:
432 afwTable.CovarianceMatrix2fKey.addFields(
433 schema=schema,
434 prefix="pm",
435 names=["ra", "dec"],
436 units=["rad/year", "rad/year"],
437 diagonalOnly=(properMotionErrDim == 2),
438 )
439 schema.addField(
440 field="pm_flag",
441 type="Flag",
442 doc="Set if proper motion or proper motion error is bad",
443 )
444
445 if addParallax:
446 schema.addField(
447 field="parallax",
448 type="Angle",
449 doc="parallax",
450 units="rad",
451 )
452 schema.addField(
453 field="parallaxErr",
454 type="Angle",
455 doc="uncertainty in parallax",
456 units="rad",
457 )
458 schema.addField(
459 field="parallax_flag",
460 type="Flag",
461 doc="Set if parallax or parallax error is bad",
462 )
463 return schema
464
465 @staticmethod
466 def _remapReferenceCatalogSchema(refCat, *, anyFilterMapsToThis=None,
467 filterMap=None, centroids=False):
468 """This function takes in a reference catalog and returns a new catalog
469 with additional columns defined from the remaining function arguments.
470
471 Parameters
472 ----------
474 Reference catalog to map to new catalog
475 anyFilterMapsToThis : `str`, optional
476 Always use this reference catalog filter.
477 Mutually exclusive with `filterMap`
478 filterMap : `dict` [`str`,`str`], optional
479 Mapping of camera filter name: reference catalog filter name.
480 centroids : `bool`, optional
481 Add centroid fields to the loaded Schema. ``loadPixelBox`` expects
482 these fields to exist.
483
484 Returns
485 -------
486 expandedCat : `lsst.afw.table.SimpleCatalog`
487 Deep copy of input reference catalog with additional columns added
488 """
489 if anyFilterMapsToThis or filterMap:
490 ReferenceObjectLoaderBase._addFluxAliases(refCat.schema, anyFilterMapsToThis, filterMap)
491
492 mapper = afwTable.SchemaMapper(refCat.schema, True)
493 mapper.addMinimalSchema(refCat.schema, True)
494 mapper.editOutputSchema().disconnectAliases()
495
496 if centroids:
497 # Add and initialize centroid and hasCentroid fields (these are
498 # added after loading to avoid wasting space in the saved catalogs).
499 # The new fields are automatically initialized to (nan, nan) and
500 # False so no need to set them explicitly.
501 mapper.editOutputSchema().addField("centroid_x", type=float, doReplace=True)
502 mapper.editOutputSchema().addField("centroid_y", type=float, doReplace=True)
503 mapper.editOutputSchema().addField("hasCentroid", type="Flag", doReplace=True)
504 mapper.editOutputSchema().getAliasMap().set("slot_Centroid", "centroid")
505
506 expandedCat = afwTable.SimpleCatalog(mapper.getOutputSchema())
507 expandedCat.setMetadata(refCat.getMetadata())
508 expandedCat.extend(refCat, mapper=mapper)
509
510 return expandedCat
511
512 @staticmethod
513 def _addFluxAliases(schema, anyFilterMapsToThis=None, filterMap=None):
514 """Add aliases for camera filter fluxes to the schema.
515
516 For each camFilter: refFilter in filterMap, adds these aliases:
517 <camFilter>_camFlux: <refFilter>_flux
518 <camFilter>_camFluxErr: <refFilter>_fluxErr, if the latter exists
519 or sets `anyFilterMapsToThis` in the schema.
520
521 Parameters
522 ----------
523 schema : `lsst.afw.table.Schema`
524 Schema for reference catalog.
525 anyFilterMapsToThis : `str`, optional
526 Always use this reference catalog filter.
527 Mutually exclusive with `filterMap`.
528 filterMap : `dict` [`str`,`str`], optional
529 Mapping of camera filter name: reference catalog filter name.
530 Mutually exclusive with `anyFilterMapsToThis`.
531
532 Raises
533 ------
534 RuntimeError
535 Raised if any required reference flux field is missing from the
536 schema.
537 """
538 # Fail on any truthy value for either of these.
539 if anyFilterMapsToThis and filterMap:
540 raise ValueError("anyFilterMapsToThis and filterMap are mutually exclusive!")
541
542 aliasMap = schema.getAliasMap()
543
544 if anyFilterMapsToThis is not None:
545 refFluxName = anyFilterMapsToThis + "_flux"
546 if refFluxName not in schema:
547 msg = f"Unknown reference filter for anyFilterMapsToThis='{refFluxName}'"
548 raise RuntimeError(msg)
549 aliasMap.set("anyFilterMapsToThis", refFluxName)
550 return # this is mutually exclusive with filterMap
551
552 def addAliasesForOneFilter(filterName, refFilterName):
553 """Add aliases for a single filter
554
555 Parameters
556 ----------
557 filterName : `str` (optional)
558 Camera filter name. The resulting alias name is
559 <filterName>_camFlux
560 refFilterName : `str`
561 Reference catalog filter name; the field
562 <refFilterName>_flux must exist.
563 """
564 camFluxName = filterName + "_camFlux"
565 refFluxName = refFilterName + "_flux"
566 if refFluxName not in schema:
567 raise RuntimeError("Unknown reference filter %s" % (refFluxName,))
568 aliasMap.set(camFluxName, refFluxName)
569 refFluxErrName = refFluxName + "Err"
570 if refFluxErrName in schema:
571 camFluxErrName = camFluxName + "Err"
572 aliasMap.set(camFluxErrName, refFluxErrName)
573
574 if filterMap is not None:
575 for filterName, refFilterName in filterMap.items():
576 addAliasesForOneFilter(filterName, refFilterName)
577
578 @staticmethod
579 def _makeBoxRegion(BBox, wcs, BBoxPadding):
580 outerLocalBBox = geom.Box2D(BBox)
581 innerLocalBBox = geom.Box2D(BBox)
582
583 # Grow the bounding box to allow for effects not fully captured by the
584 # wcs provided (which represents the current best-guess wcs solution
585 # associated with the dataset for which the calibration is to be
586 # computed using the loaded and trimmed reference catalog being defined
587 # here). These effects could include pointing errors and/or an
588 # insufficient optical distorition model for the instrument. The idea
589 # is to ensure the spherical geometric region created contains the
590 # entire region covered by the bbox.
591 # Also create an inner region that is sure to be inside the bbox.
592 outerLocalBBox.grow(BBoxPadding)
593 innerLocalBBox.grow(-1*BBoxPadding)
594
595 # Handle the case where the inner bounding box shrank to a zero sized
596 # region (which will be the case if the shrunken size of either
597 # dimension is less than or equal to zero). In this case, the inner
598 # bounding box is set to the original input bounding box. This is
599 # probably not the best way to handle an empty inner bounding box, but
600 # it is what the calling code currently expects.
601 if innerLocalBBox.getDimensions() == geom.Extent2D(0, 0):
602 innerLocalBBox = geom.Box2D(BBox)
603
604 # Convert the corners of the bounding boxes to sky coordinates.
605 innerBoxCorners = innerLocalBBox.getCorners()
606 innerSphCorners = [wcs.pixelToSky(corner).getVector() for corner in innerBoxCorners]
607 innerSkyRegion = sphgeom.ConvexPolygon(innerSphCorners)
608
609 outerBoxCorners = outerLocalBBox.getCorners()
610 outerSphCorners = [wcs.pixelToSky(corner).getVector() for corner in outerBoxCorners]
611 outerSkyRegion = sphgeom.ConvexPolygon(outerSphCorners)
612
613 return innerSkyRegion, outerSkyRegion, innerSphCorners, outerSphCorners
614
615 @staticmethod
616 def _calculateCircle(bbox, wcs, pixelMargin):
617 """Compute on-sky center and radius of search region.
618
619 Parameters
620 ----------
621 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D`
622 Pixel bounding box.
624 WCS; used to convert pixel positions to sky coordinates.
625 pixelMargin : `int`
626 Padding to add to 4 all edges of the bounding box (pixels).
627
628 Returns
629 -------
630 results : `lsst.pipe.base.Struct`
631 A Struct containing:
632
633 - coord : `lsst.geom.SpherePoint`
634 ICRS center of the search region.
635 - radius : `lsst.geom.Angle`
636 Radius of the search region.
637 - bbox : `lsst.geom.Box2D`
638 Bounding box used to compute the circle.
639 """
640 bbox = geom.Box2D(bbox) # we modify the box, so use a copy
641 bbox.grow(pixelMargin)
642 coord = wcs.pixelToSky(bbox.getCenter())
643 radius = max(coord.separation(wcs.pixelToSky(pp)) for pp in bbox.getCorners())
644 return pipeBase.Struct(coord=coord, radius=radius, bbox=bbox)
645
646 @staticmethod
647 def getMetadataCircle(coord, radius, filterName, epoch=None):
648 """Return metadata about the reference catalog being loaded.
649
650 This metadata is used for reloading the catalog (e.g. for
651 reconstituting a normalized match list).
652
653 Parameters
654 ----------
655 coord : `lsst.geom.SpherePoint`
656 ICRS center of the search region.
657 radius : `lsst.geom.Angle`
658 Radius of the search region.
659 filterName : `str`
660 Name of the camera filter.
661 epoch : `astropy.time.Time` or `None`, optional
662 Epoch to which to correct proper motion and parallax, or `None` to
663 not apply such corrections.
664
665 Returns
666 -------
668 Metadata about the catalog.
669 """
670 md = PropertyList()
671 md.add('RA', coord.getRa().asDegrees(), 'field center in degrees')
672 md.add('DEC', coord.getDec().asDegrees(), 'field center in degrees')
673 md.add('RADIUS', radius.asDegrees(), 'field radius in degrees, minimum')
674 md.add('SMATCHV', 1, 'SourceMatchVector version number')
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
678
679 def getMetadataBox(self, bbox, wcs, filterName, epoch=None,
680 bboxToSpherePadding=100):
681 """Return metadata about the load
682
683 This metadata is used for reloading the catalog (e.g., for
684 reconstituting a normalised match list).
685
686 Parameters
687 ----------
688 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D`
689 Bounding box for the pixels.
691 The WCS object associated with ``bbox``.
692 filterName : `str`
693 Name of the camera filter.
694 epoch : `astropy.time.Time` or `None`, optional
695 Epoch to which to correct proper motion and parallax, or `None` to
696 not apply such corrections.
697 bboxToSpherePadding : `int`, optional
698 Padding in pixels to account for translating a set of corners into
699 a spherical (convex) boundary that is certain to encompass the
700 enitre area covered by the box.
701
702 Returns
703 -------
705 The metadata detailing the search parameters used for this
706 dataset.
707 """
708 circle = self._calculateCircle_calculateCircle(bbox, wcs, self.configconfig.pixelMargin)
709 md = self.getMetadataCirclegetMetadataCircle(circle.coord, circle.radius, filterName, epoch=epoch)
710
711 paddedBbox = circle.bbox
712 _, _, innerCorners, outerCorners = self._makeBoxRegion_makeBoxRegion(paddedBbox, wcs, bboxToSpherePadding)
713 for box, corners in zip(("INNER", "OUTER"), (innerCorners, outerCorners)):
714 for (name, corner) in zip(("UPPER_LEFT", "UPPER_RIGHT", "LOWER_LEFT", "LOWER_RIGHT"),
715 corners):
716 md.add(f"{box}_{name}_RA", geom.SpherePoint(corner).getRa().asDegrees(), f"{box}_corner")
717 md.add(f"{box}_{name}_DEC", geom.SpherePoint(corner).getDec().asDegrees(), f"{box}_corner")
718 return md
719
720 def joinMatchListWithCatalog(self, matchCat, sourceCat):
721 """Relink an unpersisted match list to sources and reference objects.
722
723 A match list is persisted and unpersisted as a catalog of IDs
724 produced by afw.table.packMatches(), with match metadata
725 (as returned by the astrometry tasks) in the catalog's metadata
726 attribute. This method converts such a match catalog into a match
727 list, with links to source records and reference object records.
728
729 Parameters
730 ----------
731 matchCat : `lsst.afw.table.BaseCatalog`
732 Unpersisted packed match list.
733 ``matchCat.table.getMetadata()`` must contain match metadata,
734 as returned by the astrometry tasks.
735 sourceCat : `lsst.afw.table.SourceCatalog`
736 Source catalog. As a side effect, the catalog will be sorted
737 by ID.
738
739 Returns
740 -------
742 Match list.
743 """
744 return joinMatchListWithCatalogImpl(self, matchCat, sourceCat)
745
746
748 """This class facilitates loading reference catalogs with gen 3 middleware.
749
750 The QuantumGraph generation will create a list of datasets that may
751 possibly overlap a given region. These datasets are then used to construct
752 and instance of this class. The class instance should then be passed into
753 a task which needs reference catalogs. These tasks should then determine
754 the exact region of the sky reference catalogs will be loaded for, and
755 call a corresponding method to load the reference objects.
756
757 Parameters
758 ----------
759 dataIds : iterable of `lsst.daf.butler.DataCoordinate`
760 An iterable object of data IDs that point to reference catalogs
761 in a gen 3 repository.
762 refCats : iterable of `lsst.daf.butler.DeferredDatasetHandle`
763 Handles to load refCats on demand
764 log : `lsst.log.Log`, `logging.Logger` or `None`, optional
765 Logger object used to write out messages. If `None` a default
766 logger will be used.
767 """
768 def __init__(self, dataIds, refCats, log=None, config=None, **kwargs):
769 if config is None:
770 config = self.ConfigClassConfigClass()
771 super().__init__(config=config, **kwargs)
772 self.dataIdsdataIds = dataIds
773 self.refCatsrefCats = refCats
774 self.loglog = log or logging.getLogger(__name__).getChild("ReferenceObjectLoader")
775
776 def loadPixelBox(self, bbox, wcs, filterName, epoch=None,
777 bboxToSpherePadding=100):
778 """Load reference objects that are within a pixel-based rectangular
779 region.
780
781 This algorithm works by creating a spherical box whose corners
782 correspond to the WCS converted corners of the input bounding box
783 (possibly padded). It then defines a filtering function which looks at
784 the pixel position of the reference objects and accepts only those that
785 lie within the specified bounding box.
786
787 The spherical box region and filtering function are passed to the
788 generic loadRegion method which loads and filters the reference objects
789 from the datastore and returns a single catalog containing the filtered
790 set of reference objects.
791
792 Parameters
793 ----------
794 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D`
795 Box which bounds a region in pixel space.
797 Wcs object defining the pixel to sky (and inverse) transform for
798 the supplied ``bbox``.
799 filterName : `str`
800 Name of camera filter.
801 epoch : `astropy.time.Time` or `None`, optional
802 Epoch to which to correct proper motion and parallax, or `None`
803 to not apply such corrections.
804 bboxToSpherePadding : `int`, optional
805 Padding to account for translating a set of corners into a
806 spherical (convex) boundary that is certain to encompase the
807 enitre area covered by the box.
808
809 Returns
810 -------
811 output : `lsst.pipe.base.Struct`
812 Results struct with attributes:
813
814 ``refCat``
815 Catalog containing reference objects inside the specified
816 bounding box (padded by self.configconfig.pixelMargin).
817 ``fluxField``
818 Name of the field containing the flux associated with
819 ``filterName``.
820
821 Raises
822 ------
823 RuntimeError
824 Raised if no reference catalogs could be found for the specified
825 region.
826 TypeError
827 Raised if the loaded reference catalogs do not have matching
828 schemas.
829 """
830 paddedBbox = geom.Box2D(bbox)
831 paddedBbox.grow(self.configconfig.pixelMargin)
832 innerSkyRegion, outerSkyRegion, _, _ = self._makeBoxRegion_makeBoxRegion(paddedBbox, wcs, bboxToSpherePadding)
833
834 def _filterFunction(refCat, region):
835 # Perform an initial "pre filter" step based on the refCat coords
836 # and the outerSkyRegion created from the self.config.pixelMargin-
837 # paddedBbox plus an "extra" padding of bboxToSpherePadding and the
838 # raw wcs. This should ensure a large enough projected area on the
839 # sky that accounts for any projection/distortion issues, but small
840 # enough to filter out loaded reference objects that lie well
841 # beyond the projected detector of interest. This step is required
842 # due to the very local nature of the wcs available for the
843 # sky <--> pixel conversions.
844 preFiltFunc = _FilterCatalog(outerSkyRegion)
845 refCat = preFiltFunc(refCat, region)
846
847 # Add columns to the pre-filtered reference catalog relating their
848 # coordinates to equivalent pixel positions for the wcs provided
849 # and use to populate those columns.
850 refCat = self._remapReferenceCatalogSchema_remapReferenceCatalogSchema(refCat, centroids=True)
851 afwTable.updateRefCentroids(wcs, refCat)
852 # No need to filter the catalog if it is entirely contained in the
853 # region defined by the inner sky region.
854 if innerSkyRegion.contains(region):
855 return refCat
856 # Create a new reference catalog, and populate it only with records
857 # that fall inside the padded bbox.
858 filteredRefCat = type(refCat)(refCat.table)
859 centroidKey = afwTable.Point2DKey(refCat.schema['centroid'])
860 for record in refCat:
861 pixCoords = record[centroidKey]
862 if paddedBbox.contains(geom.Point2D(pixCoords)):
863 filteredRefCat.append(record)
864 return filteredRefCat
865 return self.loadRegionloadRegion(outerSkyRegion, filterName, filtFunc=_filterFunction, epoch=epoch)
866
867 def loadRegion(self, region, filterName, filtFunc=None, epoch=None):
868 """Load reference objects within a specified region.
869
870 This function loads the DataIds used to construct an instance of this
871 class which intersect or are contained within the specified region. The
872 reference catalogs which intersect but are not fully contained within
873 the input region are further filtered by the specified filter function.
874 This function returns a single source catalog containing all reference
875 objects inside the specified region.
876
877 Parameters
878 ----------
879 region : `lsst.sphgeom.Region`
880 This can be any type that is derived from `lsst.sphgeom.Region` and
881 should define the spatial region for which reference objects are to
882 be loaded.
883 filtFunc : callable or `None`, optional
884 This optional parameter should be a callable object that takes a
885 reference catalog and its corresponding region as parameters,
886 filters the catalog by some criteria and returns the filtered
887 reference catalog. If `None`, an internal filter function is used
888 which filters according to if a reference object falls within the
889 input region.
890 filterName : `str`
891 Name of camera filter.
892 epoch : `astropy.time.Time` or `None`, optional
893 Epoch to which to correct proper motion and parallax, or `None` to
894 not apply such corrections.
895
896 Returns
897 -------
898 output : `lsst.pipe.base.Struct`
899 Results struct with attributes:
900
901 ``refCat``
902 Catalog containing reference objects which intersect the
903 input region, filtered by the specified filter function.
904 ``fluxField``
905 Name of the field containing the flux associated with
906 ``filterName``.
907
908 Raises
909 ------
910 RuntimeError
911 Raised if no reference catalogs could be found for the specified
912 region.
913 TypeError
914 Raised if the loaded reference catalogs do not have matching
915 schemas.
916 """
917 regionLat = region.getBoundingBox().getLat()
918 regionLon = region.getBoundingBox().getLon()
919 self.loglog.info("Loading reference objects from %s in region bounded by "
920 "[%.8f, %.8f], [%.8f, %.8f] RA Dec",
921 # Name of refcat we're loading from is the datasetType.
922 self.refCatsrefCats[0].ref.datasetType.name,
923 regionLon.getA().asDegrees(), regionLon.getB().asDegrees(),
924 regionLat.getA().asDegrees(), regionLat.getB().asDegrees())
925 if filtFunc is None:
926 filtFunc = _FilterCatalog(region)
927 # filter out all the regions supplied by the constructor that do not overlap
928 overlapList = []
929 for dataId, refCat in zip(self.dataIdsdataIds, self.refCatsrefCats):
930 # SphGeom supports some objects intersecting others, but is not symmetric,
931 # try the intersect operation in both directions
932 try:
933 intersects = dataId.region.intersects(region)
934 except TypeError:
935 intersects = region.intersects(dataId.region)
936
937 if intersects:
938 overlapList.append((dataId, refCat))
939
940 if len(overlapList) == 0:
941 raise RuntimeError("No reference tables could be found for input region")
942
943 firstCat = overlapList[0][1].get()
944 refCat = filtFunc(firstCat, overlapList[0][0].region)
945 trimmedAmount = len(firstCat) - len(refCat)
946
947 # Load in the remaining catalogs
948 for dataId, inputRefCat in overlapList[1:]:
949 tmpCat = inputRefCat.get()
950
951 if tmpCat.schema != firstCat.schema:
952 raise TypeError("Reference catalogs have mismatching schemas")
953
954 filteredCat = filtFunc(tmpCat, dataId.region)
955 refCat.extend(filteredCat)
956 trimmedAmount += len(tmpCat) - len(filteredCat)
957
958 self.loglog.debug("Trimmed %d refCat objects lying outside padded region, leaving %d",
959 trimmedAmount, len(refCat))
960 self.loglog.info("Loaded %d reference objects", len(refCat))
961
962 # Ensure that the loaded reference catalog is continuous in memory
963 if not refCat.isContiguous():
964 refCat = refCat.copy(deep=True)
965
966 self.applyProperMotionsapplyProperMotions(refCat, epoch)
967
968 # Verify the schema is in the correct units and has the correct version; automatically convert
969 # it with a warning if this is not the case.
970 if not hasNanojanskyFluxUnits(refCat.schema) or not getFormatVersionFromRefCat(refCat) >= 1:
971 self.loglog.warning("Found version 0 reference catalog with old style units in schema.")
972 self.loglog.warning("run `meas_algorithms/bin/convert_refcat_to_nJy.py` to convert fluxes to nJy.")
973 self.loglog.warning("See RFC-575 for more details.")
974 refCat = convertToNanojansky(refCat, self.loglog)
975
976 expandedCat = self._remapReferenceCatalogSchema_remapReferenceCatalogSchema(refCat,
977 anyFilterMapsToThis=self.configconfig.anyFilterMapsToThis,
978 filterMap=self.configconfig.filterMap)
979
980 # Ensure that the returned reference catalog is continuous in memory
981 if not expandedCat.isContiguous():
982 expandedCat = expandedCat.copy(deep=True)
983
984 fluxField = getRefFluxField(expandedCat.schema, filterName)
985 return pipeBase.Struct(refCat=expandedCat, fluxField=fluxField)
986
987 def loadSkyCircle(self, ctrCoord, radius, filterName, epoch=None):
988 """Load reference objects that lie within a circular region on the sky.
989
990 This method constructs a circular region from an input center and
991 angular radius, loads reference catalogs which are contained in or
992 intersect the circle, and filters reference catalogs which intersect
993 down to objects which lie within the defined circle.
994
995 Parameters
996 ----------
997 ctrCoord : `lsst.geom.SpherePoint`
998 Point defining the center of the circular region.
999 radius : `lsst.geom.Angle`
1000 Defines the angular radius of the circular region.
1001 filterName : `str`
1002 Name of camera filter.
1003 epoch : `astropy.time.Time` or `None`, optional
1004 Epoch to which to correct proper motion and parallax, or `None` to
1005 not apply such corrections.
1006
1007 Returns
1008 -------
1009 output : `lsst.pipe.base.Struct`
1010 Results struct with attributes:
1011
1012 ``refCat``
1013 Catalog containing reference objects inside the specified
1014 search circle.
1015 ``fluxField``
1016 Name of the field containing the flux associated with
1017 ``filterName``.
1018 """
1019 centerVector = ctrCoord.getVector()
1020 sphRadius = sphgeom.Angle(radius.asRadians())
1021 circularRegion = sphgeom.Circle(centerVector, sphRadius)
1022 return self.loadRegionloadRegion(circularRegion, filterName, epoch=epoch)
1023
1024
1025def getRefFluxField(schema, filterName):
1026 """Get the name of a flux field from a schema.
1027
1028 return the alias of "anyFilterMapsToThis", if present
1029 else:
1030 return "*filterName*_camFlux" if present
1031 else return "*filterName*_flux" if present (camera filter name
1032 matches reference filter name)
1033 else throw RuntimeError
1034
1035 Parameters
1036 ----------
1037 schema : `lsst.afw.table.Schema`
1038 Reference catalog schema.
1039 filterName : `str`
1040 Name of camera filter.
1041
1042 Returns
1043 -------
1044 fluxFieldName : `str`
1045 Name of flux field.
1046
1047 Raises
1048 ------
1049 RuntimeError
1050 If an appropriate field is not found.
1051 """
1052 if not isinstance(schema, afwTable.Schema):
1053 raise RuntimeError("schema=%s is not a schema" % (schema,))
1054 try:
1055 return schema.getAliasMap().get("anyFilterMapsToThis")
1056 except LookupError:
1057 pass # try the filterMap next
1058
1059 fluxFieldList = [filterName + "_camFlux", filterName + "_flux"]
1060 for fluxField in fluxFieldList:
1061 if fluxField in schema:
1062 return fluxField
1063
1064 raise RuntimeError("Could not find flux field(s) %s" % (", ".join(fluxFieldList)))
1065
1066
1067def getRefFluxKeys(schema, filterName):
1068 """Return keys for flux and flux error.
1069
1070 Parameters
1071 ----------
1072 schema : `lsst.afw.table.Schema`
1073 Reference catalog schema.
1074 filterName : `str`
1075 Name of camera filter.
1076
1077 Returns
1078 -------
1079 keys : `tuple` of (`lsst.afw.table.Key`, `lsst.afw.table.Key`)
1080 Two keys:
1081
1082 - flux key
1083 - flux error key, if present, else None
1084
1085 Raises
1086 ------
1087 RuntimeError
1088 If flux field not found.
1089 """
1090 fluxField = getRefFluxField(schema, filterName)
1091 fluxErrField = fluxField + "Err"
1092 fluxKey = schema[fluxField].asKey()
1093 try:
1094 fluxErrKey = schema[fluxErrField].asKey()
1095 except Exception:
1096 fluxErrKey = None
1097 return (fluxKey, fluxErrKey)
1098
1099
1100class LoadReferenceObjectsTask(pipeBase.Task, ReferenceObjectLoaderBase, metaclass=abc.ABCMeta):
1101 """Abstract gen2 base class to load objects from reference catalogs.
1102 """
1103 _DefaultName = "LoadReferenceObjects"
1104
1105 def __init__(self, butler=None, *args, **kwargs):
1106 """Construct a LoadReferenceObjectsTask
1107
1108 Parameters
1109 ----------
1111 Data butler, for access reference catalogs.
1112 """
1113 pipeBase.Task.__init__(self, *args, **kwargs)
1114 self.butlerbutler = butler
1115
1116 @timeMethod
1117 def loadPixelBox(self, bbox, wcs, filterName, photoCalib=None, epoch=None):
1118 """Load reference objects that overlap a rectangular pixel region.
1119
1120 Parameters
1121 ----------
1122 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D`
1123 Bounding box for pixels.
1124 wcs : `lsst.afw.geom.SkyWcs`
1125 WCS; used to convert pixel positions to sky coordinates
1126 and vice-versa.
1127 filterName : `str`
1128 Name of filter. This can be used for flux limit comparisons.
1129 photoCalib : `None`
1130 Deprecated, only included for api compatibility.
1131 epoch : `astropy.time.Time` or `None`, optional
1132 Epoch to which to correct proper motion and parallax, or `None` to
1133 not apply such corrections.
1134
1135 Returns
1136 -------
1137 results : `lsst.pipe.base.Struct`
1138 A Struct containing the following fields:
1139 refCat : `lsst.afw.catalog.SimpleCatalog`
1140 A catalog of reference objects with the standard
1141 schema, as documented in the main doc string for
1142 `LoadReferenceObjects`.
1143 The catalog is guaranteed to be contiguous.
1144 fluxField : `str`
1145 Name of flux field for specified `filterName`.
1146
1147 Notes
1148 -----
1149 The search algorithm works by searching in a region in sky
1150 coordinates whose center is the center of the bbox and radius
1151 is large enough to just include all 4 corners of the bbox.
1152 Stars that lie outside the bbox are then trimmed from the list.
1153 """
1154 circle = self._calculateCircle_calculateCircle(bbox, wcs, self.configconfig.pixelMargin)
1155
1156 # find objects in circle
1157 self.log.info("Loading reference objects from %s using center %s and radius %s deg",
1158 self.configconfig.ref_dataset_name, circle.coord, circle.radius.asDegrees())
1159 loadRes = self.loadSkyCircleloadSkyCircle(circle.coord, circle.radius, filterName, epoch=epoch,
1160 centroids=True)
1161 refCat = loadRes.refCat
1162 numFound = len(refCat)
1163
1164 # trim objects outside bbox
1165 refCat = self._trimToBBox_trimToBBox(refCat=refCat, bbox=circle.bbox, wcs=wcs)
1166 numTrimmed = numFound - len(refCat)
1167 self.log.debug("trimmed %d out-of-bbox objects, leaving %d", numTrimmed, len(refCat))
1168 self.log.info("Loaded %d reference objects", len(refCat))
1169
1170 # make sure catalog is contiguous
1171 if not refCat.isContiguous():
1172 loadRes.refCat = refCat.copy(deep=True)
1173
1174 return loadRes
1175
1176 @abc.abstractmethod
1177 def loadSkyCircle(self, ctrCoord, radius, filterName, epoch=None, centroids=False):
1178 """Load reference objects that overlap a circular sky region.
1179
1180 Parameters
1181 ----------
1182 ctrCoord : `lsst.geom.SpherePoint`
1183 ICRS center of search region.
1184 radius : `lsst.geom.Angle`
1185 Radius of search region.
1186 filterName : `str`
1187 Name of filter. This can be used for flux limit comparisons.
1188 epoch : `astropy.time.Time` or `None`, optional
1189 Epoch to which to correct proper motion and parallax, or `None` to
1190 not apply such corrections.
1191 centroids : `bool`, optional
1192 Add centroid fields to the loaded Schema. ``loadPixelBox`` expects
1193 these fields to exist.
1194
1195 Returns
1196 -------
1197 results : `lsst.pipe.base.Struct`
1198 A Struct containing the following fields:
1199 refCat : `lsst.afw.catalog.SimpleCatalog`
1200 A catalog of reference objects with the standard
1201 schema, as documented in the main doc string for
1202 `LoadReferenceObjects`.
1203 The catalog is guaranteed to be contiguous.
1204 fluxField : `str`
1205 Name of flux field for specified `filterName`.
1206
1207 Notes
1208 -----
1209 Note that subclasses are responsible for performing the proper motion
1210 correction, since this is the lowest-level interface for retrieving
1211 the catalog.
1212 """
1213 return
1214
1215 @staticmethod
1216 def _trimToBBox(refCat, bbox, wcs):
1217 """Remove objects outside a given pixel bounding box and set
1218 centroid and hasCentroid fields.
1219
1220 Parameters
1221 ----------
1223 A catalog of objects. The schema must include fields
1224 "coord", "centroid" and "hasCentroid".
1225 The "coord" field is read.
1226 The "centroid" and "hasCentroid" fields are set.
1227 bbox : `lsst.geom.Box2D`
1228 Pixel region
1229 wcs : `lsst.afw.geom.SkyWcs`
1230 WCS; used to convert sky coordinates to pixel positions.
1231
1232 Returns
1233 -------
1235 Reference objects in the bbox, with centroid and
1236 hasCentroid fields set.
1237 """
1238 afwTable.updateRefCentroids(wcs, refCat)
1239 centroidKey = afwTable.Point2DKey(refCat.schema["centroid"])
1240 retStarCat = type(refCat)(refCat.table)
1241 for star in refCat:
1242 point = star.get(centroidKey)
1243 if bbox.contains(point):
1244 retStarCat.append(star)
1245 return retStarCat
1246
1247
1248def joinMatchListWithCatalogImpl(refObjLoader, matchCat, sourceCat):
1249 """Relink an unpersisted match list to sources and reference
1250 objects.
1251
1252 A match list is persisted and unpersisted as a catalog of IDs
1253 produced by afw.table.packMatches(), with match metadata
1254 (as returned by the astrometry tasks) in the catalog's metadata
1255 attribute. This method converts such a match catalog into a match
1256 list, with links to source records and reference object records.
1257
1258 Parameters
1259 ----------
1260 refObjLoader
1261 Reference object loader to use in getting reference objects
1262 matchCat : `lsst.afw.table.BaseCatalog`
1263 Unperisted packed match list.
1264 ``matchCat.table.getMetadata()`` must contain match metadata,
1265 as returned by the astrometry tasks.
1266 sourceCat : `lsst.afw.table.SourceCatalog`
1267 Source catalog. As a side effect, the catalog will be sorted
1268 by ID.
1269
1270 Returns
1271 -------
1273 Match list.
1274 """
1275 matchmeta = matchCat.table.getMetadata()
1276 version = matchmeta.getInt('SMATCHV')
1277 if version != 1:
1278 raise ValueError('SourceMatchVector version number is %i, not 1.' % version)
1279 filterName = matchmeta.getString('FILTER').strip()
1280 try:
1281 epoch = matchmeta.getDouble('EPOCH')
1282 except (LookupError, TypeError):
1283 epoch = None # Not present, or not correct type means it's not set
1284 if 'RADIUS' in matchmeta:
1285 # This is a circle style metadata, call loadSkyCircle
1286 ctrCoord = geom.SpherePoint(matchmeta.getDouble('RA'),
1287 matchmeta.getDouble('DEC'), geom.degrees)
1288 rad = matchmeta.getDouble('RADIUS')*geom.degrees
1289 refCat = refObjLoader.loadSkyCircle(ctrCoord, rad, filterName, epoch=epoch).refCat
1290 elif "INNER_UPPER_LEFT_RA" in matchmeta:
1291 # This is the sky box type (only triggers in the LoadReferenceObject class, not task)
1292 # Only the outer box is required to be loaded to get the maximum region, all filtering
1293 # will be done by the unpackMatches function, and no spatial filtering needs to be done
1294 # by the refObjLoader
1295 box = []
1296 for place in ("UPPER_LEFT", "UPPER_RIGHT", "LOWER_LEFT", "LOWER_RIGHT"):
1297 coord = geom.SpherePoint(matchmeta.getDouble(f"OUTER_{place}_RA"),
1298 matchmeta.getDouble(f"OUTER_{place}_DEC"),
1299 geom.degrees).getVector()
1300 box.append(coord)
1301 outerBox = sphgeom.ConvexPolygon(box)
1302 refCat = refObjLoader.loadRegion(outerBox, filterName, epoch=epoch).refCat
1303
1304 refCat.sort()
1305 sourceCat.sort()
1306 return afwTable.unpackMatches(matchCat, refCat, sourceCat)
1307
1308
1309def applyProperMotionsImpl(log, catalog, epoch):
1310 """Apply proper motion correction to a reference catalog.
1311
1312 Adjust position and position error in the ``catalog``
1313 for proper motion to the specified ``epoch``,
1314 modifying the catalog in place.
1315
1316 Parameters
1317 ----------
1318 log : `lsst.log.Log` or `logging.getLogger`
1319 Log object to write to.
1321 Catalog of positions, containing:
1322
1323 - Coordinates, retrieved by the table's coordinate key.
1324 - ``coord_raErr`` : Error in Right Ascension (rad).
1325 - ``coord_decErr`` : Error in Declination (rad).
1326 - ``pm_ra`` : Proper motion in Right Ascension (rad/yr,
1327 East positive)
1328 - ``pm_raErr`` : Error in ``pm_ra`` (rad/yr), optional.
1329 - ``pm_dec`` : Proper motion in Declination (rad/yr,
1330 North positive)
1331 - ``pm_decErr`` : Error in ``pm_dec`` (rad/yr), optional.
1332 - ``epoch`` : Mean epoch of object (an astropy.time.Time)
1333 epoch : `astropy.time.Time`
1334 Epoch to which to correct proper motion.
1335 """
1336 if "epoch" not in catalog.schema or "pm_ra" not in catalog.schema or "pm_dec" not in catalog.schema:
1337 log.warning("Proper motion correction not available from catalog")
1338 return
1339 if not catalog.isContiguous():
1340 raise RuntimeError("Catalog must be contiguous")
1341 catEpoch = astropy.time.Time(catalog["epoch"], scale="tai", format="mjd")
1342 log.info("Correcting reference catalog for proper motion to %r", epoch)
1343 # Use `epoch.tai` to make sure the time difference is in TAI
1344 timeDiffsYears = (epoch.tai - catEpoch).to(astropy.units.yr).value
1345 coordKey = catalog.table.getCoordKey()
1346 # Compute the offset of each object due to proper motion
1347 # as components of the arc of a great circle along RA and Dec
1348 pmRaRad = catalog["pm_ra"]
1349 pmDecRad = catalog["pm_dec"]
1350 offsetsRaRad = pmRaRad*timeDiffsYears
1351 offsetsDecRad = pmDecRad*timeDiffsYears
1352 # Compute the corresponding bearing and arc length of each offset
1353 # due to proper motion, and apply the offset
1354 # The factor of 1e6 for computing bearing is intended as
1355 # a reasonable scale for typical values of proper motion
1356 # in order to avoid large errors for small values of proper motion;
1357 # using the offsets is another option, but it can give
1358 # needlessly large errors for short duration
1359 offsetBearingsRad = numpy.arctan2(pmDecRad*1e6, pmRaRad*1e6)
1360 offsetAmountsRad = numpy.hypot(offsetsRaRad, offsetsDecRad)
1361 for record, bearingRad, amountRad in zip(catalog, offsetBearingsRad, offsetAmountsRad):
1362 record.set(coordKey,
1363 record.get(coordKey).offset(bearing=bearingRad*geom.radians,
1364 amount=amountRad*geom.radians))
1365 # Increase error in RA and Dec based on error in proper motion
1366 if "coord_raErr" in catalog.schema:
1367 catalog["coord_raErr"] = numpy.hypot(catalog["coord_raErr"],
1368 catalog["pm_raErr"]*timeDiffsYears)
1369 if "coord_decErr" in catalog.schema:
1370 catalog["coord_decErr"] = numpy.hypot(catalog["coord_decErr"],
1371 catalog["pm_decErr"]*timeDiffsYears)
int max
table::Key< int > type
Definition: Detector.cc:163
table::Key< int > to
A 2-dimensional celestial WCS that transform pixels to ICRS RA/Dec, using the LSST standard for pixel...
Definition: SkyWcs.h:117
A class used as a handle to a particular field in a table.
Definition: Key.h:53
Defines the fields and offsets for a table.
Definition: Schema.h:51
A mapping between the keys of two Schemas, used to copy data between them.
Definition: SchemaMapper.h:21
Custom catalog class for record/table subclasses that are guaranteed to have an ID,...
Definition: SortedCatalog.h:42
Class for storing ordered metadata with comments.
Definition: PropertyList.h:68
A class representing an angle.
Definition: Angle.h:128
A floating-point coordinate rectangle geometry.
Definition: Box.h:413
An integer coordinate rectangle.
Definition: Box.h:55
Point in an unspecified spherical coordinate system.
Definition: SpherePoint.h:57
This static class includes a variety of methods for interacting with the the logging module.
Definition: Log.h:724
def loadSkyCircle(self, ctrCoord, radius, filterName, epoch=None, centroids=False)
def loadPixelBox(self, bbox, wcs, filterName, photoCalib=None, epoch=None)
def getMetadataBox(self, bbox, wcs, filterName, epoch=None, bboxToSpherePadding=100)
def makeMinimalSchema(filterNameList, *addCentroid=False, addIsPhotometric=False, addIsResolved=False, addIsVariable=False, coordErrDim=2, addProperMotion=False, properMotionErrDim=2, addParallax=False)
def _remapReferenceCatalogSchema(refCat, *anyFilterMapsToThis=None, filterMap=None, centroids=False)
def loadSkyCircle(self, ctrCoord, radius, filterName, epoch=None)
def __init__(self, dataIds, refCats, log=None, config=None, **kwargs)
def loadRegion(self, region, filterName, filtFunc=None, epoch=None)
def loadPixelBox(self, bbox, wcs, filterName, epoch=None, bboxToSpherePadding=100)
Angle represents an angle in radians.
Definition: Angle.h:43
Circle is a circular region on the unit sphere that contains its boundary.
Definition: Circle.h:46
ConvexPolygon is a closed convex polygon on the unit sphere.
Definition: ConvexPolygon.h:57
Region is a minimal interface for 2-dimensional regions on the unit sphere.
Definition: Region.h:79
daf::base::PropertySet * set
Definition: fits.cc:912
bool strip
Definition: fits.cc:911
std::vector< Match< typename Cat1::Record, typename Cat2::Record > > unpackMatches(BaseCatalog const &matches, Cat1 const &cat1, Cat2 const &cat2)
Reconstruct a MatchVector from a BaseCatalog representation of the matches and a pair of catalogs.
Definition: Match.cc:455
void updateRefCentroids(geom::SkyWcs const &wcs, ReferenceCollection &refList)
Update centroids in a collection of reference objects.
Definition: wcsUtils.cc:72
def convertToNanojansky(catalog, log, doConvert=True)
def joinMatchListWithCatalogImpl(refObjLoader, matchCat, sourceCat)
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:174
A description of a field in a table.
Definition: Field.h:24