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

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

172

173

174

175

176

177

178

179

180

181

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

222

223

224

225

226

227

228

229

230

231

232

233

234

235

236

237

238

239

240

241

242

243

244

245

246

247

# LSST Data Management System 

# Copyright 2016 AURA/LSST. 

# 

# This product includes software developed by the 

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

# 

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

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

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

# (at your option) any later version. 

# 

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

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

# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 

# GNU General Public License for more details. 

# 

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

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

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

"""Analytic single-visit photometric error model. 

""" 

 

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

 

import astropy.units as u 

import numpy as np 

from scipy.optimize import curve_fit 

 

from lsst.verify import Blob, Datum 

 

 

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

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

 

The model is described in the LSST Overview paper: 

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

 

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

 

.. math: 

 

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

 

Where random uncertainty is (Eq. 5): 

 

.. math:: 

 

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

 

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

 

Parameters 

---------- 

mag : `astropy.unit.Quantity` 

Source magnitude. 

sigmaSq : `astropy.unit.Quantity` 

Systematic photometric error (magnitude). 

gamma : `astropy.unit.Quantity` 

Proxy for sky brightness and readout noise (dimensionless). 

m5 : `astropy.unit.Quantity` 

5-sigma depth (magnitude). 

 

Returns 

------- 

sigma : `astropy.units.Quantity` 

Photometric error for a single visit. 

""" 

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

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

sigmaSq = sigmaSys**2 + sigmaRandSq 

return np.sqrt(sigmaSq) 

 

 

def fitPhotErrModel(mag, mag_err): 

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

 

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

 

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

 

Parameters 

---------- 

mag : `astropy.units.Quantity` 

Magnitude. 

mag_err : `astropy.units.Quantity` 

Magnitude uncertainty or variation. 

 

Returns 

------- 

params : `dict` 

Fitted model. Fields are: 

 

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

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

`astropy.units.Quantity`). 

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

 

See also 

-------- 

photErrModel 

""" 

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

# mag = mag * u.mag 

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

# mag_err = mag_err * u.mag 

 

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

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

# 24.35 * u.mag] # m5 

 

# fit_params, fit_param_covariance = curve_fit( 

# photErrModel, mag, mag_err, p0=p0) 

 

# params = { 

# 'sigmaSys': fit_params[0], 

# 'gamma': fit_params[1], 

# 'm5': fit_params[2], 

# } 

if isinstance(mag, u.Quantity): 

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

if isinstance(mag_err, u.Quantity): 

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

 

p0 = [0.01, # sigmaSys (mag) 

0.039, # gamma ('') 

24.35] # m5 (mag) 

 

try: 

fit_params, fit_param_covariance = curve_fit( 

photErrModel, mag, mag_err, p0=p0) 

sigmaSys, gamma, m5 = fit_params 

except (RuntimeError, ValueError) as e: 

print("fitPhotErrorModel fitting failed with") 

print(e) 

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

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

 

params = { 

'sigmaSys': sigmaSys * u.mag, 

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

'm5': m5 * u.mag, 

} 

return params 

 

 

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

matchRef=500): 

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

 

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

(Eq 4, 5): 

 

.. math:: 

 

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

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

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

 

Parameters 

---------- 

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

A dataset containing matched statistics for stars across multiple 

visits. 

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

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

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

Median reference astrometric scatter (millimagnitudes by default). 

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

Should match at least matchRef stars. 

 

Returns 

---------- 

blob : `lsst.verify.Blob` 

Blob with datums: 

 

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

- ``sigmaSys``: Systematic error floor. 

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

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

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

 

Notes 

----- 

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

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

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

""" 

blob = Blob('PhotometricErrorModel') 

 

# FIXME add a description field to blobs? 

# _doc['doc'] \ 

# = "Photometric uncertainty model from " \ 

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

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

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

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

 

if not isinstance(medianRef, u.Quantity): 

medianRef = medianRef * u.mmag 

if not isinstance(brightSnr, u.Quantity): 

brightSnr = brightSnr * u.Unit('') 

_compute(blob, 

matchedMultiVisitDataset['snr'].quantity, 

matchedMultiVisitDataset['mag'].quantity, 

matchedMultiVisitDataset['magerr'].quantity, 

matchedMultiVisitDataset['magrms'].quantity, 

matchedMultiVisitDataset['dist'].quantity, 

len(matchedMultiVisitDataset.goodMatches), 

brightSnr, 

medianRef, 

matchRef) 

return blob 

 

 

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

brightSnr, medianRef, matchRef): 

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

label='Bright SNR', 

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

'model') 

 

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

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

label='RMS', 

description='RMS photometric scatter for good stars') 

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

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

 

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

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

label='sigma(sys)', 

description='Systematic error floor') 

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

label='gamma', 

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

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

label='m5', 

description='5-sigma depth') 

 

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

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

'reference : {1:.3f}' 

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

if nMatch < matchRef: 

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

'(should be > %d)' 

print(msg.format(nMatch, matchRef))