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 

35import lsst.afw.table as afwTable 

36import lsst.pex.config as pexConfig 

37import lsst.pex.exceptions as pexExceptions 

38import lsst.pipe.base as pipeBase 

39import lsst.pex.exceptions as pexExcept 

40import lsst.log 

41from lsst import geom 

42from lsst import sphgeom 

43from lsst.daf.base import PropertyList 

44 

45 

46def 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 

57def 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()): 61 ↛ 62line 61 didn't jump to line 62, because the condition on line 61 was never true

62 return False 

63 return True 

64 

65 

66def 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: 81 ↛ 82line 81 didn't jump to line 82, because the condition on line 81 was never true

82 return None 

83 try: 

84 return md.getScalar("REFCAT_FORMAT_VERSION") 

85 except KeyError: 

86 return None 

87 

88 

89def 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 

151class _FilterCatalog: 

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

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

154 the class. 

155 

156 Parameters 

157 ---------- 

158 region : `lsst.sphgeom.Region` 

159 The spatial region which all objects should lie within 

160 """ 

161 def __init__(self, region): 

162 self.region = region 

163 

164 def __call__(self, refCat, catRegion): 

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

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

167 

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

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

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

171 

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

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

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

175 this catalog is then returned. 

176 

177 Parameters 

178 --------- 

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

180 SourceCatalog to be filtered. 

181 catRegion : `lsst.sphgeom.Region` 

182 Region in which the catalog was created 

183 """ 

184 if catRegion.isWithin(self.region): 

185 # no filtering needed, region completely contains refcat 

186 return refCat 

187 

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

189 for record in refCat: 

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

191 filteredRefCat.append(record) 

192 return filteredRefCat 

193 

194 

195class ReferenceObjectLoader: 

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

197 

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

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

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

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

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

203 call a corresponding method to load the reference objects. 

204 """ 

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

206 """ Constructs an instance of ReferenceObjectLoader 

207 

208 Parameters 

209 ---------- 

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

211 An iterable object of DataSetRefs which point to reference catalogs 

212 in a gen 3 repository 

213 refCats : Iterable of `lsst.daf.butler.DeferedDatasetHandle` 

214 Handles to load refCats on demand 

215 log : `lsst.log.Log` 

216 Logger object used to write out messages. If `None` (default) the default 

217 lsst logger will be used 

218 

219 """ 

220 self.dataIds = dataIds 

221 self.refCats = refCats 

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

223 self.config = config 

224 

225 @staticmethod 

226 def _makeBoxRegion(BBox, wcs, BBoxPadding): 

227 outerLocalBBox = geom.Box2D(BBox) 

228 innerLocalBBox = geom.Box2D(BBox) 

229 

230 # Grow the bounding box to make sure the spherical geometry bbox will contain 

231 # the same region, as non-padded boxes may contain different regions because of optical distortion. 

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

233 outerLocalBBox.grow(BBoxPadding) 

234 innerLocalBBox.grow(-1*BBoxPadding) 

235 

236 # Handle the fact that the inner bounding box shrunk to a zero sized region in at least one 

237 # dimension, in which case all reference catalogs must be checked fully against the input 

238 # bounding box 

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

240 innerSkyRegion = sphgeom.Box() 

241 else: 

242 innerBoxCorners = innerLocalBBox.getCorners() 

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

244 

245 innerSkyRegion = sphgeom.ConvexPolygon(innerSphCorners) 

246 

247 # Convert the corners of the box to sky coordinates 

248 outerBoxCorners = outerLocalBBox.getCorners() 

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

250 

251 outerSkyRegion = sphgeom.ConvexPolygon(outerSphCorners) 

252 

253 return innerSkyRegion, outerSkyRegion, innerSphCorners, outerSphCorners 

254 

255 def loadPixelBox(self, bbox, wcs, filterName=None, epoch=None, photoCalib=None, bboxPadding=100): 

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

257 

258 This algorithm works by creating a spherical box whose corners correspond 

259 to the WCS converted corners of the input bounding box (possibly padded). 

260 It then defines a filtering function which will look at a reference 

261 objects pixel position and accept objects that lie within the specified 

262 bounding box. 

263 

264 The spherical box region and filtering function are passed to the generic 

265 loadRegion method which will load and filter the reference objects from 

266 the datastore and return a single catalog containing all reference objects 

267 

268 Parameters 

269 ---------- 

270 bbox : `lsst.geom.box2I` 

271 Box which bounds a region in pixel space 

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

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

274 of pixels of the supplied bbox 

275 filterName : `str` 

276 Name of camera filter, or None or blank for the default filter 

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

278 Epoch to which to correct proper motion and parallax, 

279 or None to not apply such corrections. 

280 photoCalib : None 

281 Deprecated and ignored, only included for api compatibility 

282 bboxPadding : `int` 

283 Number describing how much to pad the input bbox by (in pixels), defaults 

284 to 100. This parameter is necessary because optical distortions in telescopes 

285 can cause a rectangular pixel grid to map into a non "rectangular" spherical 

286 region in sky coordinates. This padding is used to create a spherical 

287 "rectangle", which will for sure enclose the input box. This padding is only 

288 used to determine if the reference catalog for a sky patch will be loaded from 

289 the data store, this function will filter out objects which lie within the 

290 padded region but fall outside the input bounding box region. 

291 

292 Returns 

293 ------- 

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

295 Catalog containing reference objects inside the specified bounding box 

296 

297 Raises 

298 ------ 

299 `lsst.pex.exception.RuntimeError` 

300 Raised if no reference catalogs could be found for the specified region 

301 

302 `lsst.pex.exception.TypeError` 

303 Raised if the loaded reference catalogs do not have matching schemas 

304 """ 

