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# LSST Data Management System 

2# Copyright 2016 AURA/LSST. 

3# 

4# This product includes software developed by the 

5# LSST Project (http://www.lsst.org/). 

6# 

7# This program is free software: you can redistribute it and/or modify 

8# it under the terms of the GNU General Public License as published by 

9# the Free Software Foundation, either version 3 of the License, or 

10# (at your option) any later version. 

11# 

12# This program is distributed in the hope that it will be useful, 

13# but WITHOUT ANY WARRANTY; without even the implied warranty of 

14# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 

15# GNU General Public License for more details. 

16# 

17# You should have received a copy of the LSST License Statement and 

18# the GNU General Public License along with this program. If not, 

19# see <https://www.lsstcorp.org/LegalNotices/>. 

20"""Analytic single-visit photometric error model. 

21""" 

22 

23__all__ = ['photErrModel', 'fitPhotErrModel', 'build_photometric_error_model'] 

24 

25import astropy.units as u 

26import numpy as np 

27from scipy.optimize import curve_fit 

28 

29from lsst.verify import Blob, Datum 

30 

31 

32def photErrModel(mag, sigmaSys, gamma, m5, **kwargs): 

33 r"""Model of photometric error for a single visit. 

34 

35 The model is described in the LSST Overview paper: 

36 http://arxiv.org/abs/0805.2366v4. 

37 

38 Photometric error of a single visit (Eq. 4): 

39 

40 .. math: 

41 

42 \sigma_1^2 = \sigma_\mathrm{sys}^2 + \sigma_\mathrm{rand}^2 

43 

44 Where random uncertainty is (Eq. 5): 

45 

46 .. math:: 

47 

48 \sigma_\mathrm{rand}^2 = (0.04 - \gamma) * x + \gamma * x^2~[\mathrm{mag}^2] 

49 

50 and $x = 10^(0.4*(m-m_5))$. 

51 

52 Parameters 

53 ---------- 

54 mag : `astropy.unit.Quantity` 

55 Source magnitude. 

56 sigmaSq : `astropy.unit.Quantity` 

57 Systematic photometric error (magnitude). 

58 gamma : `astropy.unit.Quantity` 

59 Proxy for sky brightness and readout noise (dimensionless). 

60 m5 : `astropy.unit.Quantity` 

61 5-sigma depth (magnitude). 

62 

63 Returns 

64 ------- 

65 sigma : `astropy.units.Quantity` 

66 Photometric error for a single visit. 

67 """ 

68 x = 10**(0.4*(mag - m5)) 

69 sigmaRandSq = (0.04 - gamma) * x + gamma * x**2 

70 sigmaSq = sigmaSys**2 + sigmaRandSq 

71 return np.sqrt(sigmaSq) 

72 

73 

74def fitPhotErrModel(mag, mag_err): 

75 """Fit photometric error model from the LSST Overview paper: 

76 

77 http://arxiv.org/abs/0805.2366v4 

78 

79 The fit is performed with `scipy.optimize.curvefit`. 

80 

81 Parameters 

82 ---------- 

83 mag : `astropy.units.Quantity` 

84 Magnitude. 

85 mag_err : `astropy.units.Quantity` 

86 Magnitude uncertainty or variation. 

87 

88 Returns 

89 ------- 

90 params : `dict` 

91 Fitted model. Fields are: 

92 

93 - `sigmaSys`: systematic error (magnitude, `astropy.units.Quantity`). 

94 - `gamma`: Proxy for sky brightness and readout noise (dimensionless, 

95 `astropy.units.Quantity`). 

96 - `m5`: 5-sigma limiting depth (magnitude, `astropy.units.Quantity`). 

97 

98 See also 

99 -------- 

100 photErrModel 

101 """ 

102 # if not isinstance(mag, u.Quantity): 

103 # mag = mag * u.mag 

104 # if not isinstance(mag_err, u.Quantity): 

105 # mag_err = mag_err * u.mag 

106 

107 # p0 = [0.01 * u.mag, # sigmaSys 

108 # 0.039 * u.Unit(''), # gamma 

109 # 24.35 * u.mag] # m5 

110 

111 # fit_params, fit_param_covariance = curve_fit( 

112 # photErrModel, mag, mag_err, p0=p0) 

113 

114 # params = { 

115 # 'sigmaSys': fit_params[0], 

116 # 'gamma': fit_params[1], 

117 # 'm5': fit_params[2], 

118 # } 

119 if isinstance(mag, u.Quantity): 

120 mag = mag.to(u.mag).value 

121 if isinstance(mag_err, u.Quantity): 

122 mag_err = mag_err.to(u.mag).value 

123 

124 p0 = [0.01, # sigmaSys (mag) 

125 0.039, # gamma ('') 

126 24.35] # m5 (mag) 

127 

128 try: 

129 fit_params, fit_param_covariance = curve_fit( 

130 photErrModel, mag, mag_err, p0=p0) 

131 sigmaSys, gamma, m5 = fit_params 

132 except (RuntimeError, ValueError) as e: 

133 print("fitPhotErrorModel fitting failed with") 

134 print(e) 

135 print("sigmaSys, gamma, m5 are being set to NaN.") 

136 sigmaSys, gamma, m5 = np.nan, np.nan, np.nan 

137 

