Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

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 

27import abc 

28import itertools 

29 

30import astropy.time 

31import astropy.units 

32import numpy 

33 

34import lsst.geom as geom 

35import lsst.afw.table as afwTable 

36import lsst.pex.config as pexConfig 

37import lsst.pipe.base as pipeBase 

38import lsst.log 

39from lsst import sphgeom 

40from lsst.daf.base import PropertyList 

41 

42 

43def isOldFluxField(name, units): 

44 """Return True if this name/units combination corresponds to an 

45 "old-style" reference catalog flux field. 

46 """ 

47 unitsCheck = units != 'nJy' # (units == 'Jy' or units == '' or units == '?') 

48 isFlux = name.endswith('_flux') 

49 isFluxSigma = name.endswith('_fluxSigma') 

50 isFluxErr = name.endswith('_fluxErr') 

51 return (isFlux or isFluxSigma or isFluxErr) and unitsCheck 

52 

53 

54def hasNanojanskyFluxUnits(schema): 

55 """Return True if the units of all flux and fluxErr are correct (nJy). 

56 """ 

57 for field in schema: 

58 if isOldFluxField(field.field.getName(), field.field.getUnits()): 58 ↛ 59line 58 didn't jump to line 59, because the condition on line 58 was never true

59 return False 

60 return True 

61 

62 

63def getFormatVersionFromRefCat(refCat): 

64 """"Return the format version stored in a reference catalog header. 

65 

66 Parameters 

67 ---------- 

68 refCat : `lsst.afw.table.SimpleCatalog` 

69 Reference catalog to inspect. 

70 

71 Returns 

72 ------- 

73 version : `int` or `None` 

74 Format version integer, or `None` if the catalog has no metadata 

75 or the metadata does not include a "REFCAT_FORMAT_VERSION" key. 

76 """ 

77 md = refCat.getMetadata() 

78 if md is None: 78 ↛ 79line 78 didn't jump to line 79, because the condition on line 78 was never true

79 return None 

80 try: 

81 return md.getScalar("REFCAT_FORMAT_VERSION") 

82 except KeyError: 

83 return None 

84 

85 

86def convertToNanojansky(catalog, log, doConvert=True): 

87 """Convert fluxes in a catalog from jansky to nanojansky. 

88 

89 Parameters 

90 ---------- 

91 catalog : `lsst.afw.table.SimpleCatalog` 

92 The catalog to convert. 

93 log : `lsst.log.Log` 

94 Log to send messages to. 

95 doConvert : `bool`, optional 

96 Return a converted catalog, or just identify the fields that need to be converted? 

97 This supports the "write=False" mode of `bin/convert_to_nJy.py`. 

98 

99 Returns 

100 ------- 

101 catalog : `lsst.afw.table.SimpleCatalog` or None 

102 The converted catalog, or None if ``doConvert`` is False. 

103 

104 Notes 

105 ----- 

106 Support for old units in reference catalogs will be removed after the 

107 release of late calendar year 2019. 

108 Use `meas_algorithms/bin/convert_to_nJy.py` to update your reference catalog. 

109 """ 

110 # Do not share the AliasMap: for refcats, that gets created when the 

111 # catalog is read from disk and should not be propagated. 

112 mapper = lsst.afw.table.SchemaMapper(catalog.schema, shareAliasMap=False) 

113 mapper.addMinimalSchema(lsst.afw.table.SimpleTable.makeMinimalSchema()) 

114 input_fields = [] 

115 output_fields = [] 

116 for field in catalog.schema: 

117 oldName = field.field.getName() 

118 oldUnits = field.field.getUnits() 

119 if isOldFluxField(oldName, oldUnits): 

120 units = 'nJy' 

121 # remap Sigma flux fields to Err, so we can drop the alias 

122 if oldName.endswith('_fluxSigma'): 

123 name = oldName.replace('_fluxSigma', '_fluxErr') 

124 else: 

125 name = oldName 

126 newField = lsst.afw.table.Field[field.dtype](name, field.field.getDoc(), units) 

127 mapper.addMapping(field.getKey(), newField) 

128 input_fields.append(field.field) 

129 output_fields.append(newField) 

130 else: 

131 mapper.addMapping(field.getKey()) 

132 

133 fluxFieldsStr = '; '.join("(%s, '%s')" % (field.getName(), field.getUnits()) for field in input_fields) 

134 

135 if doConvert: 

136 newSchema = mapper.getOutputSchema() 

137 output = lsst.afw.table.SimpleCatalog(newSchema) 

138 output.extend(catalog, mapper=mapper) 

139 for field in output_fields: 

140 output[field.getName()] *= 1e9 

141 log.info(f"Converted refcat flux fields to nJy (name, units): {fluxFieldsStr}") 

142 return output 

143 else: 

144 log.info(f"Found old-style refcat flux fields (name, units): {fluxFieldsStr}") 

145 return None 

146 

147 

148class _FilterCatalog: 

149 """This is a private helper class which filters catalogs by 

150 row based on the row being inside the region used to initialize 

151 the class. 

152 

153 Parameters 

154 ---------- 

155 region : `lsst.sphgeom.Region` 

156 The spatial region which all objects should lie within 

157 """ 

158 def __init__(self, region): 

159 self.region = region 

160 

161 def __call__(self, refCat, catRegion): 

162 """This call method on an instance of this class takes in a reference 

163 catalog, and the region from which the catalog was generated. 

164 

165 If the catalog region is entirely contained within the region used to 

166 initialize this class, then all the entries in the catalog must be 

167 within the region and so the whole catalog is returned. 

168 

169 If the catalog region is not entirely contained, then the location for 

170 each record is tested against the region used to initialize the class. 

171 Records which fall inside this region are added to a new catalog, and 

172 this catalog is then returned. 

173 

174 Parameters 

175 --------- 

176 refCat : `lsst.afw.table.SourceCatalog` 

177 SourceCatalog to be filtered. 

178 catRegion : `lsst.sphgeom.Region` 

179 Region in which the catalog was created 

180 """ 

181 if catRegion.isWithin(self.region): 

182 # no filtering needed, region completely contains refcat 

183 return refCat 

184 

185 filteredRefCat = type(refCat)(refCat.table) 

186 for record in refCat: 

187 if self.region.contains(record.getCoord().getVector()): 

188 filteredRefCat.append(record) 

189 return filteredRefCat 

190 

191 

192class ReferenceObjectLoaderBase: 

193 """Base class for reference object loaders, to facilitate gen2/gen3 code 

194 sharing. 

195 """ 

196 def applyProperMotions(self, catalog, epoch): 

197 """Apply proper motion correction to a reference catalog. 

198 

199 Adjust position and position error in the ``catalog`` 

200 for proper motion to the specified ``epoch``, 

201 modifying the catalog in place. 

202 

203 Parameters 

204 ---------- 

205 catalog : `lsst.afw.table.SimpleCatalog` 

206 Catalog of positions, containing at least these fields: 

207 

208 - Coordinates, retrieved by the table's coordinate key. 

209 - ``coord_raErr`` : Error in Right Ascension (rad). 

210 - ``coord_decErr`` : Error in Declination (rad). 

211 - ``pm_ra`` : Proper motion in Right Ascension (rad/yr, 

212 East positive) 

213 - ``pm_raErr`` : Error in ``pm_ra`` (rad/yr), optional. 

214 - ``pm_dec`` : Proper motion in Declination (rad/yr, 

215 North positive) 

216 - ``pm_decErr`` : Error in ``pm_dec`` (rad/yr), optional. 

217 - ``epoch`` : Mean epoch of object (an astropy.time.Time) 

218 epoch : `astropy.time.Time` 

219 Epoch to which to correct proper motion. 

220 If None, do not apply PM corrections or raise if 

221 ``config.requireProperMotion`` is True. 

222 

223 Raises 

224 ------ 

225 RuntimeError 

226 Raised if ``config.requireProperMotion`` is set but we cannot 

227 apply the proper motion correction for some reason. 

228 """ 

229 if epoch is None: 229 ↛ 237line 229 didn't jump to line 237, because the condition on line 229 was never false

230 if self.config.requireProperMotion: 230 ↛ 231line 230 didn't jump to line 231, because the condition on line 230 was never true