305 innerSkyRegion, outerSkyRegion, _, _ = self._makeBoxRegion(bbox, wcs, bboxPadding) 

306 

307 def _filterFunction(refCat, region): 

308 # Add columns to the reference catalog relating to center positions and use afwTable 

309 # to populate those columns 

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

311 afwTable.updateRefCentroids(wcs, refCat) 

312 # no need to filter the catalog if it is sure that it is entirely contained in the region 

313 # defined by given bbox 

314 if innerSkyRegion.contains(region): 

315 return refCat 

316 # Create a new reference catalog, and populate it with records which fall inside the bbox 

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

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

319 for record in refCat: 

320 pixCoords = record[centroidKey] 

321 if bbox.contains(geom.Point2I(pixCoords)): 

322 filteredRefCat.append(record) 

323 return filteredRefCat 

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

325 

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

327 """ Load reference objects within a specified region 

328 

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

330 which intersect or are contained within the specified region. The reference 

331 catalogs which intersect but are not fully contained within the input region are 

332 further filtered by the specified filter function. This function will return a 

333 single source catalog containing all reference objects inside the specified region. 

334 

335 Parameters 

336 ---------- 

337 region : `lsst.sphgeom.Region` 

338 This can be any type that is derived from `lsst.sphgeom.region` and should 

339 define the spatial region for which reference objects are to be loaded. 

340 filtFunc : callable 

341 This optional parameter should be a callable object that takes a reference 

342 catalog and its corresponding region as parameters, filters the catalog by 

343 some criteria and returns the filtered reference catalog. If the value is 

344 left as the default (None) than an internal filter function is used which 

345 filters according to if a reference object falls within the input region. 

346 filterName : `str` 

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

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

349 Epoch to which to correct proper motion and parallax, 

350 or None to not apply such corrections. 

351 

352 Returns 

353 ------- 

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

355 Catalog containing reference objects which intersect the input region, 

356 filtered by the specified filter function 

357 

358 Raises 

359 ------ 

360 `lsst.pex.exception.RuntimeError` 

361 Raised if no reference catalogs could be found for the specified region 

362 

363 `lsst.pex.exception.TypeError` 

364 Raised if the loaded reference catalogs do not have matching schemas 

365 

366 """ 

367 regionBounding = region.getBoundingBox() 

368 self.log.info("Loading reference objects from region bounded by {}, {} lat lon".format( 

369 regionBounding.getLat(), regionBounding.getLon())) 

370 if filtFunc is None: 

371 filtFunc = _FilterCatalog(region) 

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

373 overlapList = [] 

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

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

376 # try the intersect operation in both directions 

377 try: 

378 intersects = dataId.region.intersects(region) 

379 except TypeError: 

380 intersects = region.intersects(dataId.region) 

381 

382 if intersects: 

383 overlapList.append((dataId, refCat)) 

384 

385 if len(overlapList) == 0: 

386 raise pexExceptions.RuntimeError("No reference tables could be found for input region") 

387 

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

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

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

391 

392 # Load in the remaining catalogs 

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

394 tmpCat = inputRefCat.get() 

395 

396 if tmpCat.schema != firstCat.schema: 

397 raise pexExceptions.TypeError("Reference catalogs have mismatching schemas") 

398 

399 filteredCat = filtFunc(tmpCat, dataId.region) 

400 refCat.extend(filteredCat) 

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

402 

403 self.log.debug(f"Trimmed {trimmedAmount} out of region objects, leaving {len(refCat)}") 

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

405 

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

407 if not refCat.isContiguous(): 

408 refCat = refCat.copy(deep=True) 

409 

410 if epoch is not None and "pm_ra" in refCat.schema: 

411 # check for a catalog in a non-standard format 

412 if isinstance(refCat.schema["pm_ra"].asKey(), lsst.afw.table.KeyAngle): 

413 applyProperMotionsImpl(self.log, refCat, epoch) 

414 else: 

415 self.log.warn("Catalog pm_ra field is not an Angle; not applying proper motion") 

416 

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

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

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

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

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

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

423 refCat = convertToNanojansky(refCat, self.log) 

424 

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

426 

427 # Add flux aliases 

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

429 

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

431 if not expandedCat.isContiguous(): 

432 expandedCat = expandedCat.copy(deep=True) 

433 

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

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

436 

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

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

439 

440 This method constructs a circular region from an input center and angular radius, 

441 loads reference catalogs which are contained in or intersect the circle, and 

442 filters reference catalogs which intersect down to objects which lie within 

443 the defined circle. 

444 

445 Parameters 

