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