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=None) 

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 # Calling str on the epoch is a workaround for code that never worked anyway and was not used 

547 # see DM-17438 

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

549 return md 

550 

551 @staticmethod 

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

553 """Return metadata about the load 

554 

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

556 a normalized match list.) 

557 

558 Parameters 

559 ---------- 

560 coord : `lsst.geom.SpherePoint` 

561 ICRS center of a circle 

562 radius : `lsst.geom.angle` 

563 radius of a circle 

564 filterName : `str` or None 

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

566 photoCalib : None 

567 Deprecated, only included for api compatibility 

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

569 Epoch to which to correct proper motion and parallax, 

570 or None to not apply such corrections. 

571 

572 Returns 

573 ------- 

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

575 """ 

576 md = PropertyList() 

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

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

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

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

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

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

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

584 return md 

585 

586 @staticmethod 

587 def addFluxAliases(refCat, defaultFilter, filterReferenceMap): 

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

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

590 

591 Parameters 

592 ---------- 

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

594 Catalog of reference objects 

595 defaultFilter : `str` 

596 Name of the default reference filter 

597 filterReferenceMap : `dict` of `str` 

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

599 correspond to the name of the reference filter. 

600 

601 Returns 

602 ------- 

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

604 Reference catalog with columns added to track reference filters 

605 

606 Raises 

607 ------ 

608 `RuntimeError` 

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

610 reference filter map. 

611 """ 

612 refCat = ReferenceObjectLoader.remapReferenceCatalogSchema(refCat, 

613 filterNameList=filterReferenceMap.keys()) 

614 aliasMap = refCat.schema.getAliasMap() 

615 if filterReferenceMap is None: 

616 filterReferenceMap = {} 

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

618 filterReferenceMap.items()): 

619 if refFilterName: 

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

621 refFluxName = refFilterName + "_flux" 

622 if refFluxName not in refCat.schema: 

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

624 aliasMap.set(camFluxName, refFluxName) 

625 

626 refFluxErrName = refFluxName + "Err" 

627 camFluxErrName = camFluxName + "Err" 

628 aliasMap.set(camFluxErrName, refFluxErrName) 

629 

630 return refCat 

631 

632 @staticmethod 

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

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

635 columns defined the remaining function arguments. 

636 

637 Parameters 

638 ---------- 

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

640 Reference catalog to map to new catalog 

641 

642 Returns 

643 ------- 

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

645 Deep copy of input reference catalog with additional columns added 

646 """ 

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

648 mapper.addMinimalSchema(refCat.schema, True) 

649 mapper.editOutputSchema().disconnectAliases() 

650 if filterNameList: 

651 for filterName in filterNameList: 

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

653 type=numpy.float64, 

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

655 units="Jy" 

656 ) 

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

658 type=numpy.float64, 

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

660 units="Jy" 

661 ) 

662 

663 if position: 

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

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

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

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

668 

669 if photometric: 

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

671 type="Flag", 

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

673 "calibration", 

674 ) 

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

676 type="Flag", 

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

678 ) 

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

680 type="Flag", 

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

682 ) 

683 

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

685 expandedCat.setMetadata(refCat.getMetadata()) 

686 expandedCat.extend(refCat, mapper=mapper) 

687 

688 return expandedCat 

689 

690 

691def getRefFluxField(schema, filterName=None): 

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

693 

694 return the alias of "anyFilterMapsToThis", if present 

695 else if filterName is specified: 

696 return "*filterName*_camFlux" if present 

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

698 matches reference filter name) 

699 else throw RuntimeError 

700 else: 

701 return "camFlux", if present, 

702 else throw RuntimeError 

703 

704 Parameters 

705 ---------- 

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

707 Reference catalog schema. 

708 filterName : `str`, optional 

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

710 set in the refcat loader config. 

711 

712 Returns 

713 ------- 

714 fluxFieldName : `str` 

715 Name of flux field. 

716 

717 Raises 

718 ------ 

719 RuntimeError 

720 If an appropriate field is not found. 

721 """ 

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

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

724 try: 

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

726 except LookupError: 

727 pass # try the filterMap next 

728 

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

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

731 else: 

732 fluxFieldList = ["camFlux"] 

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

734 if fluxField in schema: 

735 return fluxField 

736 

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

738 

739 

740def getRefFluxKeys(schema, filterName=None): 

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

742 

743 Parameters 

744 ---------- 

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

746 Reference catalog schema. 

747 filterName : `str` 

748 Name of camera filter. 

749 

750 Returns 

751 ------- 

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

753 Two keys: 

754 

755 - flux key 

756 - flux error key, if present, else None 

757 

758 Raises 

759 ------ 

760 RuntimeError 

761 If flux field not found. 

762 """ 