446 ---------- 

447 ctrCoord : `lsst.geom.SpherePoint` 

448 Point defining the center of the circular region 

449 radius : `lsst.geom.Angle` 

450 Defines the angular radius of the circular region 

451 filterName : `str` 

452 Name of camera filter, or None or blank for the default filter 

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

454 Epoch to which to correct proper motion and parallax, 

455 or None to not apply such corrections. 

456 

457 Returns 

458 ------- 

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

460 Catalog containing reference objects inside the specified bounding box 

461 

462 Raises 

463 ------ 

464 `lsst.pex.exception.RuntimeError` 

465 Raised if no reference catalogs could be found for the specified region 

466 

467 `lsst.pex.exception.TypeError` 

468 Raised if the loaded reference catalogs do not have matching schemas 

469 

470 """ 

471 centerVector = ctrCoord.getVector() 

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

473 circularRegion = sphgeom.Circle(centerVector, sphRadius) 

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

475 

476 def joinMatchListWithCatalog(self, matchCat, sourceCat): 

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

478 objects. 

479 

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

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

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

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

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

485 

486 Parameters 

487 ---------- 

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

489 Unpersisted packed match list. 

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

491 as returned by the astrometry tasks. 

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

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

494 by ID. 

495 

496 Returns 

497 ------- 

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

499 Match list. 

500 """ 

501 return joinMatchListWithCatalogImpl(self, matchCat, sourceCat) 

502 

503 @classmethod 

504 def getMetadataBox(cls, bbox, wcs, filterName=None, photoCalib=None, epoch=None, bboxPadding=100): 

505 """Return metadata about the load 

506 

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

508 reconstituting a normalised match list.) 

509 

510 Parameters 

511 ---------- 

512 bbox : `lsst.geom.Box2I` 

513 Bounding bos for the pixels 

514 wcs : `lsst.afw.geom.SkyWcs 

515 WCS object 

516 filterName : `str` or None 

517 filterName of the camera filter, or None or blank for the default filter 

518 photoCalib : None 

519 Deprecated, only included for api compatibility 

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

521 Epoch to which to correct proper motion and parallax, 

522 or None to not apply such corrections. 

523 bboxPadding : `int` 

524 Number describing how much to pad the input bbox by (in pixels), defaults 

525 to 100. This parameter is necessary because optical distortions in telescopes 

526 can cause a rectangular pixel grid to map into a non "rectangular" spherical 

527 region in sky coordinates. This padding is used to create a spherical 

528 "rectangle", which will for sure enclose the input box. This padding is only 

529 used to determine if the reference catalog for a sky patch will be loaded from 

530 the data store, this function will filter out objects which lie within the 

531 padded region but fall outside the input bounding box region. 

532 Returns 

533 ------- 

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

535 """ 

536 _, _, innerCorners, outerCorners = cls._makeBoxRegion(bbox, wcs, bboxPadding) 

537 md = PropertyList() 

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

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

540 corners): 

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

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

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

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

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

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

547 return md 

548 

549 @staticmethod 

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

551 """Return metadata about the load 

552 

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

554 a normalized match list.) 

555 

556 Parameters 

557 ---------- 

558 coord : `lsst.geom.SpherePoint` 

559 ICRS center of a circle 

560 radius : `lsst.geom.angle` 

561 radius of a circle 

562 filterName : `str` or None 

563 filterName of the camera filter, or None or blank for the default filter 

564 photoCalib : None 

565 Deprecated, only included for api compatibility 

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

567 Epoch to which to correct proper motion and parallax, 

568 or None to not apply such corrections. 

569 

570 Returns 

571 ------- 

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

573 """ 

574 md = PropertyList() 

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

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

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

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

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

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

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

582 return md 

583 

584 @staticmethod 

585 def addFluxAliases(refCat, defaultFilter, filterReferenceMap): 

586 """This function creates a new catalog containing the information of the input refCat 

587 as well as added flux columns and aliases between camera and reference flux. 

588 

589 Parameters 

590 ---------- 

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

592 Catalog of reference objects 

593 defaultFilter : `str` 

594 Name of the default reference filter 

595 filterReferenceMap : `dict` of `str` 

596 Dictionary with keys corresponding to a filter name, and values which 

597 correspond to the name of the reference filter. 

598 

599 Returns 

600 ------- 

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

602 Reference catalog with columns added to track reference filters 

603 

604 Raises 

605 ------ 

606 `RuntimeError` 

607 If specified reference filter name is not a filter specifed as a key in the 

608 reference filter map. 

609 """ 

610 refCat = ReferenceObjectLoader.remapReferenceCatalogSchema(refCat, 

611 filterNameList=filterReferenceMap.keys()) 

612 aliasMap = refCat.schema.getAliasMap() 

613 if filterReferenceMap is None: 

614 filterReferenceMap = {} 

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

616 filterReferenceMap.items()): 

617 if refFilterName: 

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

619 refFluxName = refFilterName + "_flux" 

620 if refFluxName not in refCat.schema: 

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

622 aliasMap.set(camFluxName, refFluxName) 

623 

624 refFluxErrName = refFluxName + "Err" 