231 raise RuntimeError("requireProperMotion=True but epoch not provided to loader.") 

232 else: 

233 self.log.debug("No epoch provided: not applying proper motion corrections to refcat.") 

234 return 

235 

236 # Warn/raise for a catalog in an incorrect format, if epoch was specified. 

237 if ("pm_ra" in catalog.schema 

238 and not isinstance(catalog.schema["pm_ra"].asKey(), lsst.afw.table.KeyAngle)): 

239 if self.config.requireProperMotion: 

240 raise RuntimeError("requireProperMotion=True but refcat pm_ra field is not an Angle.") 

241 else: 

242 self.log.warn("Reference catalog pm_ra field is not an Angle; cannot apply proper motion.") 

243 return 

244 

245 if ("epoch" not in catalog.schema or "pm_ra" not in catalog.schema): 

246 if self.config.requireProperMotion: 

247 raise RuntimeError("requireProperMotion=True but PM data not available from catalog.") 

248 else: 

249 self.log.warn("Proper motion correction not available for this reference catalog.") 

250 return 

251 

252 applyProperMotionsImpl(self.log, catalog, epoch) 

253 

254 

255class ReferenceObjectLoader(ReferenceObjectLoaderBase): 

256 """This class facilitates loading reference catalogs with gen 3 middleware 

257 

258 The middleware preflight solver will create a list of datarefs that may 

259 possibly overlap a given region. These datarefs are then used to construct 

260 and instance of this class. The class instance should then be passed into 

261 a task which needs reference catalogs. These tasks should then determine 

262 the exact region of the sky reference catalogs will be loaded for, and 

263 call a corresponding method to load the reference objects. 

264 """ 

265 def __init__(self, dataIds, refCats, config, log=None): 

266 """ Constructs an instance of ReferenceObjectLoader 

267 

268 Parameters 

269 ---------- 

270 dataIds : iterable of `lsst.daf.butler.DataIds` 

271 An iterable object of DataSetRefs which point to reference catalogs 

272 in a gen 3 repository. 

273 refCats : iterable of `lsst.daf.butler.DeferedDatasetHandle` 

274 Handles to load refCats on demand 

275 config : `lsst.pex.config.configurableField` 

276 Configuration for the loader. 

277 log : `lsst.log.Log` or `None`, optional 

278 Logger object used to write out messages. If `None` the default 

279 lsst logger will be used. 

280 """ 

281 self.dataIds = dataIds 

282 self.refCats = refCats 

283 self.log = log or lsst.log.Log.getDefaultLogger() 

284 self.config = config 

285 

286 @staticmethod 

287 def _makeBoxRegion(BBox, wcs, BBoxPadding): 

288 outerLocalBBox = geom.Box2D(BBox) 

289 innerLocalBBox = geom.Box2D(BBox) 

290 

291 # Grow the bounding box to allow for effects not fully captured by the 

292 # wcs provided (which represents the current best-guess wcs solution 

293 # associated with the dataset for which the calibration is to be 

294 # computed using the loaded and trimmed reference catalog being defined 

295 # here). These effects could include pointing errors and/or an 

296 # insufficient optical distorition model for the instrument. The idea 

297 # is to ensure the spherical geometric region created contains the 

298 # entire region covered by the bbox. 

299 # Also create an inner region that is sure to be inside the bbox. 

300 outerLocalBBox.grow(BBoxPadding) 

301 innerLocalBBox.grow(-1*BBoxPadding) 

302 

303 # Handle the case where the inner bounding box shrank to a zero sized 

304 # region (which will be the case if the shrunken size of either 

305 # dimension is less than or equal to zero). In this case, the inner 

306 # bounding box is set to the original input bounding box. This is 

307 # probably not the best way to handle an empty inner bounding box, but 

308 # it is what the calling code currently expects. 

309 if innerLocalBBox.getDimensions() == geom.Extent2D(0, 0): 

310 innerLocalBBox = geom.Box2D(BBox) 

311 

312 # Convert the corners of the bounding boxes to sky coordinates. 

313 innerBoxCorners = innerLocalBBox.getCorners() 

314 innerSphCorners = [wcs.pixelToSky(corner).getVector() for corner in innerBoxCorners] 

315 innerSkyRegion = sphgeom.ConvexPolygon(innerSphCorners) 

316 

317 outerBoxCorners = outerLocalBBox.getCorners() 

318 outerSphCorners = [wcs.pixelToSky(corner).getVector() for corner in outerBoxCorners] 

319 outerSkyRegion = sphgeom.ConvexPolygon(outerSphCorners) 

320 

321 return innerSkyRegion, outerSkyRegion, innerSphCorners, outerSphCorners 

322 

323 def loadPixelBox(self, bbox, wcs, filterName=None, epoch=None, photoCalib=None, 

324 bboxToSpherePadding=100): 

325 """Load reference objects that are within a pixel-based rectangular 

326 region. 

327 

328 This algorithm works by creating a spherical box whose corners 

329 correspond to the WCS converted corners of the input bounding box 

330 (possibly padded). It then defines a filtering function which looks at 

331 the pixel position of the reference objects and accepts only those that 

332 lie within the specified bounding box. 

333 

334 The spherical box region and filtering function are passed to the 

335 generic loadRegion method which loads and filters the reference objects 

336 from the datastore and returns a single catalog containing the filtered 

337 set of reference objects. 

338 

339 Parameters 

340 ---------- 

341 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D` 

342 Box which bounds a region in pixel space. 

343 wcs : `lsst.afw.geom.SkyWcs` 

344 Wcs object defining the pixel to sky (and inverse) transform for 

345 the supplied ``bbox``. 

346 filterName : `str` or `None`, optional 

347 Name of camera filter, or `None` or blank for the default filter. 

348 epoch : `astropy.time.Time` or `None`, optional 

349 Epoch to which to correct proper motion and parallax, or `None` 

350 to not apply such corrections. 

351 photoCalib : `None` 

352 Deprecated and ignored, only included for api compatibility. 

353 bboxToSpherePadding : `int`, optional 

354 Padding to account for translating a set of corners into a 

355 spherical (convex) boundary that is certain to encompase the 

356 enitre area covered by the box. 

357 

358 Returns 

359 ------- 

360 referenceCatalog : `lsst.afw.table.SimpleCatalog` 

361 Catalog containing reference objects inside the specified bounding 

362 box (padded by self.config.pixelMargin). 

363 

364 Raises 

365 ------ 

366 RuntimeError 

367 Raised if no reference catalogs could be found for the specified 

368 region. 

369 TypeError 

370 Raised if the loaded reference catalogs do not have matching 

371 schemas. 

372 """ 

373 paddedBbox = geom.Box2D(bbox) 

374 paddedBbox.grow(self.config.pixelMargin) 

375 innerSkyRegion, outerSkyRegion, _, _ = self._makeBoxRegion(paddedBbox, wcs, bboxToSpherePadding) 

376 

377 def _filterFunction(refCat, region): 

378 # Add columns to the reference catalog relating their coordinates to 

379 # equivalent pixel positions for the wcs provided and use afwTable 

380 # to populate those columns. 

381 refCat = self.remapReferenceCatalogSchema(refCat, position=True) 

382 afwTable.updateRefCentroids(wcs, refCat) 

383 # No need to filter the catalog if it is entirely contained in the 

384 # region defined by the inner sky region. 

385 if innerSkyRegion.contains(region): 

386 return refCat 

387 # Create a new reference catalog, and populate it only with records 

388 # that fall inside the padded bbox. 

389 filteredRefCat = type(refCat)(refCat.table) 

390 centroidKey = afwTable.Point2DKey(refCat.schema['centroid']) 

391 for record in refCat: 

392 pixCoords = record[centroidKey] 

393 if paddedBbox.contains(geom.Point2D(pixCoords)): 

394 filteredRefCat.append(record) 

395 return filteredRefCat 

396 return self.loadRegion(outerSkyRegion, filtFunc=_filterFunction, epoch=epoch, filterName=filterName) 

397 

398 def loadRegion(self, region, filtFunc=None, filterName=None, epoch=None): 