763 fluxField = getRefFluxField(schema, filterName) 

764 fluxErrField = fluxField + "Err" 

765 fluxKey = schema[fluxField].asKey() 

766 try: 

767 fluxErrKey = schema[fluxErrField].asKey() 

768 except Exception: 

769 fluxErrKey = None 

770 return (fluxKey, fluxErrKey) 

771 

772 

773class LoadReferenceObjectsConfig(pexConfig.Config): 

774 pixelMargin = pexConfig.RangeField( 

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

776 dtype=int, 

777 default=300, 

778 min=0, 

779 ) 

780 defaultFilter = pexConfig.Field( 

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

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

783 dtype=str, 

784 default="", 

785 ) 

786 anyFilterMapsToThis = pexConfig.Field( 

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

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

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

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

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

792 dtype=str, 

793 default=None, 

794 optional=True 

795 ) 

796 filterMap = pexConfig.DictField( 

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

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

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

800 keytype=str, 

801 itemtype=str, 

802 default={}, 

803 ) 

804 requireProperMotion = pexConfig.Field( 

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

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

807 dtype=bool, 

808 default=False, 

809 ) 

810 

811 def validate(self): 

812 super().validate() 

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

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

815 raise pexConfig.FieldValidationError(LoadReferenceObjectsConfig.anyFilterMapsToThis, 

816 self, msg) 

817 

818 

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

820## @addtogroup LSST_task_documentation 

821## @{ 

822## @page LoadReferenceObjectsTask 

823## @ref LoadReferenceObjectsTask_ "LoadReferenceObjectsTask" 

824## @copybrief LoadReferenceObjectsTask 

825## @} 

826 

827 

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

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

830 

831 @anchor LoadReferenceObjectsTask_ 

832 

833 @section meas_algorithms_loadReferenceObjects_Contents Contents 

834 

835 - @ref meas_algorithms_loadReferenceObjects_Purpose 

836 - @ref meas_algorithms_loadReferenceObjects_Initialize 

837 - @ref meas_algorithms_loadReferenceObjects_IO 

838 - @ref meas_algorithms_loadReferenceObjects_Schema 

839 - @ref meas_algorithms_loadReferenceObjects_Config 

840 

841 @section meas_algorithms_loadReferenceObjects_Purpose Description 

842 

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

844 in a particular region of the sky. 

845 

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

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

848 

849 @section meas_algorithms_loadReferenceObjects_Initialize Task initialisation 

850 

851 @copydoc \_\_init\_\_ 

852 

853 @section meas_algorithms_loadReferenceObjects_IO Invoking the Task 

854 

855 @copydoc loadPixelBox 

856 

857 @section meas_algorithms_loadReferenceObjects_Schema Schema of the reference object catalog 

858 

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

860 (other fields may also be present). 

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

862 See also makeMinimalSchema. 

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

884 nan if unknown 

885 

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

887 information is available for at least some objects. 

888 Numeric values are `nan` if unknown: 

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

890 

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

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

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

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

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

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

897 

898 - parallax: parallax (rad, an Angle) 

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

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

901 

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

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

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

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

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

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

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

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

910 

911 @section meas_algorithms_loadReferenceObjects_Config Configuration parameters 

912 

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

914 Most subclasses will add configuration variables. 

915 """ 

916 ConfigClass = LoadReferenceObjectsConfig 

917 _DefaultName = "LoadReferenceObjects" 

918 

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

920 """Construct a LoadReferenceObjectsTask 

921 

922 Parameters 

923 ---------- 

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

925 Data butler, for access reference catalogs. 

926 """ 

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

928 self.butler = butler 

929 

930 @pipeBase.timeMethod 

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

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

933 

934 Parameters 

935 ---------- 

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

937 Bounding box for pixels. 

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

939 WCS; used to convert pixel positions to sky coordinates 

940 and vice-versa. 

941 filterName : `str` 

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

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

944 (which are not yet implemented). 

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

946 Calibration, or `None` if unknown. 

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

948 Epoch to which to correct proper motion and parallax, 

949 or None to not apply such corrections. 

950 

951 Returns 

952 ------- 

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

954 A Struct containing the following fields: 

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

956 A catalog of reference objects with the standard 

957 schema, as documented in the main doc string for 

958 `LoadReferenceObjects`. 

959 The catalog is guaranteed to be contiguous. 

960 fluxField : `str` 

961 Name of flux field for specified `filterName`. 

962 

963 Notes 

964 ----- 

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

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

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

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

969 """ 

970 circle = self._calculateCircle(bbox, wcs) 

971 

972 # find objects in circle 

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

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