625 camFluxErrName = camFluxName + "Err" 

626 aliasMap.set(camFluxErrName, refFluxErrName) 

627 

628 return refCat 

629 

630 @staticmethod 

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

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

633 columns defined the remaining function arguments. 

634 

635 Parameters 

636 ---------- 

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

638 Reference catalog to map to new catalog 

639 

640 Returns 

641 ------- 

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

643 Deep copy of input reference catalog with additional columns added 

644 """ 

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

646 mapper.addMinimalSchema(refCat.schema, True) 

647 mapper.editOutputSchema().disconnectAliases() 

648 if filterNameList: 

649 for filterName in filterNameList: 

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

651 type=numpy.float64, 

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

653 units="Jy" 

654 ) 

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

656 type=numpy.float64, 

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

658 units="Jy" 

659 ) 

660 

661 if position: 

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

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

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

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

666 

667 if photometric: 

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

669 type="Flag", 

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

671 "calibration", 

672 ) 

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

674 type="Flag", 

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

676 ) 

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

678 type="Flag", 

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

680 ) 

681 

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

683 expandedCat.setMetadata(refCat.getMetadata()) 

684 expandedCat.extend(refCat, mapper=mapper) 

685 

686 return expandedCat 

687 

688 

689def getRefFluxField(schema, filterName=None): 

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

691 

692 return the alias of "anyFilterMapsToThis", if present 

693 else if filterName is specified: 

694 return "*filterName*_camFlux" if present 

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

696 matches reference filter name) 

697 else throw RuntimeError 

698 else: 

699 return "camFlux", if present, 

700 else throw RuntimeError 

701 

702 Parameters 

703 ---------- 

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

705 Reference catalog schema. 

706 filterName : `str`, optional 

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

708 set in the refcat loader config. 

709 

710 Returns 

711 ------- 

712 fluxFieldName : `str` 

713 Name of flux field. 

714 

715 Raises 

716 ------ 

717 RuntimeError 

718 If an appropriate field is not found. 

719 """ 

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

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

722 try: 

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

724 except LookupError: 

725 pass # try the filterMap next 

726 

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

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

729 else: 

730 fluxFieldList = ["camFlux"] 

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

732 if fluxField in schema: 

733 return fluxField 

734 

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

736 

737 

738def getRefFluxKeys(schema, filterName=None): 

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

740 

741 Parameters 

742 ---------- 

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

744 Reference catalog schema. 

745 filterName : `str` 

746 Name of camera filter. 

747 

748 Returns 

749 ------- 

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

751 Two keys: 

752 

753 - flux key 

754 - flux error key, if present, else None 

755 

756 Raises 

757 ------ 

758 RuntimeError 

759 If flux field not found. 

760 """ 

761 fluxField = getRefFluxField(schema, filterName) 

762 fluxErrField = fluxField + "Err" 

763 fluxKey = schema[fluxField].asKey() 

764 try: 

765 fluxErrKey = schema[fluxErrField].asKey() 

766 except Exception: 

767 fluxErrKey = None 

768 return (fluxKey, fluxErrKey) 

769 

770 

771class LoadReferenceObjectsConfig(pexConfig.Config): 

772 pixelMargin = pexConfig.RangeField( 

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

774 dtype=int, 

775 default=300, 

776 min=0, 

777 ) 

778 defaultFilter = pexConfig.Field( 

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

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

781 dtype=str, 

782 default="", 

783 ) 

784 anyFilterMapsToThis = pexConfig.Field( 

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

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

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

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

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

790 dtype=str, 

791 default=None, 

792 optional=True 

793 ) 

794 filterMap = pexConfig.DictField( 

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

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

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

798 keytype=str, 

799 itemtype=str, 

800 default={}, 

801 ) 

802 requireProperMotion = pexConfig.Field( 

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

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

805 dtype=bool, 

806 default=False, 

807 ) 

808 

809 def validate(self): 

810 super().validate() 

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

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

813 raise pexConfig.FieldValidationError(LoadReferenceObjectsConfig.anyFilterMapsToThis, 

814 self, msg) 

815 

816 

817# The following comment block adds a link to this task from the Task Documentation page. 

818## @addtogroup LSST_task_documentation 

819## @{ 

820## @page LoadReferenceObjectsTask 

821## @ref LoadReferenceObjectsTask_ "LoadReferenceObjectsTask" 

822## @copybrief LoadReferenceObjectsTask 

823## @} 

824 

825 

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

827 r"""!Abstract base class to load objects from reference catalogs 

828 

829 @anchor LoadReferenceObjectsTask_ 

830 

831 @section meas_algorithms_loadReferenceObjects_Contents Contents 

832 

833 - @ref meas_algorithms_loadReferenceObjects_Purpose 

834 - @ref meas_algorithms_loadReferenceObjects_Initialize 

835 - @ref meas_algorithms_loadReferenceObjects_IO 

836 - @ref meas_algorithms_loadReferenceObjects_Schema 

837 - @ref meas_algorithms_loadReferenceObjects_Config 

838 

839 @section meas_algorithms_loadReferenceObjects_Purpose Description 

840 

841 Abstract base class for tasks that load objects from a reference catalog 