138 params = { 

139 'sigmaSys': sigmaSys * u.mag, 

140 'gamma': gamma * u.Unit(''), 

141 'm5': m5 * u.mag, 

142 } 

143 return params 

144 

145 

146def build_photometric_error_model(matchedMultiVisitDataset, brightSnr=100, medianRef=100, 

147 matchRef=500): 

148 r"""Returns a serializable analytic photometry error model for a single visit. 

149 

150 This model is originally presented in http://arxiv.org/abs/0805.2366v4 

151 (Eq 4, 5): 

152 

153 .. math:: 

154 

155 \sigma_1^2 &= \sigma_\mathrm{sys}^2 + \sigma_\mathrm{rand}^2 \\ 

156 x &= 10^{0.4(m-m_5)} \\ 

157 \sigma_\mathrm{rand}^2 &= (0.04 - \gamma) x + \gamma x^2~[\mathrm{mag}^2] 

158 

159 Parameters 

160 ---------- 

161 matchedMultiVisitDataset : `lsst.valididate.drp.matchreduce.MatchedMultiVisitDataset` 

162 A dataset containing matched statistics for stars across multiple 

163 visits. 

164 brightSnr : `float` or `astropy.unit.Quantity`, optional 

165 Minimum SNR for a star to be considered "bright." 

166 medianRef : `float` or `astropy.unit.Quantity`, optional 

167 Median reference astrometric scatter (millimagnitudes by default). 

168 matchRef : `int` or `astropy.unit.Quantity`, optional 

169 Should match at least matchRef stars. 

170 

171 Returns 

172 ---------- 

173 blob : `lsst.verify.Blob` 

174 Blob with datums: 

175 

176 - ``brightSnr``: Threshold in SNR for bright sources used in this model. 

177 - ``sigmaSys``: Systematic error floor. 

178 - ``gamma``: Proxy for sky brightness and read noise. 

179 - ``m5``: 5-sigma photometric depth (magnitudes). 

180 - ``photRms``: RMS photometric scatter for 'good' stars (millimagnitudes). 

181 

182 Notes 

183 ----- 

184 The scatter and match defaults are appropriate to SDSS are stored here. 

185 For SDSS, stars with mag < 19.5 should be completely well measured. 

186 This limit is a band-dependent statement most appropriate to r. 

187 """ 

188 blob = Blob('PhotometricErrorModel') 

189 

190 # FIXME add a description field to blobs? 

191 # _doc['doc'] \ 

192 # = "Photometric uncertainty model from " \ 

193 # "http://arxiv.org/abs/0805.2366v4 (Eq 4, 5): " \ 

194 # "sigma_1^2 = sigma_sys^2 + sigma_rand^2, " \ 

195 # "sigma_rand^2 = (0.04 - gamma) * x + gamma * x^2 [mag^2] " \ 

196 # "where x = 10**(0.4*(m-m_5))" 

197 

198 if not isinstance(medianRef, u.Quantity): 

199 medianRef = medianRef * u.mmag 

200 if not isinstance(brightSnr, u.Quantity): 

201 brightSnr = brightSnr * u.Unit('') 

202 _compute(blob, 

203 matchedMultiVisitDataset['snr'].quantity, 

204 matchedMultiVisitDataset['mag'].quantity, 

205 matchedMultiVisitDataset['magerr'].quantity, 

206 matchedMultiVisitDataset['magrms'].quantity, 

207 matchedMultiVisitDataset['dist'].quantity, 

208 len(matchedMultiVisitDataset.goodMatches), 

209 brightSnr, 

210 medianRef, 

211 matchRef) 

212 return blob 

213 

214 

215def _compute(blob, snr, mag, magErr, magRms, dist, nMatch, 

216 brightSnr, medianRef, matchRef): 

217 blob['brightSnr'] = Datum(quantity=brightSnr, 

218 label='Bright SNR', 

219 description='Threshold in SNR for bright sources used in this ' 

220 'model') 

221 

222 bright = np.where(snr > blob['brightSnr'].quantity) 

223 blob['photScatter'] = Datum(quantity=np.median(magRms[bright]), 

224 label='RMS', 

225 description='RMS photometric scatter for good stars') 

226 print('Photometric scatter (median) - SNR > {0:.1f} : {1:.1f}'.format( 

227 blob['brightSnr'].quantity, blob['photScatter'].quantity.to(u.mmag))) 

228 

229 fit_params = fitPhotErrModel(mag[bright], magErr[bright]) 

230 blob['sigmaSys'] = Datum(quantity=fit_params['sigmaSys'], 

231 label='sigma(sys)', 

232 description='Systematic error floor') 

233 blob['gamma'] = Datum(quantity=fit_params['gamma'], 

234 label='gamma', 

235 description='Proxy for sky brightness and read noise') 

236 blob['m5'] = Datum(quantity=fit_params['m5'], 

237 label='m5', 

238 description='5-sigma depth') 

239 

240 if blob['photScatter'].quantity > medianRef: 

241 msg = 'Median photometric scatter {0:.3f} is larger than ' \ 

242 'reference : {1:.3f}' 

243 print(msg.format(blob['photScatter'].quantity.value, medianRef)) 

244 if nMatch < matchRef: 

245 msg = 'Number of matched sources {0:d} is too small ' \ 

246 '(should be > {0:d})' 

247 print(msg.format(nMatch, matchRef))