399 """Load reference objects within a specified region. 

400 

401 This function loads the DataIds used to construct an instance of this 

402 class which intersect or are contained within the specified region. The 

403 reference catalogs which intersect but are not fully contained within 

404 the input region are further filtered by the specified filter function. 

405 This function returns a single source catalog containing all reference 

406 objects inside the specified region. 

407 

408 Parameters 

409 ---------- 

410 region : `lsst.sphgeom.Region` 

411 This can be any type that is derived from `lsst.sphgeom.Region` and 

412 should define the spatial region for which reference objects are to 

413 be loaded. 

414 filtFunc : callable or `None`, optional 

415 This optional parameter should be a callable object that takes a 

416 reference catalog and its corresponding region as parameters, 

417 filters the catalog by some criteria and returns the filtered 

418 reference catalog. If `None`, an internal filter function is used 

419 which filters according to if a reference object falls within the 

420 input region. 

421 filterName : `str` or `None`, optional 

422 Name of camera filter, or `None` or blank for the default filter. 

423 epoch : `astropy.time.Time` or `None`, optional 

424 Epoch to which to correct proper motion and parallax, or `None` to 

425 not apply such corrections. 

426 

427 Returns 

428 ------- 

429 referenceCatalog : `lsst.afw.table.SourceCatalog` 

430 Catalog containing reference objects which intersect the input region, 

431 filtered by the specified filter function. 

432 

433 Raises 

434 ------ 

435 RuntimeError 

436 Raised if no reference catalogs could be found for the specified 

437 region. 

438 TypeError 

439 Raised if the loaded reference catalogs do not have matching 

440 schemas. 

441 """ 

442 regionLat = region.getBoundingBox().getLat() 

443 regionLon = region.getBoundingBox().getLon() 

444 self.log.info("Loading reference objects from region bounded by " 

445 "[{:.8f}, {:.8f}], [{:.8f}, {:.8f}] RA Dec". 

446 format(regionLon.getA().asDegrees(), regionLon.getB().asDegrees(), 

447 regionLat.getA().asDegrees(), regionLat.getB().asDegrees())) 

448 if filtFunc is None: 

449 filtFunc = _FilterCatalog(region) 

450 # filter out all the regions supplied by the constructor that do not overlap 

451 overlapList = [] 

452 for dataId, refCat in zip(self.dataIds, self.refCats): 

453 # SphGeom supports some objects intersecting others, but is not symmetric, 

454 # try the intersect operation in both directions 

455 try: 

456 intersects = dataId.region.intersects(region) 

457 except TypeError: 

458 intersects = region.intersects(dataId.region) 

459 

460 if intersects: 

461 overlapList.append((dataId, refCat)) 

462 

463 if len(overlapList) == 0: 

464 raise RuntimeError("No reference tables could be found for input region") 

465 

466 firstCat = overlapList[0][1].get() 

467 refCat = filtFunc(firstCat, overlapList[0][0].region) 

468 trimmedAmount = len(firstCat) - len(refCat) 

469 

470 # Load in the remaining catalogs 

471 for dataId, inputRefCat in overlapList[1:]: 

472 tmpCat = inputRefCat.get() 

473 

474 if tmpCat.schema != firstCat.schema: 

475 raise TypeError("Reference catalogs have mismatching schemas") 

476 

477 filteredCat = filtFunc(tmpCat, dataId.region) 

478 refCat.extend(filteredCat) 

479 trimmedAmount += len(tmpCat) - len(filteredCat) 

480 

481 self.log.debug(f"Trimmed {trimmedAmount} refCat objects lying outside padded region, " 

482 "leaving {len(refCat)}") 

483 self.log.info(f"Loaded {len(refCat)} reference objects") 

484 

485 # Ensure that the loaded reference catalog is continuous in memory 

486 if not refCat.isContiguous(): 

487 refCat = refCat.copy(deep=True) 

488 

489 self.applyProperMotions(refCat, epoch) 

490 

491 # Verify the schema is in the correct units and has the correct version; automatically convert 

492 # it with a warning if this is not the case. 

493 if not hasNanojanskyFluxUnits(refCat.schema) or not getFormatVersionFromRefCat(refCat) >= 1: 

494 self.log.warn("Found version 0 reference catalog with old style units in schema.") 

495 self.log.warn("run `meas_algorithms/bin/convert_refcat_to_nJy.py` to convert fluxes to nJy.") 

496 self.log.warn("See RFC-575 for more details.") 

497 refCat = convertToNanojansky(refCat, self.log) 

498 

499 expandedCat = self.remapReferenceCatalogSchema(refCat, position=True) 

500 

501 # Add flux aliases 

502 expandedCat = self.addFluxAliases(expandedCat, self.config.defaultFilter, self.config.filterMap) 

503 

504 # Ensure that the returned reference catalog is continuous in memory 

505 if not expandedCat.isContiguous(): 

506 expandedCat = expandedCat.copy(deep=True) 

507 

508 fluxField = getRefFluxField(schema=expandedCat.schema, filterName=filterName) 

509 return pipeBase.Struct(refCat=expandedCat, fluxField=fluxField) 

510 

511 def loadSkyCircle(self, ctrCoord, radius, filterName=None, epoch=None): 

512 """Load reference objects that lie within a circular region on the sky. 

513 

514 This method constructs a circular region from an input center and 

515 angular radius, loads reference catalogs which are contained in or 

516 intersect the circle, and filters reference catalogs which intersect 

517 down to objects which lie within the defined circle. 

518 

519 Parameters 

520 ---------- 

521 ctrCoord : `lsst.geom.SpherePoint` 

522 Point defining the center of the circular region. 

523 radius : `lsst.geom.Angle` 

524 Defines the angular radius of the circular region. 

525 filterName : `str` or `None`, optional 

526 Name of camera filter, or `None` or blank for the default filter. 

527 epoch : `astropy.time.Time` or `None`, optional 

528 Epoch to which to correct proper motion and parallax, or `None` to 

529 not apply such corrections. 

530 

531 Returns 

532 ------- 

533 referenceCatalog : `lsst.afw.table.SourceCatalog` 

534 Catalog containing reference objects inside the specified search 

535 circle. 

536 """ 

537 centerVector = ctrCoord.getVector() 

538 sphRadius = sphgeom.Angle(radius.asRadians()) 

539 circularRegion = sphgeom.Circle(centerVector, sphRadius) 

540 return self.loadRegion(circularRegion, filterName=filterName, epoch=epoch) 

541 

542 def joinMatchListWithCatalog(self, matchCat, sourceCat): 

543 """Relink an unpersisted match list to sources and reference objects. 

544 

545 A match list is persisted and unpersisted as a catalog of IDs 

546 produced by afw.table.packMatches(), with match metadata 

547 (as returned by the astrometry tasks) in the catalog's metadata 

548 attribute. This method converts such a match catalog into a match 

549 list, with links to source records and reference object records. 

550 

551 Parameters 

552 ---------- 

553 matchCat : `lsst.afw.table.BaseCatalog` 

554 Unpersisted packed match list. 

555 ``matchCat.table.getMetadata()`` must contain match metadata, 

556 as returned by the astrometry tasks. 

557 sourceCat : `lsst.afw.table.SourceCatalog` 

558 Source catalog. As a side effect, the catalog will be sorted 

559 by ID. 

560 

561 Returns 

562 ------- 

563 matchList : `lsst.afw.table.ReferenceMatchVector` 

564 Match list. 

565 """ 

566 return joinMatchListWithCatalogImpl(self, matchCat, sourceCat) 

567 

568 def getMetadataBox(self, bbox, wcs, filterName=None, photoCalib=None, epoch=None, 

569 bboxToSpherePadding=100): 

570 """Return metadata about the load 

571 

572 This metadata is used for reloading the catalog (e.g., for 

573 reconstituting a normalised match list.) 

574 

575 Parameters 

576 ---------- 

577 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D` 

578 Bounding box for the pixels. 

579 wcs : `lsst.afw.geom.SkyWcs` 

580 The WCS object associated with ``bbox``. 

581 filterName : `str` or `None`, optional 

582 Name of the camera filter, or `None` or blank for the default 

583 filter. 

584 photoCalib : `None` 

585 Deprecated, only included for api compatibility. 

586 epoch : `astropy.time.Time` or `None`, optional 

587 Epoch to which to correct proper motion and parallax, or `None` to 

588 not apply such corrections. 

589 bboxToSpherePadding : `int`, optional 

590 Padding to account for translating a set of corners into a 

591 spherical (convex) boundary that is certain to encompase the 

592 enitre area covered by the box. 

593 

594 Returns 

595 ------- 

596 md : `lsst.daf.base.PropertyList` 

597 The metadata detailing the search parameters used for this 

598 dataset. 

599 """ 