842 in a particular region of the sky. 

843 

844 Implementations must subclass this class, override the loadSkyCircle method, 

845 and will typically override the value of ConfigClass with a task-specific config class. 

846 

847 @section meas_algorithms_loadReferenceObjects_Initialize Task initialisation 

848 

849 @copydoc \_\_init\_\_ 

850 

851 @section meas_algorithms_loadReferenceObjects_IO Invoking the Task 

852 

853 @copydoc loadPixelBox 

854 

855 @section meas_algorithms_loadReferenceObjects_Schema Schema of the reference object catalog 

856 

857 Reference object catalogs are instances of lsst.afw.table.SimpleCatalog with the following schema 

858 (other fields may also be present). 

859 The units use astropy quantity conventions, so a 2 suffix means squared. 

860 See also makeMinimalSchema. 

861 - coord: ICRS position of star on sky (an lsst.geom.SpherePoint) 

862 - centroid: position of star on an exposure, if relevant (an lsst.afw.Point2D) 

863 - hasCentroid: is centroid usable? (a Flag) 

864 - *referenceFilterName*_flux: brightness in the specified reference catalog filter (nJy) 

865 Note: you can use astropy.units to convert from AB Magnitude to nJy: 

866 `u.Magnitude(value, u.ABmag).to_value(u.nJy)` 

867 - *referenceFilterName*_fluxErr (optional): brightness standard deviation (nJy); 

868 omitted if no data is available; possibly nan if data is available for some objects but not others 

869 - camFlux: brightness in default camera filter (nJy); omitted if defaultFilter not specified 

870 - camFluxErr: brightness standard deviation for default camera filter; 

871 omitted if defaultFilter not specified or standard deviation not available that filter 

872 - *cameraFilterName*_camFlux: brightness in specified camera filter (nJy) 

873 - *cameraFilterName*_camFluxErr (optional): brightness standard deviation 

874 in specified camera filter (nJy); omitted if no data is available; 

875 possibly nan if data is available for some objects but not others 

876 - photometric (optional): is the object usable for photometric calibration? (a Flag) 

877 - resolved (optional): is the object spatially resolved? (a Flag) 

878 - variable (optional): does the object have variable brightness? (a Flag) 

879 - coord_raErr: uncertainty in `coord` along the direction of right ascension (radian, an Angle) 

880 = uncertainty in ra * cos(dec); nan if unknown 

881 - coord_decErr: uncertainty in `coord` along the direction of declination (radian, an Angle); 

882 nan if unknown 

883 

884 The following are optional; fields should only be present if the 

885 information is available for at least some objects. 

886 Numeric values are `nan` if unknown: 

887 - epoch: date of observation as TAI MJD (day) 

888 

889 - pm_ra: proper motion along the direction of right ascension (rad/year, an Angle) = dra/dt * cos(dec) 

890 - pm_dec: proper motion along the direction of declination (rad/year, and Angle) 

891 - pm_raErr: uncertainty in `pm_ra` (rad/year) 

892 - pm_decErr: uncertainty in `pm_dec` (rad/year) 

893 - pm_ra_dec_Cov: covariance between pm_ra and pm_dec (rad2/year2) 

894 - pm_flag: set if proper motion, error or covariance is bad 

895 

896 - parallax: parallax (rad, an Angle) 

897 - parallaxErr: uncertainty in `parallax` (rad) 

898 - parallax_flag: set if parallax value or parallaxErr is bad 

899 

900 - coord_ra_pm_ra_Cov: covariance between coord_ra and pm_ra (rad2/year) 

901 - coord_ra_pm_dec_Cov: covariance between coord_ra and pm_dec (rad2/year) 

902 - coord_ra_parallax_Cov: covariance between coord_ra and parallax (rad2/year) 

903 - coord_dec_pm_ra_Cov: covariance between coord_dec and pm_ra (rad2/year) 

904 - coord_dec_pm_dec_Cov: covariance between coord_dec and pm_dec (rad2/year) 

905 - coord_dec_parallax_Cov: covariance between coord_dec and parallax (rad2/year) 

906 - pm_ra_parallax_Cov: covariance between pm_ra and parallax (rad2/year) 

907 - pm_dec_parallax_Cov: covariance between pm_dec and parallax (rad2/year) 

908 

909 @section meas_algorithms_loadReferenceObjects_Config Configuration parameters 

910 

911 See @ref LoadReferenceObjectsConfig for a base set of configuration parameters. 

912 Most subclasses will add configuration variables. 

913 """ 

914 ConfigClass = LoadReferenceObjectsConfig 

915 _DefaultName = "LoadReferenceObjects" 

916 

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

918 """Construct a LoadReferenceObjectsTask 

919 

920 Parameters 

921 ---------- 

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

923 Data butler, for access reference catalogs. 

