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 astrometric accuracy model. 

21""" 

22 

23__all__ = ['astromErrModel', 'fitAstromErrModel', 'build_astrometric_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 astromErrModel(snr, theta=1000, sigmaSys=10, C=1): 

33 """Calculate expected astrometric uncertainty based on SNR. 

34 

35 mas = C*theta/SNR + sigmaSys 

36 

37 Parameters 

38 ---------- 

39 snr : `numpy.ndarray` or `astropy.unit.Quantity` 

40 S/N of photometric measurements (dimensionless). 

41 theta : `float`, `numpy.ndarray` or `astropy.unit.Quantity`, optional 

42 Seeing (default: milliarcsec). 

43 sigmaSys : `astropy.unit.Quantity` 

44 Systematic error floor (default: milliarcsec). 

45 C : `float` 

46 Scaling factor (dimensionless) 

47 

48 Returns 

49 ------- 

50 sigma : `astropy.unit.Quantity` 

51 Expected astrometric uncertainty with the same dimensions as ``snr``. 

52 Units will be those of theta and sigmaSys. 

53 

54 Notes 

55 ----- 

56 ``theta`` and ``sigmaSys`` must be given in the same units. 

57 Typically choices might be any of arcsec, milli-arcsec, or radians. 

58 The default values are reasonable astronominal values in milliarcsec. 

59 But the only thing that matters is that they're the same. 

60 """ 

61 return C*theta/snr + sigmaSys 

62 

63 

64def fitAstromErrModel(snr, dist): 

65 """Fit model of astrometric error from the LSST Overview paper: 

66 

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

68 

69 Parameters 

70 ---------- 

71 snr : `np.ndarray` or `astropy.unit.Quantity` 

72 Signal-to-noise ratio of photometric observations (dimensionless). 

73 dist : `np.ndarray` or `astropy.unit.Quantity` 

74 Scatter in measured positions (default: millarcsec) 

75 

76 Returns 

77 ------- 

78 params : `dict` 

79 Fitted model parameters. Fields are: 

80 

81 - ``C``: Model scale factor (dimensionless). 

82 - ``theta``: Seeing (default: milliarcsec). 

83 - ``sigmaSys``: Systematic astrometric uncertainty 

84 (default: milliarcsec). 

85 """ 

86 # Note that C is fixed to 1. 

87 p0 = [1, # theta 

88 0.01] # sigmaSys 

89 if isinstance(dist, u.Quantity): 

90 dist = dist.to(u.marcsec).value 

91 if isinstance(snr, u.Quantity): 

92 snr = snr.value 

93 fit_params, fit_param_covariance = curve_fit(astromErrModel, snr, dist, 

94 p0=p0) 

95 

96 params = {'C': 1 * u.Unit(''), 

97 'theta': fit_params[0] * u.marcsec, 

98 'sigmaSys': fit_params[1] * u.marcsec} 

99 return params 

100 

101 

102def build_astrometric_error_model(matchedMultiVisitDataset, selection, medianRef=100, matchRef=500): 

103 r"""Serializable model of astrometry errors across multiple visits. 

104 

105 .. math:: 

106 

107 \mathrm{astromRms} = C \theta / \mathrm{SNR} + \sigma_\mathrm{sys} 

108 

109 Parameters 

110 ---------- 

111 matchedMultiVisitDataset : `MatchedMultiVisitDataset` 

112 A dataset containing matched statistics for stars across multiple 

113 visits. 

114 selection : `np.array` of `bool` 

115 The selection of sources to use to build the model. 

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

117 Median reference astrometric scatter (default: milliarcsecond). 

118 matchRef : int, optional 

119 Should match at least matchRef number of stars (dimensionless). 

120 

121 Returns 

122 ------- 

123 blob : `lsst.verify.Blob` 

124 Blob with datums: 

125 

126 - ``C``: Model scaling factor. 

127 - ``theta``: Seeing (milliarcsecond). 

128 - ``sigmaSys``: Systematic error floor (milliarcsecond). 

129 - ``astromRms``: Astrometric scatter (RMS) for good stars (milliarcsecond). 

130 

131 Notes 

132 ----- 

133 The scatter and match defaults appropriate to SDSS are the defaults 

134 for ``medianRef`` and ``matchRef``. 

135 

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

137 """ 

138 

139 blob = Blob('AnalyticAstrometryModel') 

140 for field in ('brightSnrMin', 'brightSnrMax'): 

141 blob[field] = matchedMultiVisitDataset[field] 

142 

143 # FIXME add description field to blobs 

144 # _doc['doc'] \ 

145 # = "Astrometric astrometry model: mas = C*theta/SNR + sigmaSys" 

146 

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

148 medianRef = medianRef * u.marcsec 

149 

150 _compute(blob, 

151 selection, 

152 matchedMultiVisitDataset['snr'].quantity, 

153 matchedMultiVisitDataset['dist'].quantity, 

154 medianRef, matchRef) 

155 return blob 

156 

157 

158def _compute(blob, bright, snr, dist, medianRef, matchRef): 

159 nMatch = len(bright) 

160 median_dist = np.median(dist) 

161 msg = 'Median value of the astrometric scatter - all magnitudes: ' \ 

162 '{0:.3f}' 

163 print(msg.format(median_dist)) 

164 

165 astromScatter = np.median(dist[bright]) 

166 msg = 'Astrometric scatter (median) - snr > {0:.1f} : {1:.1f}' 

167 print(msg.format(blob['brightSnrMin'].quantity, astromScatter)) 

168 

169 fit_params = fitAstromErrModel(snr[bright], dist[bright]) 

170 

171 if astromScatter > medianRef: 

172 msg = 'Median astrometric scatter {0:.1f} is larger than ' \ 

173 'reference : {1:.1f}' 

174 print(msg.format(astromScatter, medianRef)) 

175 if nMatch < matchRef: 

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

177 '(should be > {1:d})' 

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

179 

180 blob['C'] = Datum(quantity=fit_params['C'], 

181 description='Scaling factor') 

182 blob['theta'] = Datum(quantity=fit_params['theta'], 

183 label='theta', 

184 description='Seeing') 

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

186 label='sigma(sys)', 

187 description='Systematic error floor') 

188 blob['astromRms'] = Datum(quantity=astromScatter, 

189 label='RMS', 

190 description='Astrometric scatter (RMS) for good stars')