600 paddedBbox = geom.Box2D(bbox) 

601 paddedBbox.grow(self.config.pixelMargin) 

602 _, _, innerCorners, outerCorners = self._makeBoxRegion(paddedBbox, wcs, bboxToSpherePadding) 

603 md = PropertyList() 

604 for box, corners in zip(("INNER", "OUTER"), (innerCorners, outerCorners)): 

605 for (name, corner) in zip(("UPPER_LEFT", "UPPER_RIGHT", "LOWER_LEFT", "LOWER_RIGHT"), 

606 corners): 

607 md.add(f"{box}_{name}_RA", geom.SpherePoint(corner).getRa().asDegrees(), f"{box}_corner") 

608 md.add(f"{box}_{name}_DEC", geom.SpherePoint(corner).getDec().asDegrees(), f"{box}_corner") 

609 md.add("SMATCHV", 1, 'SourceMatchVector version number') 

610 filterName = "UNKNOWN" if filterName is None else str(filterName) 

611 md.add('FILTER', filterName, 'filter name for photometric data') 

612 md.add('EPOCH', "NONE" if epoch is None else epoch.mjd, 'Epoch (TAI MJD) for catalog') 

613 return md 

614 

615 @staticmethod 

616 def getMetadataCircle(coord, radius, filterName, photoCalib=None, epoch=None): 

617 """Return metadata about the load. 

618 

619 This metadata is used for reloading the catalog (e.g. for 

620 reconstituting a normalized match list.) 

621 

622 Parameters 

623 ---------- 

624 coord : `lsst.geom.SpherePoint` 

625 ICRS center of the search region. 

626 radius : `lsst.geom.Angle` 

627 Radius of the search region. 

628 filterName : `str` or `None` 

629 Name of the camera filter, or `None` or blank for the default 

630 filter. 

631 photoCalib : `None` 

632 Deprecated, only included for api compatibility. 

633 epoch : `astropy.time.Time` or `None`, optional 

634 Epoch to which to correct proper motion and parallax, or `None` to 

635 not apply such corrections. 

636 

637 Returns 

638 ------- 

639 md : `lsst.daf.base.PropertyList` 

640 """ 

641 md = PropertyList() 

642 md.add('RA', coord.getRa().asDegrees(), 'field center in degrees') 

643 md.add('DEC', coord.getDec().asDegrees(), 'field center in degrees') 

644 md.add('RADIUS', radius.asDegrees(), 'field radius in degrees, minimum') 

645 md.add('SMATCHV', 1, 'SourceMatchVector version number') 

646 filterName = "UNKNOWN" if filterName is None else str(filterName) 

647 md.add('FILTER', filterName, 'filter name for photometric data') 

648 md.add('EPOCH', "NONE" if epoch is None else epoch.mjd, 'Epoch (TAI MJD) for catalog') 

649 return md 

650 

651 @staticmethod 

652 def addFluxAliases(refCat, defaultFilter, filterReferenceMap): 

653 """Add flux columns and aliases for camera to reference mapping. 

654 

655 Creates a new catalog containing the information of the input refCat 

656 as well as added flux columns and aliases between camera and reference 

657 fluxes. 

658 

659 Parameters 

660 ---------- 

661 refCat : `lsst.afw.table.SimpleCatalog` 

662 Catalog of reference objects 

663 defaultFilter : `str` 

664 Name of the default reference filter 

665 filterReferenceMap : `dict` of `str` 

666 Dictionary with keys corresponding to a filter name and values 

667 which correspond to the name of the reference filter. 

668 

669 Returns 

670 ------- 

671 refCat : `lsst.afw.table.SimpleCatalog` 

672 Reference catalog with columns added to track reference filters. 

673 

674 Raises 

675 ------ 

676 `RuntimeError` 

677 If the specified reference filter name is not specifed as a 

678 key in the reference filter map. 

679 """ 

680 refCat = ReferenceObjectLoader.remapReferenceCatalogSchema(refCat, 

681 filterNameList=filterReferenceMap.keys()) 

682 aliasMap = refCat.schema.getAliasMap() 

683 if filterReferenceMap is None: 

684 filterReferenceMap = {} 

685 for filterName, refFilterName in itertools.chain([(None, defaultFilter)], 

686 filterReferenceMap.items()): 

687 if refFilterName: 

688 camFluxName = filterName + "_camFlux" if filterName is not None else "camFlux" 

689 refFluxName = refFilterName + "_flux" 

690 if refFluxName not in refCat.schema: 

691 raise RuntimeError("Unknown reference filter %s" % (refFluxName,)) 

692 aliasMap.set(camFluxName, refFluxName) 

693 

694 refFluxErrName = refFluxName + "Err" 

695 camFluxErrName = camFluxName + "Err" 

696 aliasMap.set(camFluxErrName, refFluxErrName) 

697 

698 return refCat 

699 

700 @staticmethod 

701 def remapReferenceCatalogSchema(refCat, *, filterNameList=None, position=False, photometric=False): 

702 """This function takes in a reference catalog and creates a new catalog with additional 

703 columns defined the remaining function arguments. 

704 

705 Parameters 

706 ---------- 

707 refCat : `lsst.afw.table.SimpleCatalog` 

708 Reference catalog to map to new catalog 

709 

710 Returns 

711 ------- 

712 expandedCat : `lsst.afw.table.SimpleCatalog` 

713 Deep copy of input reference catalog with additional columns added 

714 """ 

715 mapper = afwTable.SchemaMapper(refCat.schema, True) 

716 mapper.addMinimalSchema(refCat.schema, True) 

717 mapper.editOutputSchema().disconnectAliases() 

718 if filterNameList: 

719 for filterName in filterNameList: 

720 mapper.editOutputSchema().addField(f"{filterName}_flux", 

721 type=numpy.float64, 

722 doc=f"flux in filter {filterName}", 

723 units="Jy" 

724 ) 

725 mapper.editOutputSchema().addField(f"{filterName}_fluxErr", 

726 type=numpy.float64, 

727 doc=f"flux uncertanty in filter {filterName}", 

728 units="Jy" 

729 ) 

730 

731 if position: 

732 mapper.editOutputSchema().addField("centroid_x", type=float, doReplace=True) 

733 mapper.editOutputSchema().addField("centroid_y", type=float, doReplace=True) 

734 mapper.editOutputSchema().addField("hasCentroid", type="Flag", doReplace=True) 

735 mapper.editOutputSchema().getAliasMap().set("slot_Centroid", "centroid") 

736 

737 if photometric: 

738 mapper.editOutputSchema().addField("photometric", 

739 type="Flag", 

740 doc="set if the object can be used for photometric" 

741 "calibration", 

742 ) 

743 mapper.editOutputSchema().addField("resolved", 

744 type="Flag", 

745 doc="set if the object is spatially resolved" 

746 ) 

747 mapper.editOutputSchema().addField("variable", 

748 type="Flag", 

749 doc="set if the object has variable brightness" 

750 ) 

751 

752 expandedCat = afwTable.SimpleCatalog(mapper.getOutputSchema()) 

753 expandedCat.setMetadata(refCat.getMetadata()) 

754 expandedCat.extend(refCat, mapper=mapper) 

755 

756 return expandedCat 

757 

758 

759def getRefFluxField(schema, filterName=None): 

760 """Get the name of a flux field from a schema. 

761 

762 return the alias of "anyFilterMapsToThis", if present 

763 else if filterName is specified: 

764 return "*filterName*_camFlux" if present 

765 else return "*filterName*_flux" if present (camera filter name 

766 matches reference filter name) 

767 else throw RuntimeError 

768 else: 

769 return "camFlux", if present, 

770 else throw RuntimeError 

771 

772 Parameters 

773 ---------- 

774 schema : `lsst.afw.table.Schema` 

775 Reference catalog schema. 

776 filterName : `str`, optional 

777 Name of camera filter. If not specified, ``defaultFilter`` needs to be 

778 set in the refcat loader config. 

779 

780 Returns 

781 ------- 

782 fluxFieldName : `str` 

783 Name of flux field. 

784 

785 Raises 

786 ------ 

787 RuntimeError 

788 If an appropriate field is not found. 

789 """ 

