LSSTApplications  17.0+124,17.0+14,17.0+73,18.0.0+37,18.0.0+80,18.0.0-4-g68ffd23+4,18.1.0-1-g0001055+12,18.1.0-1-g03d53ef+5,18.1.0-1-g1349e88+55,18.1.0-1-g2505f39+44,18.1.0-1-g5315e5e+4,18.1.0-1-g5e4b7ea+14,18.1.0-1-g7e8fceb+4,18.1.0-1-g85f8cd4+48,18.1.0-1-g8ff0b9f+4,18.1.0-1-ga2c679d+1,18.1.0-1-gd55f500+35,18.1.0-10-gb58edde+2,18.1.0-11-g0997b02+4,18.1.0-13-gfe4edf0b+12,18.1.0-14-g259bd21+21,18.1.0-19-gdb69f3f+2,18.1.0-2-g5f9922c+24,18.1.0-2-gd3b74e5+11,18.1.0-2-gfbf3545+32,18.1.0-26-g728bddb4+5,18.1.0-27-g6ff7ca9+2,18.1.0-3-g52aa583+25,18.1.0-3-g8ea57af+9,18.1.0-3-gb69f684+42,18.1.0-3-gfcaddf3+6,18.1.0-32-gd8786685a,18.1.0-4-gf3f9b77+6,18.1.0-5-g1dd662b+2,18.1.0-5-g6dbcb01+41,18.1.0-6-gae77429+3,18.1.0-7-g9d75d83+9,18.1.0-7-gae09a6d+30,18.1.0-9-gc381ef5+4,w.2019.45
LSSTDataManagementBasePackage
Classes | Functions
lsst.meas.algorithms.loadReferenceObjects Namespace Reference

Classes

class  _FilterCatalog
 
class  LoadReferenceObjectsConfig
 
class  LoadReferenceObjectsTask
 Abstract base class to load objects from reference catalogs. More...
 
class  ReferenceObjectLoader
 

Functions

def isOldFluxField (name, units)
 
def hasNanojanskyFluxUnits (schema)
 
def getFormatVersionFromRefCat (refCat)
 
def convertToNanojansky (catalog, log, doConvert=True)
 
def getRefFluxField (schema, filterName=None)
 
def getRefFluxKeys (schema, filterName=None)
 
def joinMatchListWithCatalogImpl (refObjLoader, matchCat, sourceCat)
 
def applyProperMotionsImpl (log, catalog, epoch)
 

Function Documentation

◆ applyProperMotionsImpl()

def lsst.meas.algorithms.loadReferenceObjects.applyProperMotionsImpl (   log,
  catalog,
  epoch 
)
Apply proper motion correction to a reference catalog.

Adjust position and position error in the ``catalog``
for proper motion to the specified ``epoch``,
modifying the catalong in place.

Parameters
----------
log : `lsst.log.log`
    log object to write to
catalog : `lsst.afw.table.SimpleCatalog`
    Catalog of positions, containing:

    - Coordinates, retrieved by the table's coordinate key.
    - ``coord_raErr`` : Error in Right Ascension (rad).
    - ``coord_decErr`` : Error in Declination (rad).
    - ``pm_ra`` : Proper motion in Right Ascension (rad/yr,
        East positive)
    - ``pm_raErr`` : Error in ``pm_ra`` (rad/yr), optional.
    - ``pm_dec`` : Proper motion in Declination (rad/yr,
        North positive)
    - ``pm_decErr`` : Error in ``pm_dec`` (rad/yr), optional.
    - ``epoch`` : Mean epoch of object (an astropy.time.Time)
epoch : `astropy.time.Time` (optional)
    Epoch to which to correct proper motion and parallax,
    or None to not apply such corrections.

Definition at line 1456 of file loadReferenceObjects.py.

