24 __all__ = [
"getRefFluxField",
"getRefFluxKeys",
"LoadReferenceObjectsTask",
"LoadReferenceObjectsConfig",
25 "ReferenceObjectLoader"]
42 from lsst
import sphgeom
47 """Return True if this name/units combination corresponds to an
48 "old-style" reference catalog flux field.
50 unitsCheck = units !=
'nJy'
51 isFlux = name.endswith(
'_flux')
52 isFluxSigma = name.endswith(
'_fluxSigma')
53 isFluxErr = name.endswith(
'_fluxErr')
54 return (isFlux
or isFluxSigma
or isFluxErr)
and unitsCheck
58 """Return True if the units of all flux and fluxErr are correct (nJy).
67 """"Return the format version stored in a reference catalog header.
71 refCat : `lsst.afw.table.SimpleCatalog`
72 Reference catalog to inspect.
76 version : `int` or `None`
77 Format version integer, or `None` if the catalog has no metadata
78 or the metadata does not include a "REFCAT_FORMAT_VERSION" key.
80 md = refCat.getMetadata()
84 return md.getScalar(
"REFCAT_FORMAT_VERSION")
90 """Convert fluxes in a catalog from jansky to nanojansky.
94 catalog : `lsst.afw.table.SimpleCatalog`
95 The catalog to convert.
97 Log to send messages to.
98 doConvert : `bool`, optional
99 Return a converted catalog, or just identify the fields that need to be converted?
100 This supports the "write=False" mode of `bin/convert_to_nJy.py`.
104 catalog : `lsst.afw.table.SimpleCatalog` or None
105 The converted catalog, or None if ``doConvert`` is False.
109 Support for old units in reference catalogs will be removed after the
110 release of late calendar year 2019.
111 Use `meas_algorithms/bin/convert_to_nJy.py` to update your reference catalog.
119 for field
in catalog.schema:
120 oldName = field.field.getName()
121 oldUnits = field.field.getUnits()
125 if oldName.endswith(
'_fluxSigma'):
126 name = oldName.replace(
'_fluxSigma',
'_fluxErr')
130 mapper.addMapping(field.getKey(), newField)
131 input_fields.append(field.field)
132 output_fields.append(newField)
134 mapper.addMapping(field.getKey())
136 fluxFieldsStr =
'; '.join(
"(%s, '%s')" % (field.getName(), field.getUnits())
for field
in input_fields)
139 newSchema = mapper.getOutputSchema()
141 output.extend(catalog, mapper=mapper)
142 for field
in output_fields:
143 output[field.getName()] *= 1e9
144 log.info(f
"Converted refcat flux fields to nJy (name, units): {fluxFieldsStr}")
147 log.info(f
"Found old-style refcat flux fields (name, units): {fluxFieldsStr}")
152 """This is a private helper class which filters catalogs by
153 row based on the row being inside the region used to initialize
158 region : `lsst.sphgeom.Region`
159 The spatial region which all objects should lie within
165 """This call method on an instance of this class takes in a reference
166 catalog, and the region from which the catalog was generated.
168 If the catalog region is entirely contained within the region used to
169 initialize this class, then all the entries in the catalog must be
170 within the region and so the whole catalog is returned.
172 If the catalog region is not entirely contained, then the location for
173 each record is tested against the region used to initialize the class.
174 Records which fall inside this region are added to a new catalog, and
175 this catalog is then returned.
179 refCat : `lsst.afw.table.SourceCatalog`
180 SourceCatalog to be filtered.
181 catRegion : `lsst.sphgeom.Region`
182 Region in which the catalog was created
184 if catRegion.isWithin(self.
region):
188 filteredRefCat =
type(refCat)(refCat.table)
189 for record
in refCat:
191 filteredRefCat.append(record)
192 return filteredRefCat
196 """ This class facilitates loading reference catalogs with gen 3 middleware
198 The middleware preflight solver will create a list of datarefs that may
199 possibly overlap a given region. These datarefs are then used to construct
200 and instance of this class. The class instance should then be passed into
201 a task which needs reference catalogs. These tasks should then determine
202 the exact region of the sky reference catalogs will be loaded for, and
203 call a corresponding method to load the reference objects.
205 def __init__(self, dataIds, refCats, config, log=None):
206 """ Constructs an instance of ReferenceObjectLoader
210 dataIds : iterable of `lsst.daf.butler.DataIds`
211 An iterable object of DataSetRefs which point to reference catalogs
212 in a gen 3 repository
213 refCats : Iterable of `lsst.daf.butler.DeferedDatasetHandle`
214 Handles to load refCats on demand
216 Logger object used to write out messages. If `None` (default) the default
217 lsst logger will be used
226 def _makeBoxRegion(BBox, wcs, BBoxPadding):
233 outerLocalBBox.grow(BBoxPadding)
234 innerLocalBBox.grow(-1*BBoxPadding)
242 innerBoxCorners = innerLocalBBox.getCorners()
243 innerSphCorners = [wcs.pixelToSky(corner).getVector()
for corner
in innerBoxCorners]
248 outerBoxCorners = outerLocalBBox.getCorners()
249 outerSphCorners = [wcs.pixelToSky(corner).getVector()
for corner
in outerBoxCorners]
253 return innerSkyRegion, outerSkyRegion, innerSphCorners, outerSphCorners
255 def loadPixelBox(self, bbox, wcs, filterName=None, epoch=None, photoCalib=None, bboxPadding=100):
256 """Load reference objects that are within a pixel-based rectangular region
258 This algorithm works by creating a spherical box whose corners correspond
259 to the WCS converted corners of the input bounding box (possibly padded).
260 It then defines a filtering function which will look at a reference
261 objects pixel position and accept objects that lie within the specified
264 The spherical box region and filtering function are passed to the generic
265 loadRegion method which will load and filter the reference objects from
266 the datastore and return a single catalog containing all reference objects
270 bbox : `lsst.geom.box2I`
271 Box which bounds a region in pixel space
272 wcs : `lsst.afw.geom.SkyWcs`
273 Wcs object defining the pixel to sky (and inverse) transform for the space
274 of pixels of the supplied bbox
276 Name of camera filter, or None or blank for the default filter
277 epoch : `astropy.time.Time` (optional)
278 Epoch to which to correct proper motion and parallax,
279 or None to not apply such corrections.
281 Deprecated and ignored, only included for api compatibility
283 Number describing how much to pad the input bbox by (in pixels), defaults
284 to 100. This parameter is necessary because optical distortions in telescopes
285 can cause a rectangular pixel grid to map into a non "rectangular" spherical
286 region in sky coordinates. This padding is used to create a spherical
287 "rectangle", which will for sure enclose the input box. This padding is only
288 used to determine if the reference catalog for a sky patch will be loaded from
289 the data store, this function will filter out objects which lie within the
290 padded region but fall outside the input bounding box region.
294 referenceCatalog : `lsst.afw.table.SimpleCatalog`
295 Catalog containing reference objects inside the specified bounding box
299 `lsst.pex.exception.RuntimeError`
300 Raised if no reference catalogs could be found for the specified region
302 `lsst.pex.exception.TypeError`
303 Raised if the loaded reference catalogs do not have matching schemas
305 innerSkyRegion, outerSkyRegion, _, _ = self.
_makeBoxRegion(bbox, wcs, bboxPadding)
307 def _filterFunction(refCat, region):
314 if innerSkyRegion.contains(region):
317 filteredRefCat =
type(refCat)(refCat.table)
319 for record
in refCat:
320 pixCoords = record[centroidKey]
322 filteredRefCat.append(record)
323 return filteredRefCat
324 return self.
loadRegion(outerSkyRegion, filtFunc=_filterFunction, epoch=epoch, filterName=filterName)
326 def loadRegion(self, region, filtFunc=None, filterName=None, epoch=None):
327 """ Load reference objects within a specified region
329 This function loads the DataIds used to construct an instance of this class
330 which intersect or are contained within the specified region. The reference
331 catalogs which intersect but are not fully contained within the input region are
332 further filtered by the specified filter function. This function will return a
333 single source catalog containing all reference objects inside the specified region.
337 region : `lsst.sphgeom.Region`
338 This can be any type that is derived from `lsst.sphgeom.region` and should
339 define the spatial region for which reference objects are to be loaded.
341 This optional parameter should be a callable object that takes a reference
342 catalog and its corresponding region as parameters, filters the catalog by
343 some criteria and returns the filtered reference catalog. If the value is
344 left as the default (None) than an internal filter function is used which
345 filters according to if a reference object falls within the input region.
347 Name of camera filter, or None or blank for the default filter
348 epoch : `astropy.time.Time` (optional)
349 Epoch to which to correct proper motion and parallax,
350 or None to not apply such corrections.
354 referenceCatalog : `lsst.afw.table.SourceCatalog`
355 Catalog containing reference objects which intersect the input region,
356 filtered by the specified filter function
360 `lsst.pex.exception.RuntimeError`
361 Raised if no reference catalogs could be found for the specified region
363 `lsst.pex.exception.TypeError`
364 Raised if the loaded reference catalogs do not have matching schemas
367 regionBounding = region.getBoundingBox()
368 self.
log.
info(
"Loading reference objects from region bounded by {}, {} lat lon".
format(
369 regionBounding.getLat(), regionBounding.getLon()))
378 intersects = dataId.region.intersects(region)
380 intersects = region.intersects(dataId.region)
383 overlapList.append((dataId, refCat))
385 if len(overlapList) == 0:
388 firstCat = overlapList[0][1].get()
389 refCat = filtFunc(firstCat, overlapList[0][0].region)
390 trimmedAmount = len(firstCat) - len(refCat)
393 for dataId, inputRefCat
in overlapList[1:]:
394 tmpCat = inputRefCat.get()
396 if tmpCat.schema != firstCat.schema:
399 filteredCat = filtFunc(tmpCat, dataId.region)
400 refCat.extend(filteredCat)
401 trimmedAmount += len(tmpCat) - len(filteredCat)
403 self.
log.
debug(f
"Trimmed {trimmedAmount} out of region objects, leaving {len(refCat)}")
404 self.
log.
info(f
"Loaded {len(refCat)} reference objects")
407 if not refCat.isContiguous():
408 refCat = refCat.copy(deep=
True)
410 if epoch
is not None and "pm_ra" in refCat.schema:
412 if isinstance(refCat.schema[
"pm_ra"].asKey(), lsst.afw.table.KeyAngle):
415 self.
log.
warn(
"Catalog pm_ra field is not an Angle; not applying proper motion")
420 self.
log.
warn(
"Found version 0 reference catalog with old style units in schema.")
421 self.
log.
warn(
"run `meas_algorithms/bin/convert_refcat_to_nJy.py` to convert fluxes to nJy.")
422 self.
log.
warn(
"See RFC-575 for more details.")
431 if not expandedCat.isContiguous():
432 expandedCat = expandedCat.copy(deep=
True)
434 fluxField =
getRefFluxField(schema=expandedCat.schema, filterName=filterName)
435 return pipeBase.Struct(refCat=expandedCat, fluxField=fluxField)
438 """Load reference objects that lie within a circular region on the sky
440 This method constructs a circular region from an input center and angular radius,
441 loads reference catalogs which are contained in or intersect the circle, and
442 filters reference catalogs which intersect down to objects which lie within
447 ctrCoord : `lsst.geom.SpherePoint`
448 Point defining the center of the circular region
449 radius : `lsst.geom.Angle`
450 Defines the angular radius of the circular region
452 Name of camera filter, or None or blank for the default filter
453 epoch : `astropy.time.Time` (optional)
454 Epoch to which to correct proper motion and parallax,
455 or None to not apply such corrections.
459 referenceCatalog : `lsst.afw.table.SourceCatalog`
460 Catalog containing reference objects inside the specified bounding box
464 `lsst.pex.exception.RuntimeError`
465 Raised if no reference catalogs could be found for the specified region
467 `lsst.pex.exception.TypeError`
468 Raised if the loaded reference catalogs do not have matching schemas
471 centerVector = ctrCoord.getVector()
474 return self.
loadRegion(circularRegion, filterName=filterName, epoch=
None)
477 """Relink an unpersisted match list to sources and reference
480 A match list is persisted and unpersisted as a catalog of IDs
481 produced by afw.table.packMatches(), with match metadata
482 (as returned by the astrometry tasks) in the catalog's metadata
483 attribute. This method converts such a match catalog into a match
484 list, with links to source records and reference object records.
488 matchCat : `lsst.afw.table.BaseCatalog`
489 Unpersisted packed match list.
490 ``matchCat.table.getMetadata()`` must contain match metadata,
491 as returned by the astrometry tasks.
492 sourceCat : `lsst.afw.table.SourceCatalog`
493 Source catalog. As a side effect, the catalog will be sorted
498 matchList : `lsst.afw.table.ReferenceMatchVector`
504 def getMetadataBox(cls, bbox, wcs, filterName=None, photoCalib=None, epoch=None, bboxPadding=100):
505 """Return metadata about the load
507 This metadata is used for reloading the catalog (e.g., for
508 reconstituting a normalised match list.)
512 bbox : `lsst.geom.Box2I`
513 Bounding bos for the pixels
514 wcs : `lsst.afw.geom.SkyWcs
516 filterName : `str` or None
517 filterName of the camera filter, or None or blank for the default filter
519 Deprecated, only included for api compatibility
520 epoch : `astropy.time.Time` (optional)
521 Epoch to which to correct proper motion and parallax,
522 or None to not apply such corrections.
524 Number describing how much to pad the input bbox by (in pixels), defaults
525 to 100. This parameter is necessary because optical distortions in telescopes
526 can cause a rectangular pixel grid to map into a non "rectangular" spherical
527 region in sky coordinates. This padding is used to create a spherical
528 "rectangle", which will for sure enclose the input box. This padding is only
529 used to determine if the reference catalog for a sky patch will be loaded from
530 the data store, this function will filter out objects which lie within the
531 padded region but fall outside the input bounding box region.
534 md : `lsst.daf.base.PropertyList`
536 _, _, innerCorners, outerCorners = cls.
_makeBoxRegion(bbox, wcs, bboxPadding)
538 for box, corners
in zip((
"INNER",
"OUTER"), (innerCorners, outerCorners)):
539 for (name, corner)
in zip((
"UPPER_LEFT",
"UPPER_RIGHT",
"LOWER_LEFT",
"LOWER_RIGHT"),
541 md.add(f
"{box}_{name}_RA",
geom.SpherePoint(corner).getRa().asDegrees(), f
"{box}_corner")
542 md.add(f
"{box}_{name}_DEC",
geom.SpherePoint(corner).getDec().asDegrees(), f
"{box}_corner")
543 md.add(
"SMATCHV", 1,
'SourceMatchVector version number')
544 filterName =
"UNKNOWN" if filterName
is None else str(filterName)
545 md.add(
'FILTER', filterName,
'filter name for photometric data')
548 md.add(
'EPOCH',
"NONE" if epoch
is None else str(epoch),
'Epoch (TAI MJD) for catalog')
553 """Return metadata about the load
555 This metadata is used for reloading the catalog (e.g. for reconstituting
556 a normalized match list.)
560 coord : `lsst.geom.SpherePoint`
561 ICRS center of a circle
562 radius : `lsst.geom.angle`
564 filterName : `str` or None
565 filterName of the camera filter, or None or blank for the default filter
567 Deprecated, only included for api compatibility
568 epoch : `astropy.time.Time` (optional)
569 Epoch to which to correct proper motion and parallax,
570 or None to not apply such corrections.
574 md : `lsst.daf.base.PropertyList`
577 md.add(
'RA', coord.getRa().asDegrees(),
'field center in degrees')
578 md.add(
'DEC', coord.getDec().asDegrees(),
'field center in degrees')
579 md.add(
'RADIUS', radius.asDegrees(),
'field radius in degrees, minimum')
580 md.add(
'SMATCHV', 1,
'SourceMatchVector version number')
581 filterName =
"UNKNOWN" if filterName
is None else str(filterName)
582 md.add(
'FILTER', filterName,
'filter name for photometric data')
583 md.add(
'EPOCH',
"NONE" if epoch
is None else epoch,
'Epoch (TAI MJD) for catalog')
588 """This function creates a new catalog containing the information of the input refCat
589 as well as added flux columns and aliases between camera and reference flux.
593 refCat : `lsst.afw.table.SimpleCatalog`
594 Catalog of reference objects
595 defaultFilter : `str`
596 Name of the default reference filter
597 filterReferenceMap : `dict` of `str`
598 Dictionary with keys corresponding to a filter name, and values which
599 correspond to the name of the reference filter.
603 refCat : `lsst.afw.table.SimpleCatalog`
604 Reference catalog with columns added to track reference filters
609 If specified reference filter name is not a filter specifed as a key in the
610 reference filter map.
612 refCat = ReferenceObjectLoader.remapReferenceCatalogSchema(refCat,
613 filterNameList=filterReferenceMap.keys())
614 aliasMap = refCat.schema.getAliasMap()
615 if filterReferenceMap
is None:
616 filterReferenceMap = {}
617 for filterName, refFilterName
in itertools.chain([(
None, defaultFilter)],
618 filterReferenceMap.items()):
620 camFluxName = filterName +
"_camFlux" if filterName
is not None else "camFlux"
621 refFluxName = refFilterName +
"_flux"
622 if refFluxName
not in refCat.schema:
623 raise RuntimeError(
"Unknown reference filter %s" % (refFluxName,))
624 aliasMap.set(camFluxName, refFluxName)
626 refFluxErrName = refFluxName +
"Err"
627 camFluxErrName = camFluxName +
"Err"
628 aliasMap.set(camFluxErrName, refFluxErrName)
634 """This function takes in a reference catalog and creates a new catalog with additional
635 columns defined the remaining function arguments.
639 refCat : `lsst.afw.table.SimpleCatalog`
640 Reference catalog to map to new catalog
644 expandedCat : `lsst.afw.table.SimpleCatalog`
645 Deep copy of input reference catalog with additional columns added
648 mapper.addMinimalSchema(refCat.schema,
True)
649 mapper.editOutputSchema().disconnectAliases()
651 for filterName
in filterNameList:
652 mapper.editOutputSchema().addField(f
"{filterName}_flux",
654 doc=f
"flux in filter {filterName}",
657 mapper.editOutputSchema().addField(f
"{filterName}_fluxErr",
659 doc=f
"flux uncertanty in filter {filterName}",
664 mapper.editOutputSchema().addField(
"centroid_x", type=float, doReplace=
True)
665 mapper.editOutputSchema().addField(
"centroid_y", type=float, doReplace=
True)
666 mapper.editOutputSchema().addField(
"hasCentroid", type=
"Flag", doReplace=
True)
667 mapper.editOutputSchema().getAliasMap().
set(
"slot_Centroid",
"centroid")
670 mapper.editOutputSchema().addField(
"photometric",
672 doc=
"set if the object can be used for photometric"
675 mapper.editOutputSchema().addField(
"resolved",
677 doc=
"set if the object is spatially resolved"
679 mapper.editOutputSchema().addField(
"variable",
681 doc=
"set if the object has variable brightness"
685 expandedCat.setMetadata(refCat.getMetadata())
686 expandedCat.extend(refCat, mapper=mapper)
692 """Get the name of a flux field from a schema.
694 return the alias of "anyFilterMapsToThis", if present
695 else if filterName is specified:
696 return "*filterName*_camFlux" if present
697 else return "*filterName*_flux" if present (camera filter name
698 matches reference filter name)
699 else throw RuntimeError
701 return "camFlux", if present,
702 else throw RuntimeError
706 schema : `lsst.afw.table.Schema`
707 Reference catalog schema.
708 filterName : `str`, optional
709 Name of camera filter. If not specified, ``defaultFilter`` needs to be
710 set in the refcat loader config.
714 fluxFieldName : `str`
720 If an appropriate field is not found.
723 raise RuntimeError(
"schema=%s is not a schema" % (schema,))
725 return schema.getAliasMap().get(
"anyFilterMapsToThis")
730 fluxFieldList = [filterName +
"_camFlux", filterName +
"_flux"]
732 fluxFieldList = [
"camFlux"]
733 for fluxField
in fluxFieldList:
734 if fluxField
in schema:
737 raise RuntimeError(
"Could not find flux field(s) %s" % (
", ".join(fluxFieldList)))
741 """Return keys for flux and flux error.
745 schema : `lsst.afw.table.Schema`
746 Reference catalog schema.
748 Name of camera filter.
752 keys : `tuple` of (`lsst.afw.table.Key`, `lsst.afw.table.Key`)
756 - flux error key, if present, else None
761 If flux field not found.
764 fluxErrField = fluxField +
"Err"
765 fluxKey = schema[fluxField].asKey()
767 fluxErrKey = schema[fluxErrField].asKey()
770 return (fluxKey, fluxErrKey)
774 pixelMargin = pexConfig.RangeField(
775 doc=
"Padding to add to 4 all edges of the bounding box (pixels)",
780 defaultFilter = pexConfig.Field(
781 doc=(
"Default reference catalog filter to use if filter not specified in exposure;"
782 " if blank then filter must be specified in exposure."),
786 anyFilterMapsToThis = pexConfig.Field(
787 doc=(
"Always use this reference catalog filter, no matter whether or what filter name is "
788 "supplied to the loader. Effectively a trivial filterMap: map all filter names to this filter."
789 " This can be set for purely-astrometric catalogs (e.g. Gaia DR2) where there is only one "
790 "reasonable choice for every camera filter->refcat mapping, but not for refcats used for "
791 "photometry, which need a filterMap and/or colorterms/transmission corrections."),
796 filterMap = pexConfig.DictField(
797 doc=(
"Mapping of camera filter name: reference catalog filter name; "
798 "each reference filter must exist in the refcat."
799 " Note that this does not perform any bandpass corrections: it is just a lookup."),
804 requireProperMotion = pexConfig.Field(
805 doc=
"Require that the fields needed to correct proper motion "
806 "(epoch, pm_ra and pm_dec) are present?",
814 msg =
"`filterMap` and `anyFilterMapsToThis` are mutually exclusive"
815 raise pexConfig.FieldValidationError(LoadReferenceObjectsConfig.anyFilterMapsToThis,
829 r"""!Abstract base class to load objects from reference catalogs
831 @anchor LoadReferenceObjectsTask_
833 @section meas_algorithms_loadReferenceObjects_Contents Contents
835 - @ref meas_algorithms_loadReferenceObjects_Purpose
836 - @ref meas_algorithms_loadReferenceObjects_Initialize
837 - @ref meas_algorithms_loadReferenceObjects_IO
838 - @ref meas_algorithms_loadReferenceObjects_Schema
839 - @ref meas_algorithms_loadReferenceObjects_Config
841 @section meas_algorithms_loadReferenceObjects_Purpose Description
843 Abstract base class for tasks that load objects from a reference catalog
844 in a particular region of the sky.
846 Implementations must subclass this class, override the loadSkyCircle method,
847 and will typically override the value of ConfigClass with a task-specific config class.
849 @section meas_algorithms_loadReferenceObjects_Initialize Task initialisation
851 @copydoc \_\_init\_\_
853 @section meas_algorithms_loadReferenceObjects_IO Invoking the Task
855 @copydoc loadPixelBox
857 @section meas_algorithms_loadReferenceObjects_Schema Schema of the reference object catalog
859 Reference object catalogs are instances of lsst.afw.table.SimpleCatalog with the following schema
860 (other fields may also be present).
861 The units use astropy quantity conventions, so a 2 suffix means squared.
862 See also makeMinimalSchema.
863 - coord: ICRS position of star on sky (an lsst.geom.SpherePoint)
864 - centroid: position of star on an exposure, if relevant (an lsst.afw.Point2D)
865 - hasCentroid: is centroid usable? (a Flag)
866 - *referenceFilterName*_flux: brightness in the specified reference catalog filter (nJy)
867 Note: you can use astropy.units to convert from AB Magnitude to nJy:
868 `u.Magnitude(value, u.ABmag).to_value(u.nJy)`
869 - *referenceFilterName*_fluxErr (optional): brightness standard deviation (nJy);
870 omitted if no data is available; possibly nan if data is available for some objects but not others
871 - camFlux: brightness in default camera filter (nJy); omitted if defaultFilter not specified
872 - camFluxErr: brightness standard deviation for default camera filter;
873 omitted if defaultFilter not specified or standard deviation not available that filter
874 - *cameraFilterName*_camFlux: brightness in specified camera filter (nJy)
875 - *cameraFilterName*_camFluxErr (optional): brightness standard deviation
876 in specified camera filter (nJy); omitted if no data is available;
877 possibly nan if data is available for some objects but not others
878 - photometric (optional): is the object usable for photometric calibration? (a Flag)
879 - resolved (optional): is the object spatially resolved? (a Flag)
880 - variable (optional): does the object have variable brightness? (a Flag)
881 - coord_raErr: uncertainty in `coord` along the direction of right ascension (radian, an Angle)
882 = uncertainty in ra * cos(dec); nan if unknown
883 - coord_decErr: uncertainty in `coord` along the direction of declination (radian, an Angle);
886 The following are optional; fields should only be present if the
887 information is available for at least some objects.
888 Numeric values are `nan` if unknown:
889 - epoch: date of observation as TAI MJD (day)
891 - pm_ra: proper motion along the direction of right ascension (rad/year, an Angle) = dra/dt * cos(dec)
892 - pm_dec: proper motion along the direction of declination (rad/year, and Angle)
893 - pm_raErr: uncertainty in `pm_ra` (rad/year)
894 - pm_decErr: uncertainty in `pm_dec` (rad/year)
895 - pm_ra_dec_Cov: covariance between pm_ra and pm_dec (rad2/year2)
896 - pm_flag: set if proper motion, error or covariance is bad
898 - parallax: parallax (rad, an Angle)
899 - parallaxErr: uncertainty in `parallax` (rad)
900 - parallax_flag: set if parallax value or parallaxErr is bad
902 - coord_ra_pm_ra_Cov: covariance between coord_ra and pm_ra (rad2/year)
903 - coord_ra_pm_dec_Cov: covariance between coord_ra and pm_dec (rad2/year)
904 - coord_ra_parallax_Cov: covariance between coord_ra and parallax (rad2/year)
905 - coord_dec_pm_ra_Cov: covariance between coord_dec and pm_ra (rad2/year)
906 - coord_dec_pm_dec_Cov: covariance between coord_dec and pm_dec (rad2/year)
907 - coord_dec_parallax_Cov: covariance between coord_dec and parallax (rad2/year)
908 - pm_ra_parallax_Cov: covariance between pm_ra and parallax (rad2/year)
909 - pm_dec_parallax_Cov: covariance between pm_dec and parallax (rad2/year)
911 @section meas_algorithms_loadReferenceObjects_Config Configuration parameters
913 See @ref LoadReferenceObjectsConfig for a base set of configuration parameters.
914 Most subclasses will add configuration variables.
916 ConfigClass = LoadReferenceObjectsConfig
917 _DefaultName =
"LoadReferenceObjects"
920 """Construct a LoadReferenceObjectsTask
924 butler : `lsst.daf.persistence.Butler`
925 Data butler, for access reference catalogs.
927 pipeBase.Task.__init__(self, *args, **kwargs)
931 def loadPixelBox(self, bbox, wcs, filterName=None, photoCalib=None, epoch=None):
932 """Load reference objects that overlap a rectangular pixel region.
936 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D`
937 Bounding box for pixels.
938 wcs : `lsst.afw.geom.SkyWcs`
939 WCS; used to convert pixel positions to sky coordinates
942 Name of filter, or `None` or `""` for the default filter.
943 This is used for flux values in case we have flux limits
944 (which are not yet implemented).
945 photoCalib : `lsst.afw.image.PhotoCalib` (optional)
946 Calibration, or `None` if unknown.
947 epoch : `astropy.time.Time` (optional)
948 Epoch to which to correct proper motion and parallax,
949 or None to not apply such corrections.
953 results : `lsst.pipe.base.Struct`
954 A Struct containing the following fields:
955 refCat : `lsst.afw.catalog.SimpleCatalog`
956 A catalog of reference objects with the standard
957 schema, as documented in the main doc string for
958 `LoadReferenceObjects`.
959 The catalog is guaranteed to be contiguous.
961 Name of flux field for specified `filterName`.
965 The search algorithm works by searching in a region in sky
966 coordinates whose center is the center of the bbox and radius
967 is large enough to just include all 4 corners of the bbox.
968 Stars that lie outside the bbox are then trimmed from the list.
973 self.log.
info(
"Loading reference objects using center %s and radius %s deg" %
974 (circle.coord, circle.radius.asDegrees()))
975 loadRes = self.
loadSkyCircle(circle.coord, circle.radius, filterName, centroids=
True)
976 refCat = loadRes.refCat
977 numFound = len(refCat)
980 refCat = self.
_trimToBBox(refCat=refCat, bbox=circle.bbox, wcs=wcs)
981 numTrimmed = numFound - len(refCat)
982 self.log.
debug(
"trimmed %d out-of-bbox objects, leaving %d", numTrimmed, len(refCat))
983 self.log.
info(
"Loaded %d reference objects", len(refCat))
986 if not refCat.isContiguous():
987 loadRes.refCat = refCat.copy(deep=
True)
992 def loadSkyCircle(self, ctrCoord, radius, filterName=None, epoch=None, centroids=False):
993 """Load reference objects that overlap a circular sky region.
997 ctrCoord : `lsst.geom.SpherePoint`
998 ICRS center of search region.
999 radius : `lsst.geom.Angle`
1000 Radius of search region.
1001 filterName : `str` (optional)
1002 Name of filter, or `None` or `""` for the default filter.
1003 This is used for flux values in case we have flux limits
1004 (which are not yet implemented).
1005 epoch : `astropy.time.Time` (optional)
1006 Epoch to which to correct proper motion and parallax,
1007 or None to not apply such corrections.
1008 centroids : `bool` (optional)
1009 Add centroid fields to the loaded Schema. ``loadPixelBox`` expects
1010 these fields to exist.
1014 results : `lsst.pipe.base.Struct`
1015 A Struct containing the following fields:
1016 refCat : `lsst.afw.catalog.SimpleCatalog`
1017 A catalog of reference objects with the standard
1018 schema, as documented in the main doc string for
1019 `LoadReferenceObjects`.
1020 The catalog is guaranteed to be contiguous.
1022 Name of flux field for specified `filterName`.
1026 Note that subclasses are responsible for performing the proper motion
1027 correction, since this is the lowest-level interface for retrieving
1033 def _trimToBBox(refCat, bbox, wcs):
1034 """Remove objects outside a given pixel bounding box and set
1035 centroid and hasCentroid fields.
1039 refCat : `lsst.afw.table.SimpleCatalog`
1040 A catalog of objects. The schema must include fields
1041 "coord", "centroid" and "hasCentroid".
1042 The "coord" field is read.
1043 The "centroid" and "hasCentroid" fields are set.
1044 bbox : `lsst.geom.Box2D`
1046 wcs : `lsst.afw.geom.SkyWcs`
1047 WCS; used to convert sky coordinates to pixel positions.
1049 @return a catalog of reference objects in bbox, with centroid and hasCentroid fields set
1053 retStarCat =
type(refCat)(refCat.table)
1055 point = star.get(centroidKey)
1056 if bbox.contains(point):
1057 retStarCat.append(star)
1060 def _addFluxAliases(self, schema):
1061 """Add aliases for camera filter fluxes to the schema.
1063 If self.config.defaultFilter then adds these aliases:
1064 camFlux: <defaultFilter>_flux
1065 camFluxErr: <defaultFilter>_fluxErr, if the latter exists
1067 For each camFilter: refFilter in self.config.filterMap adds these aliases:
1068 <camFilter>_camFlux: <refFilter>_flux
1069 <camFilter>_camFluxErr: <refFilter>_fluxErr, if the latter exists
1073 schema : `lsst.afw.table.Schema`
1074 Schema for reference catalog.
1079 If any reference flux field is missing from the schema.
1081 aliasMap = schema.getAliasMap()
1083 if self.config.anyFilterMapsToThis
is not None:
1084 refFluxName = self.config.anyFilterMapsToThis +
"_flux"
1085 if refFluxName
not in schema:
1086 msg = f
"Unknown reference filter for anyFilterMapsToThis='{refFluxName}'"
1087 raise RuntimeError(msg)
1088 aliasMap.set(
"anyFilterMapsToThis", refFluxName)
1091 def addAliasesForOneFilter(filterName, refFilterName):
1092 """Add aliases for a single filter
1096 filterName : `str` (optional)
1097 Camera filter name. The resulting alias name is
1098 <filterName>_camFlux, or simply "camFlux" if `filterName`
1100 refFilterName : `str`
1101 Reference catalog filter name; the field
1102 <refFilterName>_flux must exist.
1104 camFluxName = filterName +
"_camFlux" if filterName
is not None else "camFlux"
1105 refFluxName = refFilterName +
"_flux"
1106 if refFluxName
not in schema:
1107 raise RuntimeError(
"Unknown reference filter %s" % (refFluxName,))
1108 aliasMap.set(camFluxName, refFluxName)
1109 refFluxErrName = refFluxName +
"Err"
1110 if refFluxErrName
in schema:
1111 camFluxErrName = camFluxName +
"Err"
1112 aliasMap.set(camFluxErrName, refFluxErrName)
1114 if self.config.defaultFilter:
1115 addAliasesForOneFilter(
None, self.config.defaultFilter)
1117 for filterName, refFilterName
in self.config.filterMap.items():
1118 addAliasesForOneFilter(filterName, refFilterName)
1122 addIsPhotometric=False, addIsResolved=False,
1123 addIsVariable=False, coordErrDim=2,
1124 addProperMotion=False, properMotionErrDim=2,
1126 """Make a standard schema for reference object catalogs.
1130 filterNameList : `list` of `str`
1131 List of filter names. Used to create <filterName>_flux fields.
1132 addIsPhotometric : `bool`
1133 If True then add field "photometric".
1134 addIsResolved : `bool`
1135 If True then add field "resolved".
1136 addIsVariable : `bool`
1137 If True then add field "variable".
1139 Number of coord error fields; must be one of 0, 2, 3:
1141 - If 2 or 3: add fields "coord_raErr" and "coord_decErr".
1142 - If 3: also add field "coord_radecErr".
1143 addProperMotion : `bool`
1144 If True add fields "epoch", "pm_ra", "pm_dec" and "pm_flag".
1145 properMotionErrDim : `int`
1146 Number of proper motion error fields; must be one of 0, 2, 3;
1147 ignored if addProperMotion false:
1148 - If 2 or 3: add fields "pm_raErr" and "pm_decErr".
1149 - If 3: also add field "pm_radecErr".
1150 addParallax : `bool`
1151 If True add fields "epoch", "parallax", "parallaxErr"
1152 and "parallax_flag".
1156 schema : `lsst.afw.table.Schema`
1157 Schema for reference catalog, an
1158 `lsst.afw.table.SimpleCatalog`.
1162 Reference catalogs support additional covariances, such as
1163 covariance between RA and proper motion in declination,
1164 that are not supported by this method, but can be added after
1165 calling this method.
1167 schema = afwTable.SimpleTable.makeMinimalSchema()
1169 afwTable.Point2DKey.addFields(
1172 "centroid on an exposure, if relevant",
1176 field=
"hasCentroid",
1178 doc=
"is position known?",
1180 for filterName
in filterNameList:
1182 field=
"%s_flux" % (filterName,),
1184 doc=
"flux in filter %s" % (filterName,),
1187 for filterName
in filterNameList:
1189 field=
"%s_fluxErr" % (filterName,),
1191 doc=
"flux uncertainty in filter %s" % (filterName,),
1194 if addIsPhotometric:
1196 field=
"photometric",
1198 doc=
"set if the object can be used for photometric calibration",
1204 doc=
"set if the object is spatially resolved",
1210 doc=
"set if the object has variable brightness",
1212 if coordErrDim
not in (0, 2, 3):
1213 raise ValueError(
"coordErrDim={}; must be (0, 2, 3)".
format(coordErrDim))
1215 afwTable.CovarianceMatrix2fKey.addFields(
1218 names=[
"ra",
"dec"],
1219 units=[
"rad",
"rad"],
1220 diagonalOnly=(coordErrDim == 2),
1223 if addProperMotion
or addParallax:
1227 doc=
"date of observation (TAI, MJD)",
1235 doc=
"proper motion in the right ascension direction = dra/dt * cos(dec)",
1241 doc=
"proper motion in the declination direction",
1244 if properMotionErrDim
not in (0, 2, 3):
1245 raise ValueError(
"properMotionErrDim={}; must be (0, 2, 3)".
format(properMotionErrDim))
1246 if properMotionErrDim > 0:
1247 afwTable.CovarianceMatrix2fKey.addFields(
1250 names=[
"ra",
"dec"],
1251 units=[
"rad/year",
"rad/year"],
1252 diagonalOnly=(properMotionErrDim == 2),
1257 doc=
"Set if proper motion or proper motion error is bad",
1268 field=
"parallaxErr",
1270 doc=
"uncertainty in parallax",
1274 field=
"parallax_flag",
1276 doc=
"Set if parallax or parallax error is bad",
1280 def _calculateCircle(self, bbox, wcs):
1281 """Compute on-sky center and radius of search region.
1285 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D`
1287 wcs : `lsst.afw.geom.SkyWcs`
1288 WCS; used to convert pixel positions to sky coordinates.
1292 results : `lsst.pipe.base.Struct`
1293 A Struct containing:
1295 - coord : `lsst.geom.SpherePoint`
1296 ICRS center of the search region.
1297 - radius : `lsst.geom.Angle`
1298 Radius of the search region.
1299 - bbox : `lsst.geom.Box2D`
1300 Bounding box used to compute the circle.
1303 bbox.grow(self.config.pixelMargin)
1304 coord = wcs.pixelToSky(bbox.getCenter())
1305 radius =
max(coord.separation(wcs.pixelToSky(pp))
for pp
in bbox.getCorners())
1306 return pipeBase.Struct(coord=coord, radius=radius, bbox=bbox)
1309 """Return metadata about the load.
1311 This metadata is used for reloading the catalog (e.g., for
1312 reconstituting a normalised match list.
1316 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D`
1318 wcs : `lsst.afw.geom.SkyWcs`
1319 WCS; used to convert pixel positions to sky coordinates.
1321 Name of camera filter, or `None` or `""` for the default
1323 photoCalib : `lsst.afw.image.PhotoCalib` (optional)
1324 Calibration, or `None` if unknown.
1325 epoch : `astropy.time.Time` (optional)
1326 Epoch to which to correct proper motion and parallax,
1327 or None to not apply such corrections.
1331 metadata : lsst.daf.base.PropertyList
1332 Metadata about the load.
1335 return self.
getMetadataCircle(circle.coord, circle.radius, filterName, photoCalib)
1338 """Return metadata about the load.
1340 This metadata is used for reloading the catalog (e.g., for
1341 reconstituting a normalised match list.
1345 coord : `lsst.geom.SpherePoint`
1346 ICRS center of the search region.
1347 radius : `lsst.geom.Angle`
1348 Radius of the search region.
1350 Name of camera filter, or `None` or `""` for the default
1352 photoCalib : `lsst.afw.image.PhotoCalib` (optional)
1353 Calibration, or `None` if unknown.
1354 epoch : `astropy.time.Time` (optional)
1355 Epoch to which to correct proper motion and parallax,
1356 or None to not apply such corrections.
1360 metadata : lsst.daf.base.PropertyList
1361 Metadata about the load
1364 md.add(
'RA', coord.getRa().asDegrees(),
'field center in degrees')
1365 md.add(
'DEC', coord.getDec().asDegrees(),
'field center in degrees')
1366 md.add(
'RADIUS', radius.asDegrees(),
'field radius in degrees, minimum')
1367 md.add(
'SMATCHV', 1,
'SourceMatchVector version number')
1368 filterName =
"UNKNOWN" if filterName
is None else str(filterName)
1369 md.add(
'FILTER', filterName,
'filter name for photometric data')
1370 md.add(
'EPOCH',
"NONE" if epoch
is None else epoch,
'Epoch (TAI MJD) for catalog')
1374 """Relink an unpersisted match list to sources and reference
1377 A match list is persisted and unpersisted as a catalog of IDs
1378 produced by afw.table.packMatches(), with match metadata
1379 (as returned by the astrometry tasks) in the catalog's metadata
1380 attribute. This method converts such a match catalog into a match
1381 list, with links to source records and reference object records.
1385 matchCat : `lsst.afw.table.BaseCatalog`
1386 Unperisted packed match list.
1387 ``matchCat.table.getMetadata()`` must contain match metadata,
1388 as returned by the astrometry tasks.
1389 sourceCat : `lsst.afw.table.SourceCatalog`
1390 Source catalog. As a side effect, the catalog will be sorted
1395 matchList : `lsst.afw.table.ReferenceMatchVector`
1401 """Apply proper motion correction to a reference catalog.
1403 Adjust position and position error in the ``catalog``
1404 for proper motion to the specified ``epoch``,
1405 modifying the catalog in place.
1409 catalog : `lsst.afw.table.SimpleCatalog`
1410 Catalog of positions, containing:
1412 - Coordinates, retrieved by the table's coordinate key.
1413 - ``coord_raErr`` : Error in Right Ascension (rad).
1414 - ``coord_decErr`` : Error in Declination (rad).
1415 - ``pm_ra`` : Proper motion in Right Ascension (rad/yr,
1417 - ``pm_raErr`` : Error in ``pm_ra`` (rad/yr), optional.
1418 - ``pm_dec`` : Proper motion in Declination (rad/yr,
1420 - ``pm_decErr`` : Error in ``pm_dec`` (rad/yr), optional.
1421 - ``epoch`` : Mean epoch of object (an astropy.time.Time)
1422 epoch : `astropy.time.Time`
1423 Epoch to which to correct proper motion,
1425 if (
"epoch" not in catalog.schema
or "pm_ra" not in catalog.schema
or "pm_dec" not in catalog.schema):
1426 if self.config.requireProperMotion:
1427 raise RuntimeError(
"Proper motion correction required but not available from catalog")
1428 self.log.
warn(
"Proper motion correction not available from catalog")
1434 """Relink an unpersisted match list to sources and reference
1437 A match list is persisted and unpersisted as a catalog of IDs
1438 produced by afw.table.packMatches(), with match metadata
1439 (as returned by the astrometry tasks) in the catalog's metadata
1440 attribute. This method converts such a match catalog into a match
1441 list, with links to source records and reference object records.
1446 Reference object loader to use in getting reference objects
1447 matchCat : `lsst.afw.table.BaseCatalog`
1448 Unperisted packed match list.
1449 ``matchCat.table.getMetadata()`` must contain match metadata,
1450 as returned by the astrometry tasks.
1451 sourceCat : `lsst.afw.table.SourceCatalog`
1452 Source catalog. As a side effect, the catalog will be sorted
1457 matchList : `lsst.afw.table.ReferenceMatchVector`
1460 matchmeta = matchCat.table.getMetadata()
1461 version = matchmeta.getInt(
'SMATCHV')
1463 raise ValueError(
'SourceMatchVector version number is %i, not 1.' % version)
1464 filterName = matchmeta.getString(
'FILTER').
strip()
1466 epoch = matchmeta.getDouble(
'EPOCH')
1469 if 'RADIUS' in matchmeta:
1472 matchmeta.getDouble(
'DEC'), lsst.geom.degrees)
1473 rad = matchmeta.getDouble(
'RADIUS') * lsst.geom.degrees
1474 refCat = refObjLoader.loadSkyCircle(ctrCoord, rad, filterName, epoch=epoch).refCat
1475 elif "INNER_UPPER_LEFT_RA" in matchmeta:
1481 for place
in (
"UPPER_LEFT",
"UPPER_RIGHT",
"LOWER_LEFT",
"LOWER_RIGHT"):
1483 matchmeta.getDouble(f
"OUTER_{place}_DEC"),
1484 lsst.geom.degrees).getVector()
1487 refCat = refObjLoader.loadRegion(outerBox, filterName=filterName, epoch=epoch).refCat
1495 """Apply proper motion correction to a reference catalog.
1497 Adjust position and position error in the ``catalog``
1498 for proper motion to the specified ``epoch``,
1499 modifying the catalog in place.
1503 log : `lsst.log.Log`
1504 Log object to write to.
1505 catalog : `lsst.afw.table.SimpleCatalog`
1506 Catalog of positions, containing:
1508 - Coordinates, retrieved by the table's coordinate key.
1509 - ``coord_raErr`` : Error in Right Ascension (rad).
1510 - ``coord_decErr`` : Error in Declination (rad).
1511 - ``pm_ra`` : Proper motion in Right Ascension (rad/yr,
1513 - ``pm_raErr`` : Error in ``pm_ra`` (rad/yr), optional.
1514 - ``pm_dec`` : Proper motion in Declination (rad/yr,
1516 - ``pm_decErr`` : Error in ``pm_dec`` (rad/yr), optional.
1517 - ``epoch`` : Mean epoch of object (an astropy.time.Time)
1518 epoch : `astropy.time.Time`
1519 Epoch to which to correct proper motion.
1521 if "epoch" not in catalog.schema
or "pm_ra" not in catalog.schema
or "pm_dec" not in catalog.schema:
1522 log.warn(
"Proper motion correction not available from catalog")
1524 if not catalog.isContiguous():
1525 raise RuntimeError(
"Catalog must be contiguous")
1526 catEpoch = astropy.time.Time(catalog[
"epoch"], scale=
"tai", format=
"mjd")
1527 log.info(
"Correcting reference catalog for proper motion to %r", epoch)
1529 timeDiffsYears = (epoch.tai - catEpoch).
to(astropy.units.yr).value
1530 coordKey = catalog.table.getCoordKey()
1533 pmRaRad = catalog[
"pm_ra"]
1534 pmDecRad = catalog[
"pm_dec"]
1535 offsetsRaRad = pmRaRad*timeDiffsYears
1536 offsetsDecRad = pmDecRad*timeDiffsYears
1544 offsetBearingsRad = numpy.arctan2(pmDecRad*1e6, pmRaRad*1e6)
1545 offsetAmountsRad = numpy.hypot(offsetsRaRad, offsetsDecRad)
1546 for record, bearingRad, amountRad
in zip(catalog, offsetBearingsRad, offsetAmountsRad):
1547 record.set(coordKey,
1548 record.get(coordKey).offset(bearing=bearingRad*lsst.geom.radians,
1549 amount=amountRad*lsst.geom.radians))
1551 if "coord_raErr" in catalog.schema:
1552 catalog[
"coord_raErr"] = numpy.hypot(catalog[
"coord_raErr"],
1553 catalog[
"pm_raErr"]*timeDiffsYears)
1554 if "coord_decErr" in catalog.schema:
1555 catalog[
"coord_decErr"] = numpy.hypot(catalog[
"coord_decErr"],
1556 catalog[
"pm_decErr"]*timeDiffsYears)