790 if not isinstance(schema, afwTable.Schema): 790 ↛ 791line 790 didn't jump to line 791, because the condition on line 790 was never true

791 raise RuntimeError("schema=%s is not a schema" % (schema,)) 

792 try: 

793 return schema.getAliasMap().get("anyFilterMapsToThis") 

794 except LookupError: 

795 pass # try the filterMap next 

796 

797 if filterName: 797 ↛ 800line 797 didn't jump to line 800, because the condition on line 797 was never false

798 fluxFieldList = [filterName + "_camFlux", filterName + "_flux"] 

799 else: 

800 fluxFieldList = ["camFlux"] 

801 for fluxField in fluxFieldList: 801 ↛ 805line 801 didn't jump to line 805, because the loop on line 801 didn't complete

802 if fluxField in schema: 

803 return fluxField 

804 

805 raise RuntimeError("Could not find flux field(s) %s" % (", ".join(fluxFieldList))) 

806 

807 

808def getRefFluxKeys(schema, filterName=None): 

809 """Return keys for flux and flux error. 

810 

811 Parameters 

812 ---------- 

813 schema : `lsst.afw.table.Schema` 

814 Reference catalog schema. 

815 filterName : `str` 

816 Name of camera filter. 

817 

818 Returns 

819 ------- 

820 keys : `tuple` of (`lsst.afw.table.Key`, `lsst.afw.table.Key`) 

821 Two keys: 

822 

823 - flux key 

824 - flux error key, if present, else None 

825 

826 Raises 

827 ------ 

828 RuntimeError 

829 If flux field not found. 

830 """ 

831 fluxField = getRefFluxField(schema, filterName) 

832 fluxErrField = fluxField + "Err" 

833 fluxKey = schema[fluxField].asKey() 

834 try: 

835 fluxErrKey = schema[fluxErrField].asKey() 

836 except Exception: 

837 fluxErrKey = None 

838 return (fluxKey, fluxErrKey) 

839 

840 

841class LoadReferenceObjectsConfig(pexConfig.Config): 

842 pixelMargin = pexConfig.RangeField( 

843 doc="Padding to add to 4 all edges of the bounding box (pixels)", 

844 dtype=int, 

845 default=300, 

846 min=0, 

847 ) 

848 defaultFilter = pexConfig.Field( 

849 doc=("Default reference catalog filter to use if filter not specified in exposure;" 

850 " if blank then filter must be specified in exposure."), 

851 dtype=str, 

852 default="", 

853 deprecated="defaultFilter is deprecated by RFC-716. Will be removed after v22." 

854 ) 

855 anyFilterMapsToThis = pexConfig.Field( 

856 doc=("Always use this reference catalog filter, no matter whether or what filter name is " 

857 "supplied to the loader. Effectively a trivial filterMap: map all filter names to this filter." 

858 " This can be set for purely-astrometric catalogs (e.g. Gaia DR2) where there is only one " 

859 "reasonable choice for every camera filter->refcat mapping, but not for refcats used for " 

860 "photometry, which need a filterMap and/or colorterms/transmission corrections."), 

861 dtype=str, 

862 default=None, 

863 optional=True 

864 ) 

865 filterMap = pexConfig.DictField( 

866 doc=("Mapping of camera filter name: reference catalog filter name; " 

867 "each reference filter must exist in the refcat." 

868 " Note that this does not perform any bandpass corrections: it is just a lookup."), 

869 keytype=str, 

870 itemtype=str, 

871 default={}, 

872 ) 

873 requireProperMotion = pexConfig.Field( 

874 doc="Require that the fields needed to correct proper motion " 

875 "(epoch, pm_ra and pm_dec) are present?", 

876 dtype=bool, 

877 default=False, 

878 ) 

879 

880 def validate(self): 

881 super().validate() 

882 if self.filterMap != {} and self.anyFilterMapsToThis is not None: 

883 msg = "`filterMap` and `anyFilterMapsToThis` are mutually exclusive" 

884 raise pexConfig.FieldValidationError(LoadReferenceObjectsConfig.anyFilterMapsToThis, 

885 self, msg) 

886 

887 

888class LoadReferenceObjectsTask(pipeBase.Task, ReferenceObjectLoaderBase, metaclass=abc.ABCMeta): 

889 """Abstract base class to load objects from reference catalogs. 

890 """ 

891 ConfigClass = LoadReferenceObjectsConfig 

892 _DefaultName = "LoadReferenceObjects" 

893 

894 def __init__(self, butler=None, *args, **kwargs): 

895 """Construct a LoadReferenceObjectsTask 

896 

897 Parameters 

898 ---------- 

899 butler : `lsst.daf.persistence.Butler` 

900 Data butler, for access reference catalogs. 

901 """ 

902 pipeBase.Task.__init__(self, *args, **kwargs) 

903 self.butler = butler 

904 

905 @pipeBase.timeMethod 

906 def loadPixelBox(self, bbox, wcs, filterName=None, photoCalib=None, epoch=None): 

907 """Load reference objects that overlap a rectangular pixel region. 

908 

909 Parameters 

910 ---------- 

911 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D` 

912 Bounding box for pixels. 

913 wcs : `lsst.afw.geom.SkyWcs` 

914 WCS; used to convert pixel positions to sky coordinates 

915 and vice-versa. 

916 filterName : `str` or `None`, optional 

917 Name of filter, or `None` or `""` for the default filter. 

918 This is used for flux values in case we have flux limits 

919 (which are not yet implemented). 

920 photoCalib : `None` 

921 Deprecated, only included for api compatibility. 

922 epoch : `astropy.time.Time` or `None`, optional 

923 Epoch to which to correct proper motion and parallax, or `None` to 

924 not apply such corrections. 

925 

926 Returns 

927 ------- 

928 results : `lsst.pipe.base.Struct` 

929 A Struct containing the following fields: 

930 refCat : `lsst.afw.catalog.SimpleCatalog` 

931 A catalog of reference objects with the standard 

932 schema, as documented in the main doc string for 

933 `LoadReferenceObjects`. 

934 The catalog is guaranteed to be contiguous. 

935 fluxField : `str` 

936 Name of flux field for specified `filterName`. 

937 

938 Notes 

939 ----- 

940 The search algorithm works by searching in a region in sky 

941 coordinates whose center is the center of the bbox and radius 

942 is large enough to just include all 4 corners of the bbox. 

943 Stars that lie outside the bbox are then trimmed from the list. 

944 """ 

945 circle = self._calculateCircle(bbox, wcs) 

946 

947 # find objects in circle 

948 self.log.info("Loading reference objects using center %s and radius %s deg" % 

949 (circle.coord, circle.radius.asDegrees())) 

950 loadRes = self.loadSkyCircle(circle.coord, circle.radius, filterName=filterName, epoch=epoch, 

951 centroids=True) 

952 refCat = loadRes.refCat 

953 numFound = len(refCat) 

954 

955 # trim objects outside bbox 

956 refCat = self._trimToBBox(refCat=refCat, bbox=circle.bbox, wcs=wcs) 

957 numTrimmed = numFound - len(refCat) 

958 self.log.debug("trimmed %d out-of-bbox objects, leaving %d", numTrimmed, len(refCat)) 

959 self.log.info("Loaded %d reference objects", len(refCat)) 

960 

961 # make sure catalog is contiguous 

962 if not refCat.isContiguous(): 

963 loadRes.refCat = refCat.copy(deep=True) 

964 

965 return loadRes 

966 

967 @abc.abstractmethod 

968 def loadSkyCircle(self, ctrCoord, radius, filterName=None, epoch=None, centroids=False): 