975 loadRes = self.loadSkyCircle(circle.coord, circle.radius, filterName, centroids=True) 

976 refCat = loadRes.refCat 

977 numFound = len(refCat) 

978 

979 # trim objects outside bbox 

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

981 numTrimmed = numFound - len(refCat) 

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

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

984 

985 # make sure catalog is contiguous 

986 if not refCat.isContiguous(): 

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

988 

989 return loadRes 

990 

991 @abc.abstractmethod 

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

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

994 

995 Parameters 

996 ---------- 

997 ctrCoord : `lsst.geom.SpherePoint` 

998 ICRS center of search region. 

999 radius : `lsst.geom.Angle` 

1000 Radius of search region. 

1001 filterName : `str` (optional) 

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

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

1004 (which are not yet implemented). 

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

1006 Epoch to which to correct proper motion and parallax, 

1007 or None to not apply such corrections. 

1008 centroids : `bool` (optional) 

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

1010 these fields to exist. 

1011 

1012 Returns 

1013 ------- 

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

1015 A Struct containing the following fields: 

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

1017 A catalog of reference objects with the standard 

1018 schema, as documented in the main doc string for 

1019 `LoadReferenceObjects`. 

1020 The catalog is guaranteed to be contiguous. 

1021 fluxField : `str` 

1022 Name of flux field for specified `filterName`. 

1023 

1024 Notes 

1025 ----- 

1026 Note that subclasses are responsible for performing the proper motion 

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

1028 the catalog. 

1029 """ 

1030 return 

1031 

1032 @staticmethod 

1033 def _trimToBBox(refCat, bbox, wcs): 

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

1035 centroid and hasCentroid fields. 

1036 

1037 Parameters 

1038 ---------- 

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

1040 A catalog of objects. The schema must include fields 

1041 "coord", "centroid" and "hasCentroid". 

1042 The "coord" field is read. 

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

1044 bbox : `lsst.geom.Box2D` 

1045 Pixel region 

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

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

1048 

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

1050 """ 

1051 afwTable.updateRefCentroids(wcs, refCat) 

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

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

1054 for star in refCat: 

1055 point = star.get(centroidKey) 

1056 if bbox.contains(point): 

1057 retStarCat.append(star) 

1058 return retStarCat 

1059 

1060 def _addFluxAliases(self, schema): 

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

1062 

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

1064 camFlux: <defaultFilter>_flux 

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

1066 

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

1068 <camFilter>_camFlux: <refFilter>_flux 

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

1070 

1071 Parameters 

1072 ---------- 

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

1074 Schema for reference catalog. 

1075 

1076 Raises 

1077 ------ 

1078 RuntimeError 

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

1080 """ 

1081 aliasMap = schema.getAliasMap() 

1082 

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

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

1085 if refFluxName not in schema: 

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

1087 raise RuntimeError(msg) 

1088 aliasMap.set("anyFilterMapsToThis", refFluxName) 

1089 return # this is mutually exclusive with filterMap 

1090 

1091 def addAliasesForOneFilter(filterName, refFilterName): 

1092 """Add aliases for a single filter 

1093 

1094 Parameters 

1095 ---------- 

1096 filterName : `str` (optional) 

1097 Camera filter name. The resulting alias name is 

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

1099 is `None` or `""`. 

1100 refFilterName : `str` 

1101 Reference catalog filter name; the field 

1102 <refFilterName>_flux must exist. 

1103 """ 

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

1105 refFluxName = refFilterName + "_flux" 

1106 if refFluxName not in schema: 

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

1108 aliasMap.set(camFluxName, refFluxName) 

1109 refFluxErrName = refFluxName + "Err" 

1110 if refFluxErrName in schema: 

1111 camFluxErrName = camFluxName + "Err" 

1112 aliasMap.set(camFluxErrName, refFluxErrName) 

1113 

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

1115 addAliasesForOneFilter(None, self.config.defaultFilter) 

1116 

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

1118 addAliasesForOneFilter(filterName, refFilterName) 

1119 

1120 @staticmethod 

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

1122 addIsPhotometric=False, addIsResolved=False, 

1123 addIsVariable=False, coordErrDim=2, 

1124 addProperMotion=False, properMotionErrDim=2, 

1125 addParallax=False): 

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

1127 

1128 Parameters 

1129 ---------- 

1130 filterNameList : `list` of `str` 

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

1132 addIsPhotometric : `bool` 

1133 If True then add field "photometric". 

1134 addIsResolved : `bool` 

1135 If True then add field "resolved". 

1136 addIsVariable : `bool` 

1137 If True then add field "variable". 

1138 coordErrDim : `int` 

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

1140 

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

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

1143 addProperMotion : `bool` 

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

1145 properMotionErrDim : `int` 

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

1147 ignored if addProperMotion false: 

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

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

1150 addParallax : `bool` 

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

1152 and "parallax_flag". 

1153 

1154 Returns 

1155 ------- 

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

1157 Schema for reference catalog, an 

1158 `lsst.afw.table.SimpleCatalog`. 

1159 

1160 Notes 

1161 ----- 

1162 Reference catalogs support additional covariances, such as 

1163 covariance between RA and proper motion in declination, 

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

1165 calling this method. 

1166 """ 