1456 def applyProperMotionsImpl(log, catalog, epoch):
1457  """Apply proper motion correction to a reference catalog.
1458 
1459  Adjust position and position error in the ``catalog``
1460  for proper motion to the specified ``epoch``,
1461  modifying the catalong in place.
1462 
1463  Parameters
1464  ----------
1465  log : `lsst.log.log`
1466  log object to write to
1467  catalog : `lsst.afw.table.SimpleCatalog`
1468  Catalog of positions, containing:
1469 
1470  - Coordinates, retrieved by the table's coordinate key.
1471  - ``coord_raErr`` : Error in Right Ascension (rad).
1472  - ``coord_decErr`` : Error in Declination (rad).
1473  - ``pm_ra`` : Proper motion in Right Ascension (rad/yr,
1474  East positive)
1475  - ``pm_raErr`` : Error in ``pm_ra`` (rad/yr), optional.
1476  - ``pm_dec`` : Proper motion in Declination (rad/yr,
1477  North positive)
1478  - ``pm_decErr`` : Error in ``pm_dec`` (rad/yr), optional.
1479  - ``epoch`` : Mean epoch of object (an astropy.time.Time)
1480  epoch : `astropy.time.Time` (optional)
1481  Epoch to which to correct proper motion and parallax,
1482  or None to not apply such corrections.
1483  """
1484  if "epoch" not in catalog.schema or "pm_ra" not in catalog.schema or "pm_dec" not in catalog.schema:
1485  log.warn("Proper motion correction not available from catalog")
1486  return
1487  if not catalog.isContiguous():
1488  raise RuntimeError("Catalog must be contiguous")
1489  catEpoch = astropy.time.Time(catalog["epoch"], scale="tai", format="mjd")
1490  log.debug("Correcting reference catalog for proper motion to %r", epoch)
1491  # Use `epoch.tai` to make sure the time difference is in TAI
1492  timeDiffsYears = (epoch.tai - catEpoch).to(astropy.units.yr).value
1493  coordKey = catalog.table.getCoordKey()
1494  # Compute the offset of each object due to proper motion
1495  # as components of the arc of a great circle along RA and Dec
1496  pmRaRad = catalog["pm_ra"]
1497  pmDecRad = catalog["pm_dec"]
1498  offsetsRaRad = pmRaRad*timeDiffsYears
1499  offsetsDecRad = pmDecRad*timeDiffsYears
1500  # Compute the corresponding bearing and arc length of each offset
1501  # due to proper motion, and apply the offset
1502  # The factor of 1e6 for computing bearing is intended as
1503  # a reasonable scale for typical values of proper motion
1504  # in order to avoid large errors for small values of proper motion;
1505  # using the offsets is another option, but it can give
1506  # needlessly large errors for short duration
1507  offsetBearingsRad = numpy.arctan2(pmDecRad*1e6, pmRaRad*1e6)
1508  offsetAmountsRad = numpy.hypot(offsetsRaRad, offsetsDecRad)
1509  for record, bearingRad, amountRad in zip(catalog, offsetBearingsRad, offsetAmountsRad):
1510  record.set(coordKey,
1511  record.get(coordKey).offset(bearing=bearingRad*lsst.geom.radians,
1512  amount=amountRad*lsst.geom.radians))
1513  # Increase error in RA and Dec based on error in proper motion
1514  if "coord_raErr" in catalog.schema:
1515  catalog["coord_raErr"] = numpy.hypot(catalog["coord_raErr"],
1516  catalog["pm_raErr"]*timeDiffsYears)
1517  if "coord_decErr" in catalog.schema:
1518  catalog["coord_decErr"] = numpy.hypot(catalog["coord_decErr"],
1519  catalog["pm_decErr"]*timeDiffsYears)
1520 
table::Key< int > to

◆ convertToNanojansky()

def lsst.meas.algorithms.loadReferenceObjects.convertToNanojansky (   catalog,
  log,
  doConvert = True 
)
Convert fluxes in a catalog from jansky to nanojansky.

Parameters
----------
catalog : `lsst.afw.table.SimpleCatalog`
    The catalog to convert.
log : `lsst.log.Log`
    Log to send messages to.
doConvert : `bool`, optional
    Return a converted catalog, or just identify the fields that need to be converted?
    This supports the "write=False" mode of `bin/convert_to_nJy.py`.

Returns
-------
catalog : `lsst.afw.table.SimpleCatalog` or None
    The converted catalog, or None if ``doConvert`` is False.

Notes
-----
Support for old units in reference catalogs will be removed after the
release of late calendar year 2019.
Use `meas_algorithms/bin/convert_to_nJy.py` to update your reference catalog.

Definition at line 89 of file loadReferenceObjects.py.