969 """Load reference objects that overlap a circular sky region. 

970 

971 Parameters 

972 ---------- 

973 ctrCoord : `lsst.geom.SpherePoint` 

974 ICRS center of search region. 

975 radius : `lsst.geom.Angle` 

976 Radius of search region. 

977 filterName : `str` or `None`, optional 

978 Name of filter, or `None` or `""` for the default filter. 

979 This is used for flux values in case we have flux limits 

980 (which are not yet implemented). 

981 epoch : `astropy.time.Time` or `None`, optional 

982 Epoch to which to correct proper motion and parallax, or `None` to 

983 not apply such corrections. 

984 centroids : `bool`, optional 

985 Add centroid fields to the loaded Schema. ``loadPixelBox`` expects 

986 these fields to exist. 

987 

988 Returns 

989 ------- 

990 results : `lsst.pipe.base.Struct` 

991 A Struct containing the following fields: 

992 refCat : `lsst.afw.catalog.SimpleCatalog` 

993 A catalog of reference objects with the standard 

994 schema, as documented in the main doc string for 

995 `LoadReferenceObjects`. 

996 The catalog is guaranteed to be contiguous. 

997 fluxField : `str` 

998 Name of flux field for specified `filterName`. 

999 

1000 Notes 

1001 ----- 

1002 Note that subclasses are responsible for performing the proper motion 

1003 correction, since this is the lowest-level interface for retrieving 

1004 the catalog. 

1005 """ 

1006 return 

1007 

1008 @staticmethod 

1009 def _trimToBBox(refCat, bbox, wcs): 

1010 """Remove objects outside a given pixel bounding box and set 

1011 centroid and hasCentroid fields. 

1012 

1013 Parameters 

1014 ---------- 

1015 refCat : `lsst.afw.table.SimpleCatalog` 

1016 A catalog of objects. The schema must include fields 

1017 "coord", "centroid" and "hasCentroid". 

1018 The "coord" field is read. 

1019 The "centroid" and "hasCentroid" fields are set. 

1020 bbox : `lsst.geom.Box2D` 

1021 Pixel region 

1022 wcs : `lsst.afw.geom.SkyWcs` 

1023 WCS; used to convert sky coordinates to pixel positions. 

1024 

1025 Returns 

1026 ------- 

1027 catalog : `lsst.afw.table.SimpleCatalog` 

1028 Reference objects in the bbox, with centroid and 

1029 hasCentroid fields set. 

1030 """ 

1031 afwTable.updateRefCentroids(wcs, refCat) 

1032 centroidKey = afwTable.Point2DKey(refCat.schema["centroid"]) 

1033 retStarCat = type(refCat)(refCat.table) 

1034 for star in refCat: 

1035 point = star.get(centroidKey) 

1036 if bbox.contains(point): 

1037 retStarCat.append(star) 

1038 return retStarCat 

1039 

1040 def _addFluxAliases(self, schema): 

1041 """Add aliases for camera filter fluxes to the schema. 

1042 

1043 If self.config.defaultFilter then adds these aliases: 

1044 camFlux: <defaultFilter>_flux 

1045 camFluxErr: <defaultFilter>_fluxErr, if the latter exists 

1046 

1047 For each camFilter: refFilter in self.config.filterMap adds these aliases: 

1048 <camFilter>_camFlux: <refFilter>_flux 

1049 <camFilter>_camFluxErr: <refFilter>_fluxErr, if the latter exists 

1050 

1051 Parameters 

1052 ---------- 

1053 schema : `lsst.afw.table.Schema` 

1054 Schema for reference catalog. 

1055 

1056 Raises 

1057 ------ 

1058 RuntimeError 

1059 If any reference flux field is missing from the schema. 

1060 """ 

1061 aliasMap = schema.getAliasMap() 

1062 

1063 if self.config.anyFilterMapsToThis is not None: 1063 ↛ 1064line 1063 didn't jump to line 1064, because the condition on line 1063 was never true

1064 refFluxName = self.config.anyFilterMapsToThis + "_flux" 

1065 if refFluxName not in schema: 

1066 msg = f"Unknown reference filter for anyFilterMapsToThis='{refFluxName}'" 

1067 raise RuntimeError(msg) 

1068 aliasMap.set("anyFilterMapsToThis", refFluxName) 

1069 return # this is mutually exclusive with filterMap 

1070 

1071 def addAliasesForOneFilter(filterName, refFilterName): 

1072 """Add aliases for a single filter 

1073 

1074 Parameters 

1075 ---------- 

1076 filterName : `str` (optional) 

1077 Camera filter name. The resulting alias name is 

1078 <filterName>_camFlux, or simply "camFlux" if `filterName` 

1079 is `None` or `""`. 

1080 refFilterName : `str` 

1081 Reference catalog filter name; the field 

1082 <refFilterName>_flux must exist. 

1083 """ 

1084 camFluxName = filterName + "_camFlux" if filterName is not None else "camFlux" 

1085 refFluxName = refFilterName + "_flux" 

1086 if refFluxName not in schema: 

1087 raise RuntimeError("Unknown reference filter %s" % (refFluxName,)) 

1088 aliasMap.set(camFluxName, refFluxName) 

1089 refFluxErrName = refFluxName + "Err" 

1090 if refFluxErrName in schema: 

1091 camFluxErrName = camFluxName + "Err" 

1092 aliasMap.set(camFluxErrName, refFluxErrName) 

1093 

1094 if self.config.defaultFilter: 1094 ↛ 1095line 1094 didn't jump to line 1095, because the condition on line 1094 was never true

1095 addAliasesForOneFilter(None, self.config.defaultFilter) 

1096 

1097 for filterName, refFilterName in self.config.filterMap.items(): 1097 ↛ 1098line 1097 didn't jump to line 1098, because the loop on line 1097 never started

1098 addAliasesForOneFilter(filterName, refFilterName) 

1099 

1100 @staticmethod 

1101 def makeMinimalSchema(filterNameList, *, addCentroid=False, 

1102 addIsPhotometric=False, addIsResolved=False, 

1103 addIsVariable=False, coordErrDim=2, 

1104 addProperMotion=False, properMotionErrDim=2, 

1105 addParallax=False): 

1106 """Make a standard schema for reference object catalogs. 

1107 

1108 Parameters 

1109 ---------- 

1110 filterNameList : `list` of `str` 

1111 List of filter names. Used to create <filterName>_flux fields. 

1112 addIsPhotometric : `bool` 

1113 If True then add field "photometric". 

1114 addIsResolved : `bool` 

1115 If True then add field "resolved". 

1116 addIsVariable : `bool` 

1117 If True then add field "variable". 

1118 coordErrDim : `int` 

1119 Number of coord error fields; must be one of 0, 2, 3: 

1120 

1121 - If 2 or 3: add fields "coord_raErr" and "coord_decErr". 

1122 - If 3: also add field "coord_radecErr". 

1123 addProperMotion : `bool` 

1124 If True add fields "epoch", "pm_ra", "pm_dec" and "pm_flag". 

1125 properMotionErrDim : `int` 

1126 Number of proper motion error fields; must be one of 0, 2, 3; 

1127 ignored if addProperMotion false: 

1128 - If 2 or 3: add fields "pm_raErr" and "pm_decErr". 

1129 - If 3: also add field "pm_radecErr". 

1130 addParallax : `bool` 

1131 If True add fields "epoch", "parallax", "parallaxErr" 

1132 and "parallax_flag". 

1133 

1134 Returns 

1135 ------- 

1136 schema : `lsst.afw.table.Schema` 

1137 Schema for reference catalog, an 

1138 `lsst.afw.table.SimpleCatalog`. 

1139 

1140 Notes 

1141 ----- 

1142 Reference catalogs support additional covariances, such as 

1143 covariance between RA and proper motion in declination, 

1144 that are not supported by this method, but can be added after 

1145 calling this method. 

1146 """ 

1147 schema = afwTable.SimpleTable.makeMinimalSchema() 

1148 if addCentroid: 1148 ↛ 1149line 1148 didn't jump to line 1149, because the condition on line 1148 was never true

1149 afwTable.Point2DKey.addFields( 

1150 schema, 

1151 "centroid", 

1152 "centroid on an exposure, if relevant", 

1153 "pixel", 

1154 ) 

1155 schema.addField( 

1156 field="hasCentroid", 

1157 type="Flag", 

1158 doc="is position known?", 

1159 ) 

1160 for filterName in filterNameList: 

