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()): 

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: 

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 if filterName is specified: 

695 return *filterName*_camFlux if present 

696 else return *filterName*_flux if present (camera filter name 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` 

707 Name of camera filter. 

708 

709 Returns 

710 ------- 

711 fluxFieldName : `str` 

712 Name of flux field. 

713 

714 Raises 

715 ------ 

716 RuntimeError 

717 If an appropriate field is not found. 

718 """ 

719 if not isinstance(schema, afwTable.Schema): 

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

721 if filterName: 

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

723 else: 

724 fluxFieldList = ["camFlux"] 

725 for fluxField in fluxFieldList: 

726 if fluxField in schema: 

727 return fluxField 

728 

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

730 

731 

732def getRefFluxKeys(schema, filterName=None): 

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

734 

735 Parameters 

736 ---------- 

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

738 Reference catalog schema. 

739 filterName : `str` 

740 Name of camera filter. 

741 

742 Returns 

743 ------- 

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

745 Two keys: 

746 

747 - flux key 

748 - flux error key, if present, else None 

749 

750 Raises 

751 ------ 

752 RuntimeError 

753 If flux field not found. 

754 """ 

755 fluxField = getRefFluxField(schema, filterName) 

756 fluxErrField = fluxField + "Err" 

757 fluxKey = schema[fluxField].asKey() 

758 try: 

759 fluxErrKey = schema[fluxErrField].asKey() 

760 except Exception: 

761 fluxErrKey = None 

762 return (fluxKey, fluxErrKey) 

763 

764 

765class LoadReferenceObjectsConfig(pexConfig.Config): 

766 pixelMargin = pexConfig.RangeField( 

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

768 dtype=int, 

769 default=300, 

770 min=0, 

771 ) 

772 defaultFilter = pexConfig.Field( 

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

774 "if blank then filter must be specified in exposure", 

775 dtype=str, 

776 default="", 

777 ) 

778 filterMap = pexConfig.DictField( 

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

780 "each reference filter must exist", 

781 keytype=str, 

782 itemtype=str, 

783 default={}, 

784 ) 

785 requireProperMotion = pexConfig.Field( 

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

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

788 dtype=bool, 

789 default=False, 

790 ) 

791 

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

793## @addtogroup LSST_task_documentation 

794## @{ 

795## @page LoadReferenceObjectsTask 

796## @ref LoadReferenceObjectsTask_ "LoadReferenceObjectsTask" 

797## @copybrief LoadReferenceObjectsTask 

798## @} 

799 

800 

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

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

803 

804 @anchor LoadReferenceObjectsTask_ 

805 

806 @section meas_algorithms_loadReferenceObjects_Contents Contents 

807 

808 - @ref meas_algorithms_loadReferenceObjects_Purpose 

809 - @ref meas_algorithms_loadReferenceObjects_Initialize 

810 - @ref meas_algorithms_loadReferenceObjects_IO 

811 - @ref meas_algorithms_loadReferenceObjects_Schema 

812 - @ref meas_algorithms_loadReferenceObjects_Config 

813 

814 @section meas_algorithms_loadReferenceObjects_Purpose Description 

815 

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

817 in a particular region of the sky. 

818 

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

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

821 

822 @section meas_algorithms_loadReferenceObjects_Initialize Task initialisation 

823 

824 @copydoc \_\_init\_\_ 

825 

826 @section meas_algorithms_loadReferenceObjects_IO Invoking the Task 

827 

828 @copydoc loadPixelBox 

829 

830 @section meas_algorithms_loadReferenceObjects_Schema Schema of the reference object catalog 

831 

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

833 (other fields may also be present). 

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

835 See also makeMinimalSchema. 

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

857 nan if unknown 

858 

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

860 information is available for at least some objects. 

861 Numeric values are `nan` if unknown: 

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

863 

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

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

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

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

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

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

870 

871 - parallax: parallax (rad, an Angle) 

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

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

874 

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

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

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

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

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

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

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

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

883 

884 @section meas_algorithms_loadReferenceObjects_Config Configuration parameters 

885 

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

887 Most subclasses will add configuration variables. 

888 """ 

889 ConfigClass = LoadReferenceObjectsConfig 

890 _DefaultName = "LoadReferenceObjects" 

891 

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

893 """Construct a LoadReferenceObjectsTask 

894 

895 Parameters 

896 ---------- 

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

898 Data butler, for access reference catalogs. 

899 """ 

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

901 self.butler = butler 

902 

903 @pipeBase.timeMethod 

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

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

906 

907 Parameters 

908 ---------- 

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

910 Bounding box for pixels. 

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

912 WCS; used to convert pixel positions to sky coordinates 

913 and vice-versa. 

914 filterName : `str` 

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

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

917 (which are not yet implemented). 

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

919 Calibration, or `None` if unknown. 

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

921 Epoch to which to correct proper motion and parallax, 

922 or None to not apply such corrections. 

923 

924 Returns 

925 ------- 

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

927 A Struct containing the following fields: 

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

929 A catalog of reference objects with the standard 

930 schema, as documented in the main doc string for 

931 `LoadReferenceObjects`. 

932 The catalog is guaranteed to be contiguous. 

933 fluxField : `str` 

934 Name of flux field for specified `filterName`. 

935 

936 Notes 

937 ----- 

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

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

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

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

942 """ 

943 circle = self._calculateCircle(bbox, wcs) 

944 

945 # find objects in circle 

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

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

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

949 refCat = loadRes.refCat 

950 numFound = len(refCat) 

951 

952 # trim objects outside bbox 

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

954 numTrimmed = numFound - len(refCat) 

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

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

957 

958 # make sure catalog is contiguous 

959 if not refCat.isContiguous(): 

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

961 

962 return loadRes 

963 

964 @abc.abstractmethod 

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

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

967 

968 Parameters 

969 ---------- 

970 ctrCoord : `lsst.geom.SpherePoint` 

971 ICRS center of search region. 

972 radius : `lsst.geom.Angle` 

973 Radius of search region. 

974 filterName : `str` (optional) 

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

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

977 (which are not yet implemented). 

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

979 Epoch to which to correct proper motion and parallax, 

980 or None to not apply such corrections. 

981 centroids : `bool` (optional) 

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

983 these fields to exist. 

984 

985 Returns 

986 ------- 

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

988 A Struct containing the following fields: 

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

990 A catalog of reference objects with the standard 

991 schema, as documented in the main doc string for 

992 `LoadReferenceObjects`. 

993 The catalog is guaranteed to be contiguous. 

994 fluxField : `str` 

995 Name of flux field for specified `filterName`. 

996 

997 Notes 

998 ----- 

999 Note that subclasses are responsible for performing the proper motion 

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

1001 the catalog. 

1002 """ 

1003 return 

1004 

1005 @staticmethod 

1006 def _trimToBBox(refCat, bbox, wcs): 

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

1008 centroid and hasCentroid fields. 

1009 

1010 Parameters 

1011 ---------- 

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

1013 A catalog of objects. The schema must include fields 

1014 "coord", "centroid" and "hasCentroid". 

1015 The "coord" field is read. 

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

1017 bbox : `lsst.geom.Box2D` 

1018 Pixel region 

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

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

1021 

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

1023 """ 

1024 afwTable.updateRefCentroids(wcs, refCat) 

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

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

1027 for star in refCat: 

1028 point = star.get(centroidKey) 

1029 if bbox.contains(point): 

1030 retStarCat.append(star) 

1031 return retStarCat 

1032 

1033 def _addFluxAliases(self, schema): 

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

1035 

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

1037 camFlux: <defaultFilter>_flux 

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

1039 

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

1041 <camFilter>_camFlux: <refFilter>_flux 

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

1043 

1044 Parameters 

1045 ---------- 

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

1047 Schema for reference catalog. 

1048 

1049 Throws 

1050 ------ 

1051 RuntimeError 

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

1053 """ 

1054 aliasMap = schema.getAliasMap() 

1055 

1056 def addAliasesForOneFilter(filterName, refFilterName): 

1057 """Add aliases for a single filter 

1058 

1059 Parameters 

1060 ---------- 

1061 filterName : `str` (optional) 

1062 Camera filter name. The resulting alias name is 

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

1064 is `None` or `""`. 

1065 refFilterName : `str` 

1066 Reference catalog filter name; the field 

1067 <refFilterName>_flux must exist. 

1068 """ 

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

1070 refFluxName = refFilterName + "_flux" 

1071 if refFluxName not in schema: 

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

1073 aliasMap.set(camFluxName, refFluxName) 

1074 refFluxErrName = refFluxName + "Err" 

1075 if refFluxErrName in schema: 

1076 camFluxErrName = camFluxName + "Err" 

1077 aliasMap.set(camFluxErrName, refFluxErrName) 

1078 

1079 if self.config.defaultFilter: 

1080 addAliasesForOneFilter(None, self.config.defaultFilter) 

1081 

1082 for filterName, refFilterName in self.config.filterMap.items(): 

1083 addAliasesForOneFilter(filterName, refFilterName) 

1084 

1085 @staticmethod 

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

1087 addIsPhotometric=False, addIsResolved=False, 

1088 addIsVariable=False, coordErrDim=2, 

1089 addProperMotion=False, properMotionErrDim=2, 

1090 addParallax=False): 

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

1092 

1093 Parameters 

1094 ---------- 

1095 filterNameList : `list` of `str` 

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

1097 addIsPhotometric : `bool` 

1098 If True then add field "photometric". 

1099 addIsResolved : `bool` 

1100 If True then add field "resolved". 

1101 addIsVariable : `bool` 

1102 If True then add field "variable". 

1103 coordErrDim : `int` 

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

1105 

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

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

1108 addProperMotion : `bool` 

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

1110 properMotionErrDim : `int` 

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

1112 ignored if addProperMotion false: 

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

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

1115 addParallax : `bool` 

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

1117 and "parallax_flag". 

1118 

1119 Returns 

1120 ------- 

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

1122 Schema for reference catalog, an 

1123 `lsst.afw.table.SimpleCatalog`. 

1124 

1125 Notes 

1126 ----- 

1127 Reference catalogs support additional covariances, such as 

1128 covariance between RA and proper motion in declination, 

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

1130 calling this method. 

1131 """ 

1132 schema = afwTable.SimpleTable.makeMinimalSchema() 

1133 if addCentroid: 

1134 afwTable.Point2DKey.addFields( 

1135 schema, 

1136 "centroid", 

1137 "centroid on an exposure, if relevant", 

1138 "pixel", 

1139 ) 

1140 schema.addField( 

1141 field="hasCentroid", 

1142 type="Flag", 

1143 doc="is position known?", 

1144 ) 

1145 for filterName in filterNameList: 

1146 schema.addField( 

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

1148 type=numpy.float64, 

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

1150 units="nJy", 

1151 ) 

1152 for filterName in filterNameList: 

1153 schema.addField( 

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

1155 type=numpy.float64, 

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

1157 units="nJy", 

1158 ) 

1159 if addIsPhotometric: 

1160 schema.addField( 

1161 field="photometric", 

1162 type="Flag", 

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

1164 ) 

1165 if addIsResolved: 

1166 schema.addField( 

1167 field="resolved", 

1168 type="Flag", 

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

1170 ) 

1171 if addIsVariable: 

1172 schema.addField( 

1173 field="variable", 

1174 type="Flag", 

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

1176 ) 

1177 if coordErrDim not in (0, 2, 3): 

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

1179 if coordErrDim > 0: 

1180 afwTable.CovarianceMatrix2fKey.addFields( 

1181 schema=schema, 

1182 prefix="coord", 

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

1184 units=["rad", "rad"], 

1185 diagonalOnly=(coordErrDim == 2), 

1186 ) 

1187 

1188 if addProperMotion or addParallax: 

1189 schema.addField( 

1190 field="epoch", 

1191 type=numpy.float64, 

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

1193 units="day", 

1194 ) 

1195 

1196 if addProperMotion: 

1197 schema.addField( 

1198 field="pm_ra", 

1199 type="Angle", 

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

1201 units="rad/year", 

1202 ) 

1203 schema.addField( 

1204 field="pm_dec", 

1205 type="Angle", 

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

1207 units="rad/year", 

1208 ) 

1209 if properMotionErrDim not in (0, 2, 3): 

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

1211 if properMotionErrDim > 0: 

1212 afwTable.CovarianceMatrix2fKey.addFields( 

1213 schema=schema, 

1214 prefix="pm", 

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

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

1217 diagonalOnly=(properMotionErrDim == 2), 

1218 ) 

1219 schema.addField( 

1220 field="pm_flag", 

1221 type="Flag", 

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

1223 ) 

1224 

1225 if addParallax: 

1226 schema.addField( 

1227 field="parallax", 

1228 type="Angle", 

1229 doc="parallax", 

1230 units="rad", 

1231 ) 

1232 schema.addField( 

1233 field="parallaxErr", 

1234 type="Angle", 

1235 doc="uncertainty in parallax", 

1236 units="rad", 

1237 ) 

1238 schema.addField( 

1239 field="parallax_flag", 

1240 type="Flag", 

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

1242 ) 

1243 return schema 

1244 

1245 def _calculateCircle(self, bbox, wcs): 

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

1247 

1248 Parameters 

1249 ---------- 

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

1251 Pixel bounding box. 

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

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

1254 

1255 Returns 

1256 ------- 

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

1258 A Struct containing: 

1259 

1260 - coord : `lsst.geom.SpherePoint` 

1261 ICRS center of the search region. 

1262 - radius : `lsst.geom.Angle` 

1263 Radius of the search region. 

1264 - bbox : `lsst.geom.Box2D` 

1265 Bounding box used to compute the circle. 

1266 """ 

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

1268 bbox.grow(self.config.pixelMargin) 

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

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

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

1272 

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

1274 """Return metadata about the load. 

1275 

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

1277 reconstituting a normalised match list. 

1278 

1279 Parameters 

1280 ---------- 

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

1282 Pixel bounding box. 

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

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

1285 filterName : `str` 

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

1287 filter. 

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

1289 Calibration, or `None` if unknown. 

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

1291 Epoch to which to correct proper motion and parallax, 

1292 or None to not apply such corrections. 

1293 

1294 Returns 

1295 ------- 

1296 metadata : lsst.daf.base.PropertyList 

1297 Metadata about the load. 

1298 """ 

1299 circle = self._calculateCircle(bbox, wcs) 

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

1301 

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

1303 """Return metadata about the load. 

1304 

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

1306 reconstituting a normalised match list. 

1307 

1308 Parameters 

1309 ---------- 

1310 coord : `lsst.geom.SpherePoint` 

1311 ICRS center of the search region. 

1312 radius : `lsst.geom.Angle` 

1313 Radius of the search region. 

1314 filterName : `str` 

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

1316 filter. 

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

1318 Calibration, or `None` if unknown. 

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

1320 Epoch to which to correct proper motion and parallax, 

1321 or None to not apply such corrections. 

1322 

1323 Returns 

1324 ------- 

1325 metadata : lsst.daf.base.PropertyList 

1326 Metadata about the load 

1327 """ 

1328 md = PropertyList() 

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

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

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

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

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

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

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

1336 return md 

1337 

1338 def joinMatchListWithCatalog(self, matchCat, sourceCat): 

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

1340 objects. 

1341 

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

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

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

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

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

1347 

1348 Parameters 

1349 ---------- 

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

1351 Unperisted packed match list. 

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

1353 as returned by the astrometry tasks. 

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

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

1356 by ID. 

1357 

1358 Returns 

1359 ------- 

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

1361 Match list. 

1362 """ 

1363 return joinMatchListWithCatalogImpl(self, matchCat, sourceCat) 

1364 

1365 def applyProperMotions(self, catalog, epoch): 

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

1367 

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

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

1370 modifying the catalong in place. 

1371 

1372 Parameters 

1373 ---------- 

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

1375 Catalog of positions, containing: 

1376 

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

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

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

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

1381 East positive) 

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

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

1384 North positive) 

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

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

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

1388 Epoch to which to correct proper motion and parallax, 

1389 or None to not apply such corrections. 

1390 """ 

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

1392 if self.config.requireProperMotion: 

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

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

1395 return 

1396 applyProperMotionsImpl(self.log, catalog, epoch) 

1397 

1398 

1399def joinMatchListWithCatalogImpl(refObjLoader, matchCat, sourceCat): 

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

1401 objects. 

1402 

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

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

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

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

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

1408 

1409 Parameters 

1410 ---------- 

1411 refObjLoader 

1412 Reference object loader to use in getting reference objects 

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

1414 Unperisted packed match list. 

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

1416 as returned by the astrometry tasks. 

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

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

1419 by ID. 

1420 

1421 Returns 

1422 ------- 

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

1424 Match list. 

1425 """ 

1426 matchmeta = matchCat.table.getMetadata() 

1427 version = matchmeta.getInt('SMATCHV') 

1428 if version != 1: 

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

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

1431 try: 

1432 epoch = matchmeta.getDouble('EPOCH') 

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

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

1435 if 'RADIUS' in matchmeta: 

1436 # This is a circle style metadata, call loadSkyCircle 

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

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

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

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

1441 elif "INNER_UPPER_LEFT_RA" in matchmeta: 

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

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

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

1445 # by the refObjLoader 

1446 box = [] 

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

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

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

1450 lsst.geom.degrees).getVector() 

1451 box.append(coord) 

1452 outerBox = sphgeom.ConvexPolygon(box) 

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

1454 

1455 refCat.sort() 

1456 sourceCat.sort() 

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

1458 

1459 

1460def applyProperMotionsImpl(log, catalog, epoch): 

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

1462 

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

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

1465 modifying the catalong in place. 

1466 

1467 Parameters 

1468 ---------- 

1469 log : `lsst.log.log` 

1470 log object to write to 

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

1472 Catalog of positions, containing: 

1473 

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

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

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

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

1478 East positive) 

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

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

1481 North positive) 

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

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

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

1485 Epoch to which to correct proper motion and parallax, 

1486 or None to not apply such corrections. 

1487 """ 

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

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

1490 return 

1491 if not catalog.isContiguous(): 

1492 raise RuntimeError("Catalog must be contiguous") 

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

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

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

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

1497 coordKey = catalog.table.getCoordKey() 

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

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

1500 pmRaRad = catalog["pm_ra"] 

1501 pmDecRad = catalog["pm_dec"] 

1502 offsetsRaRad = pmRaRad*timeDiffsYears 

1503 offsetsDecRad = pmDecRad*timeDiffsYears 

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

1505 # due to proper motion, and apply the offset 

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

1507 # a reasonable scale for typical values of proper motion 

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

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

1510 # needlessly large errors for short duration 

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

1512 offsetAmountsRad = numpy.hypot(offsetsRaRad, offsetsDecRad) 

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

1514 record.set(coordKey, 

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

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

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

1518 if "coord_raErr" in catalog.schema: 

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

1520 catalog["pm_raErr"]*timeDiffsYears) 

1521 if "coord_decErr" in catalog.schema: 

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

1523 catalog["pm_decErr"]*timeDiffsYears)