924 """ 

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

926 self.butler = butler 

927 

928 @pipeBase.timeMethod 

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

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

931 

932 Parameters 

933 ---------- 

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

935 Bounding box for pixels. 

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

937 WCS; used to convert pixel positions to sky coordinates 

938 and vice-versa. 

939 filterName : `str` 

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

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

942 (which are not yet implemented). 

943 photoCalib : `lsst.afw.image.PhotoCalib` (optional) 

944 Calibration, or `None` if unknown. 

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

946 Epoch to which to correct proper motion and parallax, 

947 or None to not apply such corrections. 

948 

949 Returns 

950 ------- 

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

952 A Struct containing the following fields: 

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

954 A catalog of reference objects with the standard 

955 schema, as documented in the main doc string for 

956 `LoadReferenceObjects`. 

957 The catalog is guaranteed to be contiguous. 

958 fluxField : `str` 

959 Name of flux field for specified `filterName`. 

960 

961 Notes 

962 ----- 

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

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

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

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

967 """ 

968 circle = self._calculateCircle(bbox, wcs) 

969 

970 # find objects in circle 

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

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

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

974 centroids=True) 

975 refCat = loadRes.refCat 

976 numFound = len(refCat) 

977 

978 # trim objects outside bbox 

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

980 numTrimmed = numFound - len(refCat) 

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

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

983 

984 # make sure catalog is contiguous 

985 if not refCat.isContiguous(): 

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

987 

988 return loadRes 

989 

990 @abc.abstractmethod 

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

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

993 

994 Parameters 

995 ---------- 

996 ctrCoord : `lsst.geom.SpherePoint` 

997 ICRS center of search region. 

998 radius : `lsst.geom.Angle` 

999 Radius of search region. 

1000 filterName : `str` (optional) 

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

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

1003 (which are not yet implemented). 

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

1005 Epoch to which to correct proper motion and parallax, 

1006 or None to not apply such corrections. 

1007 centroids : `bool` (optional) 

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

1009 these fields to exist. 

1010 

1011 Returns 

1012 ------- 

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

1014 A Struct containing the following fields: 

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

1016 A catalog of reference objects with the standard 

1017 schema, as documented in the main doc string for 

1018 `LoadReferenceObjects`. 

1019 The catalog is guaranteed to be contiguous. 

1020 fluxField : `str` 

1021 Name of flux field for specified `filterName`. 

1022 

1023 Notes 

1024 ----- 

1025 Note that subclasses are responsible for performing the proper motion 

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

1027 the catalog. 

1028 """ 

1029 return 

1030 

1031 @staticmethod 

1032 def _trimToBBox(refCat, bbox, wcs): 

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

1034 centroid and hasCentroid fields. 

1035 

1036 Parameters 

1037 ---------- 

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

1039 A catalog of objects. The schema must include fields 

1040 "coord", "centroid" and "hasCentroid". 

1041 The "coord" field is read. 

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

1043 bbox : `lsst.geom.Box2D` 

1044 Pixel region 

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

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

1047 

1048 @return a catalog of reference objects in bbox, with centroid and hasCentroid fields set 

1049 """ 

1050 afwTable.updateRefCentroids(wcs, refCat) 

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

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

1053 for star in refCat: 

1054 point = star.get(centroidKey) 

1055 if bbox.contains(point): 

1056 retStarCat.append(star) 

1057 return retStarCat 

1058 

1059 def _addFluxAliases(self, schema): 

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

1061 

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

1063 camFlux: <defaultFilter>_flux 

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

1065 

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

1067 <camFilter>_camFlux: <refFilter>_flux 

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

1069 

1070 Parameters 

1071 ---------- 

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

1073 Schema for reference catalog. 

1074 

1075 Raises 

1076 ------ 

1077 RuntimeError 

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

1079 """ 

1080 aliasMap = schema.getAliasMap() 

1081 

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

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

1084 if refFluxName not in schema: 

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

1086 raise RuntimeError(msg) 

1087 aliasMap.set("anyFilterMapsToThis", refFluxName) 

1088 return # this is mutually exclusive with filterMap 

1089 

1090 def addAliasesForOneFilter(filterName, refFilterName): 

1091 """Add aliases for a single filter 

1092 

1093 Parameters 

1094 ---------- 

1095 filterName : `str` (optional) 

1096 Camera filter name. The resulting alias name is 

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

1098 is `None` or `""`. 

1099 refFilterName : `str` 

1100 Reference catalog filter name; the field 

1101 <refFilterName>_flux must exist. 

1102 """ 

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

1104 refFluxName = refFilterName + "_flux" 

1105 if refFluxName not in schema: 

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

1107 aliasMap.set(camFluxName, refFluxName) 

1108 refFluxErrName = refFluxName + "Err" 

1109 if refFluxErrName in schema: 

1110 camFluxErrName = camFluxName + "Err" 

1111 aliasMap.set(camFluxErrName, refFluxErrName) 

1112 

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

1114 addAliasesForOneFilter(None, self.config.defaultFilter) 

1115 

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

1117 addAliasesForOneFilter(filterName, refFilterName) 

1118 

1119 @staticmethod 

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

1121 addIsPhotometric=False, addIsResolved=False, 

1122 addIsVariable=False, coordErrDim=2, 

1123 addProperMotion=False, properMotionErrDim=2, 

1124 addParallax=False): 

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

1126 

1127 Parameters 

1128 ---------- 

1129 filterNameList : `list` of `str` 

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

1131 addIsPhotometric : `bool` 

1132 If True then add field "photometric". 

1133 addIsResolved : `bool` 

1134 If True then add field "resolved". 

1135 addIsVariable : `bool` 

1136 If True then add field "variable". 

1137 coordErrDim : `int` 

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

1139 

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

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

1142 addProperMotion : `bool` 

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

1144 properMotionErrDim : `int` 

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

1146 ignored if addProperMotion false: 

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

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

1149 addParallax : `bool` 

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

1151 and "parallax_flag". 

1152 

1153 Returns 

1154 ------- 

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

1156 Schema for reference catalog, an 

1157 `lsst.afw.table.SimpleCatalog`. 

1158 

1159 Notes 

1160 ----- 

1161 Reference catalogs support additional covariances, such as 

1162 covariance between RA and proper motion in declination, 

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

1164 calling this method. 

1165 """ 