1161 schema.addField( 

1162 field="%s_flux" % (filterName,), 

1163 type=numpy.float64, 

1164 doc="flux in filter %s" % (filterName,), 

1165 units="nJy", 

1166 ) 

1167 for filterName in filterNameList: 

1168 schema.addField( 

1169 field="%s_fluxErr" % (filterName,), 

1170 type=numpy.float64, 

1171 doc="flux uncertainty in filter %s" % (filterName,), 

1172 units="nJy", 

1173 ) 

1174 if addIsPhotometric: 1174 ↛ 1175line 1174 didn't jump to line 1175, because the condition on line 1174 was never true

1175 schema.addField( 

1176 field="photometric", 

1177 type="Flag", 

1178 doc="set if the object can be used for photometric calibration", 

1179 ) 

1180 if addIsResolved: 1180 ↛ 1181line 1180 didn't jump to line 1181, because the condition on line 1180 was never true

1181 schema.addField( 

1182 field="resolved", 

1183 type="Flag", 

1184 doc="set if the object is spatially resolved", 

1185 ) 

1186 if addIsVariable: 1186 ↛ 1187line 1186 didn't jump to line 1187, because the condition on line 1186 was never true

1187 schema.addField( 

1188 field="variable", 

1189 type="Flag", 

1190 doc="set if the object has variable brightness", 

1191 ) 

1192 if coordErrDim not in (0, 2, 3): 1192 ↛ 1193line 1192 didn't jump to line 1193, because the condition on line 1192 was never true

1193 raise ValueError("coordErrDim={}; must be (0, 2, 3)".format(coordErrDim)) 

1194 if coordErrDim > 0: 

1195 afwTable.CovarianceMatrix2fKey.addFields( 

1196 schema=schema, 

1197 prefix="coord", 

1198 names=["ra", "dec"], 

1199 units=["rad", "rad"], 

1200 diagonalOnly=(coordErrDim == 2), 

1201 ) 

1202 

1203 if addProperMotion or addParallax: 

1204 schema.addField( 

1205 field="epoch", 

1206 type=numpy.float64, 

1207 doc="date of observation (TAI, MJD)", 

1208 units="day", 

1209 ) 

1210 

1211 if addProperMotion: 

1212 schema.addField( 

1213 field="pm_ra", 

1214 type="Angle", 

1215 doc="proper motion in the right ascension direction = dra/dt * cos(dec)", 

1216 units="rad/year", 

1217 ) 

1218 schema.addField( 

1219 field="pm_dec", 

1220 type="Angle", 

1221 doc="proper motion in the declination direction", 

1222 units="rad/year", 

1223 ) 

1224 if properMotionErrDim not in (0, 2, 3): 1224 ↛ 1225line 1224 didn't jump to line 1225, because the condition on line 1224 was never true

1225 raise ValueError("properMotionErrDim={}; must be (0, 2, 3)".format(properMotionErrDim)) 

1226 if properMotionErrDim > 0: 1226 ↛ 1234line 1226 didn't jump to line 1234, because the condition on line 1226 was never false

1227 afwTable.CovarianceMatrix2fKey.addFields( 

1228 schema=schema, 

1229 prefix="pm", 

1230 names=["ra", "dec"], 

1231 units=["rad/year", "rad/year"], 

1232 diagonalOnly=(properMotionErrDim == 2), 

1233 ) 

1234 schema.addField( 

1235 field="pm_flag", 

1236 type="Flag", 

1237 doc="Set if proper motion or proper motion error is bad", 

1238 ) 

1239 

1240 if addParallax: 

1241 schema.addField( 

1242 field="parallax", 

1243 type="Angle", 

1244 doc="parallax", 

1245 units="rad", 

1246 ) 

1247 schema.addField( 

1248 field="parallaxErr", 

1249 type="Angle", 

1250 doc="uncertainty in parallax", 

1251 units="rad", 

1252 ) 

1253 schema.addField( 

1254 field="parallax_flag", 

1255 type="Flag", 

1256 doc="Set if parallax or parallax error is bad", 

1257 ) 

1258 return schema 

1259 

1260 def _calculateCircle(self, bbox, wcs): 

1261 """Compute on-sky center and radius of search region. 

1262 

1263 Parameters 

1264 ---------- 

1265 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D` 

1266 Pixel bounding box. 

1267 wcs : `lsst.afw.geom.SkyWcs` 

1268 WCS; used to convert pixel positions to sky coordinates. 

1269 

1270 Returns 

1271 ------- 

1272 results : `lsst.pipe.base.Struct` 

1273 A Struct containing: 

1274 

1275 - coord : `lsst.geom.SpherePoint` 

1276 ICRS center of the search region. 

1277 - radius : `lsst.geom.Angle` 

1278 Radius of the search region. 

1279 - bbox : `lsst.geom.Box2D` 

1280 Bounding box used to compute the circle. 

1281 """ 

1282 bbox = geom.Box2D(bbox) # make sure bbox is double and that we have a copy 

1283 bbox.grow(self.config.pixelMargin) 

1284 coord = wcs.pixelToSky(bbox.getCenter()) 

1285 radius = max(coord.separation(wcs.pixelToSky(pp)) for pp in bbox.getCorners()) 

1286 return pipeBase.Struct(coord=coord, radius=radius, bbox=bbox) 

1287 

1288 def getMetadataBox(self, bbox, wcs, filterName=None, photoCalib=None, epoch=None): 

1289 """Return metadata about the load. 

1290 

1291 This metadata is used for reloading the catalog (e.g., for 

1292 reconstituting a normalised match list. 

1293 

1294 Parameters 

1295 ---------- 

1296 bbox : `lsst.geom.Box2I` or `lsst.geom.Box2D` 

1297 Pixel bounding box. 

1298 wcs : `lsst.afw.geom.SkyWcs` 

1299 WCS; used to convert pixel positions to sky coordinates. 

1300 filterName : `str` or `None`, optional 

1301 Name of camera filter, or `None` or `""` for the default 

1302 filter. 

1303 photoCalib : `None` 

1304 Deprecated, only included for api compatibility. 

1305 epoch : `astropy.time.Time` or `None`, optional 

1306 Epoch to which to correct proper motion and parallax, 

1307 or None to not apply such corrections. 

1308 

1309 Returns 

1310 ------- 

1311 metadata : `lsst.daf.base.PropertyList` 

1312 Metadata about the load. 

1313 """ 

1314 circle = self._calculateCircle(bbox, wcs) 

1315 return self.getMetadataCircle(circle.coord, circle.radius, filterName, epoch=epoch) 

1316 

1317 def getMetadataCircle(self, coord, radius, filterName, photoCalib=None, epoch=None): 

1318 """Return metadata about the load. 

1319 

1320 This metadata is used for reloading the catalog (e.g., for 

1321 reconstituting a normalised match list. 

1322 

1323 Parameters 

1324 ---------- 

1325 coord : `lsst.geom.SpherePoint` 

1326 ICRS center of the search region. 

1327 radius : `lsst.geom.Angle` 

1328 Radius of the search region. 

1329 filterName : `str` 

1330 Name of camera filter, or `None` or `""` for the default 

1331 filter. 

1332 photoCalib : `None` 

1333 Deprecated, only included for api compatibility. 

1334 epoch : `astropy.time.Time` (optional) 

1335 Epoch to which to correct proper motion and parallax, or `None` to 

1336 not apply such corrections. 

1337 

1338 Returns 

1339 ------- 

1340 metadata : lsst.daf.base.PropertyList 

1341 Metadata about the load 

1342 """ 

1343 md = PropertyList() 

1344 md.add('RA', coord.getRa().asDegrees(), 'field center in degrees') 

1345 md.add('DEC', coord.getDec().asDegrees(), 'field center in degrees') 

1346 md.add('RADIUS', radius.asDegrees(), 'field radius in degrees, minimum') 

1347 md.add('SMATCHV', 1, 'SourceMatchVector version number') 

1348 filterName = "UNKNOWN" if filterName is None else str(filterName) 

1349 md.add('FILTER', filterName, 'filter name for photometric data') 

1350 md.add('EPOCH', "NONE" if epoch is None else epoch.mjd, 'Epoch (TAI MJD) for catalog') 

1351 return md 

1352 