89 def convertToNanojansky(catalog, log, doConvert=True):
90  """Convert fluxes in a catalog from jansky to nanojansky.
91 
92  Parameters
93  ----------
94  catalog : `lsst.afw.table.SimpleCatalog`
95  The catalog to convert.
96  log : `lsst.log.Log`
97  Log to send messages to.
98  doConvert : `bool`, optional
99  Return a converted catalog, or just identify the fields that need to be converted?
100  This supports the "write=False" mode of `bin/convert_to_nJy.py`.
101 
102  Returns
103  -------
104  catalog : `lsst.afw.table.SimpleCatalog` or None
105  The converted catalog, or None if ``doConvert`` is False.
106 
107  Notes
108  -----
109  Support for old units in reference catalogs will be removed after the
110  release of late calendar year 2019.
111  Use `meas_algorithms/bin/convert_to_nJy.py` to update your reference catalog.
112  """
113  # Do not share the AliasMap: for refcats, that gets created when the
114  # catalog is read from disk and should not be propagated.
115  mapper = lsst.afw.table.SchemaMapper(catalog.schema, shareAliasMap=False)
116  mapper.addMinimalSchema(lsst.afw.table.SimpleTable.makeMinimalSchema())
117  input_fields = []
118  output_fields = []
119  for field in catalog.schema:
120  oldName = field.field.getName()
121  oldUnits = field.field.getUnits()
122  if isOldFluxField(oldName, oldUnits):
123  units = 'nJy'
124  # remap Sigma flux fields to Err, so we can drop the alias
125  if oldName.endswith('_fluxSigma'):
126  name = oldName.replace('_fluxSigma', '_fluxErr')
127  else:
128  name = oldName
129  newField = lsst.afw.table.Field[field.dtype](name, field.field.getDoc(), units)
130  mapper.addMapping(field.getKey(), newField)
131  input_fields.append(field.field)
132  output_fields.append(newField)
133  else:
134  mapper.addMapping(field.getKey())
135 
136  fluxFieldsStr = '; '.join("(%s, '%s')" % (field.getName(), field.getUnits()) for field in input_fields)
137 
138  if doConvert:
139  newSchema = mapper.getOutputSchema()
140  output = lsst.afw.table.SimpleCatalog(newSchema)
141  output.extend(catalog, mapper=mapper)
142  for field in output_fields:
143  output[field.getName()] *= 1e9
144  log.info(f"Converted refcat flux fields to nJy (name, units): {fluxFieldsStr}")
145  return output
146  else:
147  log.info(f"Found old-style refcat flux fields (name, units): {fluxFieldsStr}")
148  return None
149 
150 
static Schema makeMinimalSchema()
Return a minimal schema for Simple tables and records.
Definition: Simple.h:140
A mapping between the keys of two Schemas, used to copy data between them.
Definition: SchemaMapper.h:21
def convertToNanojansky(catalog, log, doConvert=True)
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
A description of a field in a table.
Definition: Field.h:24

◆ getFormatVersionFromRefCat()

def lsst.meas.algorithms.loadReferenceObjects.getFormatVersionFromRefCat (   refCat)
"Return the format version stored in a reference catalog header.

Parameters
----------
refCat : `lsst.afw.table.SimpleCatalog`
    Reference catalog to inspect.

Returns
-------
version : `int` or `None`
    Format version integer, or `None` if the catalog has no metadata
    or the metadata does not include a "REFCAT_FORMAT_VERSION" key.

Definition at line 66 of file loadReferenceObjects.py.

66 def getFormatVersionFromRefCat(refCat):
67  """"Return the format version stored in a reference catalog header.
68 
69  Parameters
70  ----------
71  refCat : `lsst.afw.table.SimpleCatalog`
72  Reference catalog to inspect.
73 
74  Returns
75  -------
76  version : `int` or `None`
77  Format version integer, or `None` if the catalog has no metadata
78  or the metadata does not include a "REFCAT_FORMAT_VERSION" key.
79  """
80  md = refCat.getMetadata()
81  if md is None:
82  return None
83  try:
84  return md.getScalar("REFCAT_FORMAT_VERSION")
85  except KeyError:
86  return None
87 
88 

◆ getRefFluxField()

def lsst.meas.algorithms.loadReferenceObjects.getRefFluxField (   schema,
  filterName = None 
)
Get the name of a flux field from a schema.

if filterName is specified:
    return *filterName*_camFlux if present
    else return *filterName*_flux if present (camera filter name matches reference filter name)
    else throw RuntimeError
else:
    return camFlux, if present,
    else throw RuntimeError

Parameters
----------
schema : `lsst.afw.table.Schema`
    Reference catalog schema.
filterName : `str`
    Name of camera filter.

Returns
-------
fluxFieldName : `str`
    Name of flux field.

