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