1166 schema = afwTable.SimpleTable.makeMinimalSchema() 

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

1168 afwTable.Point2DKey.addFields( 

1169 schema, 

1170 "centroid", 

1171 "centroid on an exposure, if relevant", 

1172 "pixel", 

1173 ) 

1174 schema.addField( 

1175 field="hasCentroid", 

1176 type="Flag", 

1177 doc="is position known?", 

1178 ) 

1179 for filterName in filterNameList: 

1180 schema.addField( 

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

1182 type=numpy.float64, 

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

1184 units="nJy", 

1185 ) 

1186 for filterName in filterNameList: 

1187 schema.addField( 

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

1189 type=numpy.float64, 

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

1191 units="nJy", 

1192 ) 

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

1194 schema.addField( 

1195 field="photometric", 

1196 type="Flag", 

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

1198 ) 

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

1200 schema.addField( 

1201 field="resolved", 

1202 type="Flag", 

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

1204 ) 

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

1206 schema.addField( 

1207 field="variable", 

1208 type="Flag", 

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

1210 ) 

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

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

1213 if coordErrDim > 0: 

1214 afwTable.CovarianceMatrix2fKey.addFields( 

1215 schema=schema, 

1216 prefix="coord", 

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

1218 units=["rad", "rad"], 

1219 diagonalOnly=(coordErrDim == 2), 

1220 ) 

1221 

1222 if addProperMotion or addParallax: 

1223 schema.addField( 

1224 field="epoch", 

1225 type=numpy.float64, 

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

1227 units="day", 

1228 ) 

1229 

1230 if addProperMotion: 

1231 schema.addField( 

1232 field="pm_ra", 

1233 type="Angle", 

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

1235 units="rad/year", 

1236 ) 

1237 schema.addField( 

1238 field="pm_dec", 

1239 type="Angle", 

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

1241 units="rad/year", 

1242 ) 

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

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

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

1246 afwTable.CovarianceMatrix2fKey.addFields( 

1247 schema=schema, 

1248 prefix="pm", 

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

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

1251 diagonalOnly=(properMotionErrDim == 2), 

1252 ) 

1253 schema.addField( 

1254 field="pm_flag", 

1255 type="Flag", 

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

1257 ) 

1258 

1259 if addParallax: 

1260 schema.addField( 

1261 field="parallax", 

1262 type="Angle", 

1263 doc="parallax", 

1264 units="rad", 

1265 ) 

1266 schema.addField( 

1267 field="parallaxErr", 

1268 type="Angle", 

1269 doc="uncertainty in parallax", 

1270 units="rad", 

1271 ) 

1272 schema.addField( 

1273 field="parallax_flag", 

1274 type="Flag", 

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

1276 ) 

1277 return schema 

1278 

1279 def _calculateCircle(self, bbox, wcs): 

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

1281 

1282 Parameters 

1283 ---------- 

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

1285 Pixel bounding box. 

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

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

1288 

1289 Returns 

1290 ------- 

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

1292 A Struct containing: 

1293 

1294 - coord : `lsst.geom.SpherePoint` 

1295 ICRS center of the search region. 

1296 - radius : `lsst.geom.Angle` 

1297 Radius of the search region. 

1298 - bbox : `lsst.geom.Box2D` 

1299 Bounding box used to compute the circle. 

1300 """ 

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

1302 bbox.grow(self.config.pixelMargin) 

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

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

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

1306 

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

1308 """Return metadata about the load. 

1309 

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

1311 reconstituting a normalised match list. 

1312 

1313 Parameters 

1314 ---------- 

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

1316 Pixel bounding box. 

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

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

1319 filterName : `str` 

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

1321 filter. 

1322 photoCalib : `lsst.afw.image.PhotoCalib` (optional) 

1323 Calibration, or `None` if unknown. 

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

1325 Epoch to which to correct proper motion and parallax, 

1326 or None to not apply such corrections. 

1327 

1328 Returns 

1329 ------- 

1330 metadata : lsst.daf.base.PropertyList 

1331 Metadata about the load. 

1332 """ 

1333 circle = self._calculateCircle(bbox, wcs) 

1334 return self.getMetadataCircle(circle.coord, circle.radius, filterName, photoCalib=photoCalib, 

1335 epoch=epoch) 

1336 

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

1338 """Return metadata about the load. 