1167 schema = afwTable.SimpleTable.makeMinimalSchema() 

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

1169 afwTable.Point2DKey.addFields( 

1170 schema, 

1171 "centroid", 

1172 "centroid on an exposure, if relevant", 

1173 "pixel", 

1174 ) 

1175 schema.addField( 

1176 field="hasCentroid", 

1177 type="Flag", 

1178 doc="is position known?", 

1179 ) 

1180 for filterName in filterNameList: 

1181 schema.addField( 

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

1183 type=numpy.float64, 

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

1185 units="nJy", 

1186 ) 

1187 for filterName in filterNameList: 

1188 schema.addField( 

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

1190 type=numpy.float64, 

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

1192 units="nJy", 

1193 ) 

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

1195 schema.addField( 

1196 field="photometric", 

1197 type="Flag", 

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

1199 ) 

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

1201 schema.addField( 

1202 field="resolved", 

1203 type="Flag", 

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

1205 ) 

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

1207 schema.addField( 

1208 field="variable", 

1209 type="Flag", 

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

1211 ) 

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

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

1214 if coordErrDim > 0: 

1215 afwTable.CovarianceMatrix2fKey.addFields( 

1216 schema=schema, 

1217 prefix="coord", 

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

1219 units=["rad", "rad"], 

1220 diagonalOnly=(coordErrDim == 2), 

1221 ) 

1222 

1223 if addProperMotion or addParallax: 

1224 schema.addField( 

1225 field="epoch", 

1226 type=numpy.float64, 

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

1228 units="day", 

1229 ) 

1230 

1231 if addProperMotion: 

1232 schema.addField( 

1233 field="pm_ra", 

1234 type="Angle", 

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

1236 units="rad/year", 

1237 ) 

1238 schema.addField( 

1239 field="pm_dec", 

1240 type="Angle", 

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

1242 units="rad/year", 

1243 ) 

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

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

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

1247 afwTable.CovarianceMatrix2fKey.addFields( 

1248 schema=schema, 

1249 prefix="pm", 

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

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

1252 diagonalOnly=(properMotionErrDim == 2), 

1253 ) 

1254 schema.addField( 

1255 field="pm_flag", 

1256 type="Flag", 

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

1258 ) 

1259 

1260 if addParallax: 

1261 schema.addField( 

1262 field="parallax", 

1263 type="Angle", 

1264 doc="parallax", 

1265 units="rad", 

1266 ) 

1267 schema.addField( 

1268 field="parallaxErr", 

1269 type="Angle", 

1270 doc="uncertainty in parallax", 

1271 units="rad", 

1272 ) 

1273 schema.addField( 

1274 field="parallax_flag", 

1275 type="Flag", 

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

1277 ) 

1278 return schema 

1279 

1280 def _calculateCircle(self, bbox, wcs): 

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

1282 

1283 Parameters 

1284 ---------- 

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

1286 Pixel bounding box. 

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

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

1289 

1290 Returns 

1291 ------- 

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

1293 A Struct containing: 

1294 

1295 - coord : `lsst.geom.SpherePoint` 

1296 ICRS center of the search region. 

1297 - radius : `lsst.geom.Angle` 

1298 Radius of the search region. 

1299 - bbox : `lsst.geom.Box2D` 

1300 Bounding box used to compute the circle. 

1301 """ 

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

1303 bbox.grow(self.config.pixelMargin) 

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

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

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

1307 

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

1309 """Return metadata about the load. 

1310 

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

1312 reconstituting a normalised match list. 

1313 

1314 Parameters 

1315 ---------- 

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

1317 Pixel bounding box. 

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

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

1320 filterName : `str` 

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

1322 filter. 

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

1324 Calibration, or `None` if unknown. 

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

1326 Epoch to which to correct proper motion and parallax, 

1327 or None to not apply such corrections. 

1328 

1329 Returns 

1330 ------- 

1331 metadata : lsst.daf.base.PropertyList 

1332 Metadata about the load. 

1333 """ 

1334 circle = self._calculateCircle(bbox, wcs) 

1335 return self.getMetadataCircle(circle.coord, circle.radius, filterName, photoCalib) 

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, '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)