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__ = ("MeasureApCorrConfig", "MeasureApCorrTask") 

25 

26import numpy 

27 

28import lsst.pex.config 

29from lsst.afw.image import ApCorrMap 

30from lsst.afw.math import ChebyshevBoundedField, ChebyshevBoundedFieldConfig 

31from lsst.pipe.base import Task, Struct 

32from lsst.meas.base.apCorrRegistry import getApCorrNameSet 

33 

34from .sourceSelector import sourceSelectorRegistry 

35 

36 

37class FluxKeys: 

38 """A collection of keys for a given flux measurement algorithm 

39 """ 

40 __slots__ = ("flux", "err", "flag", "used") # prevent accidentally adding fields 

41 

42 def __init__(self, name, schema): 

43 """Construct a FluxKeys 

44 

45 @parma[in] name name of flux measurement algorithm, e.g. "base_PsfFlux" 

46 @param[in,out] schema catalog schema containing the flux field 

47 read: {name}_instFlux, {name}_instFluxErr, {name}_flag 

48 added: apcorr_{name}_used 

49 """ 

50 self.flux = schema.find(name + "_instFlux").key 

51 self.err = schema.find(name + "_instFluxErr").key 

52 self.flag = schema.find(name + "_flag").key 

53 self.used = schema.addField("apcorr_" + name + "_used", type="Flag", 

54 doc="set if source was used in measuring aperture correction") 

55 

56# The following block adds links to these tasks from the Task Documentation page. 

57## @addtogroup LSST_task_documentation 

58## @{ 

59## @page measureApCorrTask 

60## @ref MeasureApCorrTask "MeasureApCorrTask" 

61## Task to measure aperture correction 

62## @} 

63 

64 

65class MeasureApCorrConfig(lsst.pex.config.Config): 

66 """!Configuration for MeasureApCorrTask 

67 """ 

68 refFluxName = lsst.pex.config.Field( 

69 doc="Field name prefix for the flux other measurements should be aperture corrected to match", 

70 dtype=str, 

71 default="slot_CalibFlux", 

72 ) 

73 sourceSelector = sourceSelectorRegistry.makeField( 

74 doc="Selector that sets the stars that aperture corrections will be measured from.", 

75 default="flagged", 

76 ) 

77 minDegreesOfFreedom = lsst.pex.config.RangeField( 

78 doc="Minimum number of degrees of freedom (# of valid data points - # of parameters);" 

79 " if this is exceeded, the order of the fit is decreased (in both dimensions), and" 

80 " if we can't decrease it enough, we'll raise ValueError.", 

81 dtype=int, 

82 default=1, 

83 min=1, 

84 ) 

85 fitConfig = lsst.pex.config.ConfigField( 

86 doc="Configuration used in fitting the aperture correction fields", 

87 dtype=ChebyshevBoundedFieldConfig, 

88 ) 

89 numIter = lsst.pex.config.Field( 

90 doc="Number of iterations for sigma clipping", 

91 dtype=int, 

92 default=4, 

93 ) 

94 numSigmaClip = lsst.pex.config.Field( 

95 doc="Number of standard devisations to clip at", 

96 dtype=float, 

97 default=3.0, 

98 ) 

99 allowFailure = lsst.pex.config.ListField( 

100 doc="Allow these measurement algorithms to fail without an exception", 

101 dtype=str, 

102 default=[], 

103 ) 

104 

105 def validate(self): 

106 lsst.pex.config.Config.validate(self) 

107 if self.sourceSelector.target.usesMatches: 

108 raise lsst.pex.config.FieldValidationError( 

109 MeasureApCorrConfig.sourceSelector, 

110 self, 

111 "Star selectors that require matches are not permitted") 

112 

113 

114class MeasureApCorrTask(Task): 

115 r"""!Task to measure aperture correction 

116 

117 @section measAlg_MeasureApCorrTask_Contents Contents 

118 

119 - @ref measAlg_MeasureApCorrTask_Purpose 

120 - @ref measAlg_MeasureApCorrTask_Config 

121 - @ref measAlg_MeasureApCorrTask_Debug 

122 

123 @section measAlg_MeasureApCorrTask_Purpose Description 

124 

125 @copybrief MeasureApCorrTask 

126 

127 This task measures aperture correction for the flux fields returned by 

128 lsst.meas.base.getApCorrNameSet() 

129 

130 The main method is @ref MeasureApCorrTask.run "run". 

131 

132 @section measAlg_MeasureApCorrTask_Config Configuration parameters 

133 

134 See @ref MeasureApCorrConfig 

135 

136 @section measAlg_MeasureApCorrTask_Debug Debug variables 

137 

138 The @link lsst.pipe.base.cmdLineTask.CmdLineTask command line task@endlink interface supports a flag 

139 `--debug` to import `debug.py` from your `$PYTHONPATH`; see @ref baseDebug for more about `debug.py`. 

140 

141 MeasureApCorrTask has a debug dictionary containing a single boolean key: 

142 <dl> 

143 <dt>display 

144 <dd>If True: will show plots as aperture corrections are fitted 

145 </dl> 

146 

147 For example, put something like: 

148 @code{.py} 

149 import lsstDebug 

150 def DebugInfo(name): 

151 di = lsstDebug.getInfo(name) # N.b. lsstDebug.Info(name) would call us recursively 

152 if name == "lsst.meas.algorithms.measureApCorr": 

153 di.display = dict( 

154 unsubtracted = 1, 

155 subtracted = 2, 

156 background = 3, 

157 ) 

158 

159 return di 

160 

161 lsstDebug.Info = DebugInfo 

162 @endcode 

163 into your `debug.py` file and run your command-line task with the `--debug` flag (or `import debug`). 

164 """ 