Raises
------
RuntimeError
    If an appropriate field is not found.

Definition at line 687 of file loadReferenceObjects.py.

687 def getRefFluxField(schema, filterName=None):
688  """Get the name of a flux field from a schema.
689 
690  if filterName is specified:
691  return *filterName*_camFlux if present
692  else return *filterName*_flux if present (camera filter name matches reference filter name)
693  else throw RuntimeError
694  else:
695  return camFlux, if present,
696  else throw RuntimeError
697 
698  Parameters
699  ----------
700  schema : `lsst.afw.table.Schema`
701  Reference catalog schema.
702  filterName : `str`
703  Name of camera filter.
704 
705  Returns
706  -------
707  fluxFieldName : `str`
708  Name of flux field.
709 
710  Raises
711  ------
712  RuntimeError
713  If an appropriate field is not found.
714  """
715  if not isinstance(schema, afwTable.Schema):
716  raise RuntimeError("schema=%s is not a schema" % (schema,))
717  if filterName:
718  fluxFieldList = [filterName + "_camFlux", filterName + "_flux"]
719  else:
720  fluxFieldList = ["camFlux"]
721  for fluxField in fluxFieldList:
722  if fluxField in schema:
723  return fluxField
724 
725  raise RuntimeError("Could not find flux field(s) %s" % (", ".join(fluxFieldList)))
726 
727 
Defines the fields and offsets for a table.
Definition: Schema.h:50

◆ getRefFluxKeys()

def lsst.meas.algorithms.loadReferenceObjects.getRefFluxKeys (   schema,
  filterName = None 
)
Return keys for flux and flux error.

Parameters
----------
schema : `lsst.afw.table.Schema`
    Reference catalog schema.
filterName : `str`
    Name of camera filter.

Returns
-------
keys : `tuple` of (`lsst.afw.table.Key`, `lsst.afw.table.Key`)
    Two keys:

    - flux key
    - flux error key, if present, else None

Raises
------
RuntimeError
    If flux field not found.

Definition at line 728 of file loadReferenceObjects.py.

728 def getRefFluxKeys(schema, filterName=None):
729  """Return keys for flux and flux error.
730 
731  Parameters
732  ----------
733  schema : `lsst.afw.table.Schema`
734  Reference catalog schema.
735  filterName : `str`
736  Name of camera filter.
737 
738  Returns
739  -------
740  keys : `tuple` of (`lsst.afw.table.Key`, `lsst.afw.table.Key`)
741  Two keys:
742 
743  - flux key
744  - flux error key, if present, else None
745 
746  Raises
747  ------
748  RuntimeError
749  If flux field not found.
750  """
751  fluxField = getRefFluxField(schema, filterName)
752  fluxErrField = fluxField + "Err"
753  fluxKey = schema[fluxField].asKey()
754  try:
755  fluxErrKey = schema[fluxErrField].asKey()
756  except Exception:
757  fluxErrKey = None
758  return (fluxKey, fluxErrKey)
759 
760 

◆ hasNanojanskyFluxUnits()

def lsst.meas.algorithms.loadReferenceObjects.hasNanojanskyFluxUnits (   schema)
Return True if the units of all flux and fluxErr are correct (nJy).

Definition at line 57 of file loadReferenceObjects.py.

57 def hasNanojanskyFluxUnits(schema):
58  """Return True if the units of all flux and fluxErr are correct (nJy).
59  """
60  for field in schema:
61  if isOldFluxField(field.field.getName(), field.field.getUnits()):
62  return False
63  return True
64 
65 

◆ isOldFluxField()

def lsst.meas.algorithms.loadReferenceObjects.isOldFluxField (   name,
  units 
)
Return True if this name/units combination corresponds to an
"old-style" reference catalog flux field.

Definition at line 46 of file loadReferenceObjects.py.

46 def isOldFluxField(name, units):
47  """Return True if this name/units combination corresponds to an
48  "old-style" reference catalog flux field.
49  """
50  unitsCheck = units != 'nJy' # (units == 'Jy' or units == '' or units == '?')
51  isFlux = name.endswith('_flux')
52  isFluxSigma = name.endswith('_fluxSigma')
53  isFluxErr = name.endswith('_fluxErr')
54  return (isFlux or isFluxSigma or isFluxErr) and unitsCheck
55 
56 

◆ joinMatchListWithCatalogImpl()