1339 

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

1341 reconstituting a normalised match list. 

1342 

1343 Parameters 

1344 ---------- 

1345 coord : `lsst.geom.SpherePoint` 

1346 ICRS center of the search region. 

1347 radius : `lsst.geom.Angle` 

1348 Radius of the search region. 

1349 filterName : `str` 

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

1351 filter. 

1352 photoCalib : `lsst.afw.image.PhotoCalib` (optional) 

1353 Calibration, or `None` if unknown. 

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

1355 Epoch to which to correct proper motion and parallax, 

1356 or None to not apply such corrections. 

1357 

1358 Returns 

1359 ------- 

1360 metadata : lsst.daf.base.PropertyList 

1361 Metadata about the load 

1362 """ 

1363 md = PropertyList() 

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

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

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

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

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

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

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

1371 return md 

1372 

1373 def joinMatchListWithCatalog(self, matchCat, sourceCat): 

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

1375 objects. 

1376 

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

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

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

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

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

1382 

1383 Parameters 

1384 ---------- 

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

1386 Unperisted packed match list. 

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

1388 as returned by the astrometry tasks. 

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

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

1391 by ID. 

1392 

1393 Returns 

1394 ------- 

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

1396 Match list. 

1397 """ 

1398 return joinMatchListWithCatalogImpl(self, matchCat, sourceCat) 

1399 

1400 def applyProperMotions(self, catalog, epoch): 

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

1402 

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

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

1405 modifying the catalog in place. 

1406 

1407 Parameters 

1408 ---------- 

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

1410 Catalog of positions, containing: 

1411 

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

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

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

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

1416 East positive) 

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

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

1419 North positive) 

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

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

1422 epoch : `astropy.time.Time` 

1423 Epoch to which to correct proper motion, 

1424 """ 

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

1426 if self.config.requireProperMotion: 

1427 raise RuntimeError("Proper motion correction required but not available from catalog") 

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

1429 return 

1430 applyProperMotionsImpl(self.log, catalog, epoch) 

1431 

1432 

1433def joinMatchListWithCatalogImpl(refObjLoader, matchCat, sourceCat): 

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

1435 objects. 

1436 

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

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

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

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

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

1442 

1443 Parameters 

1444 ---------- 

1445 refObjLoader 

1446 Reference object loader to use in getting reference objects 

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

1448 Unperisted packed match list. 

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

1450 as returned by the astrometry tasks. 

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

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

1453 by ID. 

1454 

1455 Returns 

1456 ------- 

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

1458 Match list. 

1459 """ 

1460 matchmeta = matchCat.table.getMetadata() 

1461 version = matchmeta.getInt('SMATCHV') 

1462 if version != 1: 

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

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

1465 try: 

1466 epoch = matchmeta.getDouble('EPOCH') 

1467 except (pexExcept.NotFoundError, pexExcept.TypeError): 

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

1469 if 'RADIUS' in matchmeta: 

1470 # This is a circle style metadata, call loadSkyCircle 

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

1472 matchmeta.getDouble('DEC'), lsst.geom.degrees) 

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

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

1475 elif "INNER_UPPER_LEFT_RA" in matchmeta: 

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

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

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

1479 # by the refObjLoader 

1480 box = [] 

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

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

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

1484 lsst.geom.degrees).getVector() 

1485 box.append(coord) 

1486 outerBox = sphgeom.ConvexPolygon(box) 

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

1488 

1489 refCat.sort() 

1490 sourceCat.sort() 

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

1492 

1493 

1494def applyProperMotionsImpl(log, catalog, epoch): 

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

1496 

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

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

1499 modifying the catalog in place. 

1500 

1501 Parameters 

1502 ---------- 

1503 log : `lsst.log.Log` 

1504 Log object to write to. 

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

1506 Catalog of positions, containing: 

1507 

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

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

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

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

1512 East positive) 

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

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

1515 North positive) 

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

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

1518 epoch : `astropy.time.Time` 

1519 Epoch to which to correct proper motion. 

1520 """ 

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

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

1523 return 

1524 if not catalog.isContiguous(): 

1525 raise RuntimeError("Catalog must be contiguous") 

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

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

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

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

1530 coordKey = catalog.table.getCoordKey() 

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

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

1533 pmRaRad = catalog["pm_ra"] 

1534 pmDecRad = catalog["pm_dec"] 

1535 offsetsRaRad = pmRaRad*timeDiffsYears 

1536 offsetsDecRad = pmDecRad*timeDiffsYears 

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

1538 # due to proper motion, and apply the offset 

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

1540 # a reasonable scale for typical values of proper motion 

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

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

1543 # needlessly large errors for short duration 

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

1545 offsetAmountsRad = numpy.hypot(offsetsRaRad, offsetsDecRad) 

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

1547 record.set(coordKey, 

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

1549 amount=amountRad*lsst.geom.radians)) 

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

1551 if "coord_raErr" in catalog.schema: 

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

1553 catalog["pm_raErr"]*timeDiffsYears) 

1554 if "coord_decErr" in catalog.schema: 

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

1556 catalog["pm_decErr"]*timeDiffsYears)