165 ConfigClass = MeasureApCorrConfig 

166 _DefaultName = "measureApCorr" 

167 

168 def __init__(self, schema, **kwds): 

169 """!Construct a MeasureApCorrTask 

170 

171 For every name in lsst.meas.base.getApCorrNameSet(): 

172 - If the corresponding flux fields exist in the schema: 

173 - Add a new field apcorr_{name}_used 

174 - Add an entry to the self.toCorrect dict 

175 - Otherwise silently skip the name 

176 """ 

177 Task.__init__(self, **kwds) 

178 self.refFluxKeys = FluxKeys(self.config.refFluxName, schema) 

179 self.toCorrect = {} # dict of flux field name prefix: FluxKeys instance 

180 for name in getApCorrNameSet(): 

181 try: 

182 self.toCorrect[name] = FluxKeys(name, schema) 

183 except KeyError: 

184 # if a field in the registry is missing, just ignore it. 

185 pass 

186 self.makeSubtask("sourceSelector") 

187 

188 def run(self, exposure, catalog): 

189 """!Measure aperture correction 

190 

191 @param[in] exposure Exposure aperture corrections are being measured 

192 on. The bounding box is retrieved from it, and 

193 it is passed to the sourceSelector. 

194 The output aperture correction map is *not* 

195 added to the exposure; this is left to the 

196 caller. 

197 

198 @param[in] catalog SourceCatalog containing measurements to be used 

199 to compute aperturecorrections. 

200 

201 @return an lsst.pipe.base.Struct containing: 

202 - apCorrMap: an aperture correction map (lsst.afw.image.ApCorrMap) that contains two entries 

203 for each flux field: 

204 - flux field (e.g. base_PsfFlux_instFlux): 2d model 

205 - flux sigma field (e.g. base_PsfFlux_instFluxErr): 2d model of error 

206 """ 

207 bbox = exposure.getBBox() 

208 import lsstDebug 

209 display = lsstDebug.Info(__name__).display 

210 doPause = lsstDebug.Info(__name__).doPause 

211 

212 self.log.info("Measuring aperture corrections for %d flux fields" % (len(self.toCorrect),)) 

213 # First, create a subset of the catalog that contains only selected stars 

214 # with non-flagged reference fluxes. 

215 subset1 = [record for record in self.sourceSelector.run(catalog, exposure=exposure).sourceCat 

216 if (not record.get(self.refFluxKeys.flag) 

217 and numpy.isfinite(record.get(self.refFluxKeys.flux)))] 

218 

219 apCorrMap = ApCorrMap() 

220 

221 # Outer loop over the fields we want to correct 

222 for name, keys in self.toCorrect.items(): 

223 fluxName = name + "_instFlux" 

224 fluxErrName = name + "_instFluxErr" 

225 

226 # Create a more restricted subset with only the objects where the to-be-correct flux 

227 # is not flagged. 

228 fluxes = numpy.fromiter((record.get(keys.flux) for record in subset1), float) 

229 with numpy.errstate(invalid="ignore"): # suppress NAN warnings 

230 isGood = numpy.logical_and.reduce([ 

231 numpy.fromiter((not record.get(keys.flag) for record in subset1), bool), 

232 numpy.isfinite(fluxes), 

233 fluxes > 0.0, 

234 ]) 

235 subset2 = [record for record, good in zip(subset1, isGood) if good] 

236 

237 # Check that we have enough data points that we have at least the minimum of degrees of 

238 # freedom specified in the config. 

239 if len(subset2) - 1 < self.config.minDegreesOfFreedom: 

240 if name in self.config.allowFailure: 

241 self.log.warn("Unable to measure aperture correction for '%s': " 

242 "only %d sources, but require at least %d." % 

243 (name, len(subset2), self.config.minDegreesOfFreedom+1)) 

244 continue 

245 raise RuntimeError("Unable to measure aperture correction for required algorithm '%s': " 

246 "only %d sources, but require at least %d." % 

247 (name, len(subset2), self.config.minDegreesOfFreedom+1)) 

248 

249 # If we don't have enough data points to constrain the fit, reduce the order until we do 

250 ctrl = self.config.fitConfig.makeControl() 

251 while len(subset2) - ctrl.computeSize() < self.config.minDegreesOfFreedom: 

252 if ctrl.orderX > 0: 252 ↛ 254line 252 didn't jump to line 254, because the condition on line 252 was never false

253 ctrl.orderX -= 1 