def lsst.meas.algorithms.loadReferenceObjects.joinMatchListWithCatalogImpl (   refObjLoader,
  matchCat,
  sourceCat 
)
Relink an unpersisted match list to sources and reference
objects.

A match list is persisted and unpersisted as a catalog of IDs
produced by afw.table.packMatches(), with match metadata
(as returned by the astrometry tasks) in the catalog's metadata
attribute. This method converts such a match catalog into a match
list, with links to source records and reference object records.

Parameters
----------
refObjLoader
    Reference object loader to use in getting reference objects
matchCat : `lsst.afw.table.BaseCatalog`
    Unperisted packed match list.
    ``matchCat.table.getMetadata()`` must contain match metadata,
    as returned by the astrometry tasks.
sourceCat : `lsst.afw.table.SourceCatalog`
    Source catalog. As a side effect, the catalog will be sorted
    by ID.

Returns
-------
matchList : `lsst.afw.table.ReferenceMatchVector`
    Match list.

Definition at line 1395 of file loadReferenceObjects.py.

1395 def joinMatchListWithCatalogImpl(refObjLoader, matchCat, sourceCat):
1396  """Relink an unpersisted match list to sources and reference
1397  objects.
1398 
1399  A match list is persisted and unpersisted as a catalog of IDs
1400  produced by afw.table.packMatches(), with match metadata
1401  (as returned by the astrometry tasks) in the catalog's metadata
1402  attribute. This method converts such a match catalog into a match
1403  list, with links to source records and reference object records.
1404 
1405  Parameters
1406  ----------
1407  refObjLoader
1408  Reference object loader to use in getting reference objects
1409  matchCat : `lsst.afw.table.BaseCatalog`
1410  Unperisted packed match list.
1411  ``matchCat.table.getMetadata()`` must contain match metadata,
1412  as returned by the astrometry tasks.
1413  sourceCat : `lsst.afw.table.SourceCatalog`
1414  Source catalog. As a side effect, the catalog will be sorted
1415  by ID.
1416 
1417  Returns
1418  -------
1419  matchList : `lsst.afw.table.ReferenceMatchVector`
1420  Match list.
1421  """
1422  matchmeta = matchCat.table.getMetadata()
1423  version = matchmeta.getInt('SMATCHV')
1424  if version != 1:
1425  raise ValueError('SourceMatchVector version number is %i, not 1.' % version)
1426  filterName = matchmeta.getString('FILTER').strip()
1427  try:
1428  epoch = matchmeta.getDouble('EPOCH')
1430  epoch = None # Not present, or not correct type means it's not set
1431  if 'RADIUS' in matchmeta:
1432  # This is a circle style metadata, call loadSkyCircle
1433  ctrCoord = lsst.geom.SpherePoint(matchmeta.getDouble('RA'),
1434  matchmeta.getDouble('DEC'), lsst.geom.degrees)
1435  rad = matchmeta.getDouble('RADIUS') * lsst.geom.degrees
1436  refCat = refObjLoader.loadSkyCircle(ctrCoord, rad, filterName, epoch=epoch).refCat
1437  elif "INNER_UPPER_LEFT_RA" in matchmeta:
1438  # This is the sky box type (only triggers in the LoadReferenceObject class, not task)
1439  # Only the outer box is required to be loaded to get the maximum region, all filtering
1440  # will be done by the unpackMatches function, and no spatial filtering needs to be done
1441  # by the refObjLoader
1442  box = []
1443  for place in ("UPPER_LEFT", "UPPER_RIGHT", "LOWER_LEFT", "LOWER_RIGHT"):
1444  coord = lsst.geom.SpherePoint(matchmeta.getDouble(f"OUTER_{place}_RA"),
1445  matchmeta.getDouble(f"OUTER_{place}_DEC"),
1446  lsst.geom.degrees).getVector()
1447  box.append(coord)
1448  outerBox = sphgeom.ConvexPolygon(box)
1449  refCat = refObjLoader.loadRegion(outerBox, filterName=filterName, epoch=epoch).refCat
1450 
1451  refCat.sort()
1452  sourceCat.sort()
1453  return afwTable.unpackMatches(matchCat, refCat, sourceCat)
1454 
1455 
Reports attempts to access elements using an invalid key.
Definition: Runtime.h:151
def joinMatchListWithCatalogImpl(refObjLoader, matchCat, sourceCat)
Point in an unspecified spherical coordinate system.
Definition: SpherePoint.h:57
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
bool strip
Definition: fits.cc:901