1353 def joinMatchListWithCatalog(self, matchCat, sourceCat): 

1354 """Relink an unpersisted match list to sources and reference 

1355 objects. 

1356 

1357 A match list is persisted and unpersisted as a catalog of IDs 

1358 produced by afw.table.packMatches(), with match metadata 

1359 (as returned by the astrometry tasks) in the catalog's metadata 

1360 attribute. This method converts such a match catalog into a match 

1361 list, with links to source records and reference object records. 

1362 

1363 Parameters 

1364 ---------- 

1365 matchCat : `lsst.afw.table.BaseCatalog` 

1366 Unperisted packed match list. 

1367 ``matchCat.table.getMetadata()`` must contain match metadata, 

1368 as returned by the astrometry tasks. 

1369 sourceCat : `lsst.afw.table.SourceCatalog` 

1370 Source catalog. As a side effect, the catalog will be sorted 

1371 by ID. 

1372 

1373 Returns 

1374 ------- 

1375 matchList : `lsst.afw.table.ReferenceMatchVector` 

1376 Match list. 

1377 """ 

1378 return joinMatchListWithCatalogImpl(self, matchCat, sourceCat) 

1379 

1380 

1381def joinMatchListWithCatalogImpl(refObjLoader, matchCat, sourceCat): 

1382 """Relink an unpersisted match list to sources and reference 

1383 objects. 

1384 

1385 A match list is persisted and unpersisted as a catalog of IDs 

1386 produced by afw.table.packMatches(), with match metadata 

1387 (as returned by the astrometry tasks) in the catalog's metadata 

1388 attribute. This method converts such a match catalog into a match 

1389 list, with links to source records and reference object records. 

1390 

1391 Parameters 

1392 ---------- 

1393 refObjLoader 

1394 Reference object loader to use in getting reference objects 

1395 matchCat : `lsst.afw.table.BaseCatalog` 

1396 Unperisted packed match list. 

1397 ``matchCat.table.getMetadata()`` must contain match metadata, 

1398 as returned by the astrometry tasks. 

1399 sourceCat : `lsst.afw.table.SourceCatalog` 

1400 Source catalog. As a side effect, the catalog will be sorted 

1401 by ID. 

1402 

1403 Returns 

1404 ------- 

1405 matchList : `lsst.afw.table.ReferenceMatchVector` 

1406 Match list. 

1407 """ 

1408 matchmeta = matchCat.table.getMetadata() 

1409 version = matchmeta.getInt('SMATCHV') 

1410 if version != 1: 

1411 raise ValueError('SourceMatchVector version number is %i, not 1.' % version) 

1412 filterName = matchmeta.getString('FILTER').strip() 

1413 try: 

1414 epoch = matchmeta.getDouble('EPOCH') 

1415 except (LookupError, TypeError): 

1416 epoch = None # Not present, or not correct type means it's not set 

1417 if 'RADIUS' in matchmeta: 

1418 # This is a circle style metadata, call loadSkyCircle 

1419 ctrCoord = geom.SpherePoint(matchmeta.getDouble('RA'), 

1420 matchmeta.getDouble('DEC'), geom.degrees) 

1421 rad = matchmeta.getDouble('RADIUS')*geom.degrees 

1422 refCat = refObjLoader.loadSkyCircle(ctrCoord, rad, filterName, epoch=epoch).refCat 

1423 elif "INNER_UPPER_LEFT_RA" in matchmeta: 

1424 # This is the sky box type (only triggers in the LoadReferenceObject class, not task) 

1425 # Only the outer box is required to be loaded to get the maximum region, all filtering 

1426 # will be done by the unpackMatches function, and no spatial filtering needs to be done 

1427 # by the refObjLoader 

1428 box = [] 

1429 for place in ("UPPER_LEFT", "UPPER_RIGHT", "LOWER_LEFT", "LOWER_RIGHT"): 

1430 coord = geom.SpherePoint(matchmeta.getDouble(f"OUTER_{place}_RA"), 

1431 matchmeta.getDouble(f"OUTER_{place}_DEC"), 

1432 geom.degrees).getVector() 

1433 box.append(coord) 

1434 outerBox = sphgeom.ConvexPolygon(box) 

1435 refCat = refObjLoader.loadRegion(outerBox, filterName=filterName, epoch=epoch).refCat 

1436 

1437 refCat.sort() 

1438 sourceCat.sort() 

1439 return afwTable.unpackMatches(matchCat, refCat, sourceCat) 

1440 

1441 

1442def applyProperMotionsImpl(log, catalog, epoch): 

1443 """Apply proper motion correction to a reference catalog. 

1444 

1445 Adjust position and position error in the ``catalog`` 

1446 for proper motion to the specified ``epoch``, 

1447 modifying the catalog in place. 

1448 

1449 Parameters 

1450 ---------- 

1451 log : `lsst.log.Log` 

1452 Log object to write to. 

1453 catalog : `lsst.afw.table.SimpleCatalog` 

1454 Catalog of positions, containing: 

1455 

1456 - Coordinates, retrieved by the table's coordinate key. 

1457 - ``coord_raErr`` : Error in Right Ascension (rad). 

1458 - ``coord_decErr`` : Error in Declination (rad). 

1459 - ``pm_ra`` : Proper motion in Right Ascension (rad/yr, 

1460 East positive) 

1461 - ``pm_raErr`` : Error in ``pm_ra`` (rad/yr), optional. 

1462 - ``pm_dec`` : Proper motion in Declination (rad/yr, 

1463 North positive) 

1464 - ``pm_decErr`` : Error in ``pm_dec`` (rad/yr), optional. 

1465 - ``epoch`` : Mean epoch of object (an astropy.time.Time) 

1466 epoch : `astropy.time.Time` 

1467 Epoch to which to correct proper motion. 

1468 """ 

1469 if "epoch" not in catalog.schema or "pm_ra" not in catalog.schema or "pm_dec" not in catalog.schema: 

1470 log.warn("Proper motion correction not available from catalog") 

1471 return 

1472 if not catalog.isContiguous(): 

1473 raise RuntimeError("Catalog must be contiguous") 

1474 catEpoch = astropy.time.Time(catalog["epoch"], scale="tai", format="mjd") 

1475 log.info("Correcting reference catalog for proper motion to %r", epoch) 

1476 # Use `epoch.tai` to make sure the time difference is in TAI 

1477 timeDiffsYears = (epoch.tai - catEpoch).to(astropy.units.yr).value 

1478 coordKey = catalog.table.getCoordKey() 

1479 # Compute the offset of each object due to proper motion 

1480 # as components of the arc of a great circle along RA and Dec 

1481 pmRaRad = catalog["pm_ra"] 

1482 pmDecRad = catalog["pm_dec"] 

1483 offsetsRaRad = pmRaRad*timeDiffsYears 

1484 offsetsDecRad = pmDecRad*timeDiffsYears 

1485 # Compute the corresponding bearing and arc length of each offset 

1486 # due to proper motion, and apply the offset 

1487 # The factor of 1e6 for computing bearing is intended as 

1488 # a reasonable scale for typical values of proper motion 

1489 # in order to avoid large errors for small values of proper motion; 

1490 # using the offsets is another option, but it can give 

1491 # needlessly large errors for short duration 

1492 offsetBearingsRad = numpy.arctan2(pmDecRad*1e6, pmRaRad*1e6) 

1493 offsetAmountsRad = numpy.hypot(offsetsRaRad, offsetsDecRad) 

1494 for record, bearingRad, amountRad in zip(catalog, offsetBearingsRad, offsetAmountsRad): 

1495 record.set(coordKey, 

1496 record.get(coordKey).offset(bearing=bearingRad*geom.radians, 

1497 amount=amountRad*geom.radians)) 

1498 # Increase error in RA and Dec based on error in proper motion 

1499 if "coord_raErr" in catalog.schema: 

1500 catalog["coord_raErr"] = numpy.hypot(catalog["coord_raErr"], 

1501 catalog["pm_raErr"]*timeDiffsYears) 

1502 if "coord_decErr" in catalog.schema: 

1503 catalog["coord_decErr"] = numpy.hypot(catalog["coord_decErr"], 

1504 catalog["pm_decErr"]*timeDiffsYears)