254 if ctrl.orderY > 0: 254 ↛ 251line 254 didn't jump to line 251, because the condition on line 254 was never false

255 ctrl.orderY -= 1 

256 

257 # Fill numpy arrays with positions and the ratio of the reference flux to the to-correct flux 

258 x = numpy.zeros(len(subset2), dtype=float) 

259 y = numpy.zeros(len(subset2), dtype=float) 

260 apCorrData = numpy.zeros(len(subset2), dtype=float) 

261 indices = numpy.arange(len(subset2), dtype=int) 

262 for n, record in enumerate(subset2): 

263 x[n] = record.getX() 

264 y[n] = record.getY() 

265 apCorrData[n] = record.get(self.refFluxKeys.flux)/record.get(keys.flux) 

266 

267 for _i in range(self.config.numIter): 

268 

269 # Do the fit, save it in the output map 

270 apCorrField = ChebyshevBoundedField.fit(bbox, x, y, apCorrData, ctrl) 

271 

272 if display: 272 ↛ 273line 272 didn't jump to line 273, because the condition on line 272 was never true

273 plotApCorr(bbox, x, y, apCorrData, apCorrField, "%s, iteration %d" % (name, _i), doPause) 

274 

275 # Compute errors empirically, using the RMS difference between the true reference flux and the 

276 # corrected to-be-corrected flux. 

277 apCorrDiffs = apCorrField.evaluate(x, y) 

278 apCorrDiffs -= apCorrData 

279 apCorrErr = numpy.mean(apCorrDiffs**2)**0.5 

280 

281 # Clip bad data points 

282 apCorrDiffLim = self.config.numSigmaClip * apCorrErr 

283 with numpy.errstate(invalid="ignore"): # suppress NAN warning 

284 keep = numpy.fabs(apCorrDiffs) <= apCorrDiffLim 

285 x = x[keep] 

286 y = y[keep] 

287 apCorrData = apCorrData[keep] 

288 indices = indices[keep] 

289 

290 # Final fit after clipping 

291 apCorrField = ChebyshevBoundedField.fit(bbox, x, y, apCorrData, ctrl) 

292 

293 self.log.info("Aperture correction for %s: RMS %f from %d" % 

294 (name, numpy.mean((apCorrField.evaluate(x, y) - apCorrData)**2)**0.5, len(indices))) 

295 

296 if display: 296 ↛ 297line 296 didn't jump to line 297, because the condition on line 296 was never true

297 plotApCorr(bbox, x, y, apCorrData, apCorrField, "%s, final" % (name,), doPause) 

298 

299 # Save the result in the output map 

300 # The error is constant spatially (we could imagine being 

301 # more clever, but we're not yet sure if it's worth the effort). 

302 # We save the errors as a 0th-order ChebyshevBoundedField 

303 apCorrMap[fluxName] = apCorrField 

304 apCorrErrCoefficients = numpy.array([[apCorrErr]], dtype=float) 

305 apCorrMap[fluxErrName] = ChebyshevBoundedField(bbox, apCorrErrCoefficients) 

306 

307 # Record which sources were used 

308 for i in indices: 

309 subset2[i].set(keys.used, True) 

310 

311 return Struct( 

312 apCorrMap=apCorrMap, 

313 ) 

314 

315 

316def plotApCorr(bbox, xx, yy, zzMeasure, field, title, doPause): 

317 """Plot aperture correction fit residuals 

318 

319 There are two subplots: residuals against x and y. 

320 

321 Intended for debugging. 

322 

323 @param bbox Bounding box (for bounds) 

324 @param xx x coordinates 

325 @param yy y coordinates 

326 @param zzMeasure Measured value of the aperture correction 

327 @param field Fit aperture correction field 

328 @param title Title for plot 

329 @param doPause Pause to inspect the residuals plot? If False, 

330 there will be a 4 second delay to allow for 

331 inspection of the plot before closing it and 

332 moving on. 

333 """ 

334 import matplotlib.pyplot as plt 

335 

336 zzFit = field.evaluate(xx, yy) 

337 residuals = zzMeasure - zzFit 

338 

339 fig, axes = plt.subplots(2, 1) 

340 

341 axes[0].scatter(xx, residuals, s=3, marker='o', lw=0, alpha=0.7) 

342 axes[1].scatter(yy, residuals, s=3, marker='o', lw=0, alpha=0.7) 

343 for ax in axes: 

344 ax.set_ylabel("ApCorr Fit Residual") 

345 ax.set_ylim(0.9*residuals.min(), 1.1*residuals.max()) 

346 axes[0].set_xlabel("x") 

347 axes[0].set_xlim(bbox.getMinX(), bbox.getMaxX()) 

348 axes[1].set_xlabel("y") 

349 axes[1].set_xlim(bbox.getMinY(), bbox.getMaxY()) 

350 plt.suptitle(title) 

351 

352 if not doPause: 

353 try: 

354 plt.pause(4) 

355 plt.close() 

356 except Exception: 

357 print("%s: plt.pause() failed. Please close plots when done." % __name__) 

358 plt.show() 

359 else: 

360 print("%s: Please close plots when done." % __name__) 

361 plt.show()