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

1import functools 

2import math 

3import numpy as np 

4from scipy.stats import norm 

5import astropy.units as u 

6 

7import lsst.pipe.base as pipeBase 

8from lsst.faro.utils.filtermatches import filterMatches 

9 

10__all__ = ("photRepeat", "calcPhotRepeat", "calcPhotRepeatSample", "computeWidths", 

11 "getRandomDiffRmsInMmags", "getRandomDiff") 

12 

13 

14def photRepeat(matchedCatalog, numRandomShuffles=50, randomSeed=None, **filterargs): 

15 filteredCat = filterMatches(matchedCatalog, **filterargs) 

16 magKey = filteredCat.schema.find('slot_PsfFlux_mag').key 

17 

18 # Require at least nMinPhotRepeat objects to calculate the repeatability: 

19 nMinPhotRepeat = 50 

20 if filteredCat.count > nMinPhotRepeat: 

21 phot_resid_meas = calcPhotRepeat(filteredCat, magKey, numRandomShuffles=50, randomSeed=randomSeed) 

22 return phot_resid_meas 

23 else: 

24 return {'nomeas': np.nan*u.mmag} 

25 

26 

27def calcPhotRepeat(matches, magKey, numRandomShuffles=50, randomSeed=None): 

28 """Calculate the photometric repeatability of measurements across a set 

29 of randomly selected pairs of visits. 

30 Parameters 

31 ---------- 

32 matches : `lsst.afw.table.GroupView` 

33 `~lsst.afw.table.GroupView` of sources matched between visits, 

34 from MultiMatch, provided by 

35 `lsst.faro.utils.matcher.match_catalogs`. 

36 magKey : `lsst.afw.table` schema key 

37 Magnitude column key in the ``groupView``. 

38 E.g., ``magKey = allMatches.schema.find("slot_ModelFlux_mag").key`` 

39 where ``allMatches`` is the result of 

40 `lsst.afw.table.MultiMatch.finish()`. 

41 numRandomShuffles : int 

42 Number of times to draw random pairs from the different observations. 

43 randomSeed : int 

44 Seed for random number generation when choosing samples. 

45 Returns 

46 ------- 

47 statistics : `dict` 

48 Statistics to compute model_phot_rep. Fields are: 

49 - ``model_phot_rep``: scalar `~astropy.unit.Quantity` of mean ``iqr``. 

50 This is formally the model_phot_rep metric measurement. 

51 - ``rms``: `~astropy.unit.Quantity` array in mmag of photometric 

52 repeatability RMS across ``numRandomShuffles``. 

53 Shape: ``(nRandomSamples,)``. 

54 - ``iqr``: `~astropy.unit.Quantity` array in mmag of inter-quartile 

55 range of photometric repeatability distribution. 

56 Shape: ``(nRandomSamples,)``. 

57 - ``magDiff``: `~astropy.unit.Quantity` array of magnitude differences 

58 between pairs of sources. Shape: ``(nRandomSamples, nMatches)``. 

59 - ``magMean``: `~astropy.unit.Quantity` array of mean magnitudes of 

60 each pair of sources. Shape: ``(nRandomSamples, nMatches)``. 

61 Notes 

62 ----- 

63 We calculate differences for ``numRandomShuffles`` different random 

64 realizations of the measurement pairs, to provide some estimate of the 

65 uncertainty on our RMS estimates due to the random shuffling. This 

66 estimate could be stated and calculated from a more formally derived 

67 motivation but in practice 50 should be sufficient. 

68 The LSST Science Requirements Document (LPM-17), or SRD, characterizes the 

69 photometric repeatability by putting a requirement on the median RMS of 

70 measurements of non-variable bright stars. This quantity is PA1, with a 

71 design, minimum, and stretch goals of (5, 8, 3) millimag following LPM-17 

72 as of 2011-07-06, available at http://ls.st/LPM-17. model_phot_rep is a 

73 similar quantity measured for extended sources (almost entirely galaxies), 

74 for which no requirement currently exists in the SRD. 

75 This present routine calculates this quantity in two different ways: 

76 1. RMS 

77 2. interquartile range (IQR) 

78 **The repeatability scalar measurement is the median of the IQR.** 

79 This function also returns additional quantities of interest: 

80 - the pair differences of observations of sources, 

81 - the mean magnitude of each source 

82 Examples 

83 -------- 

84 Normally ``calcPhotRepeat`` is called by `measurePhotRepeat`, using 

85 data from `metric_pipeline_utils.matcher.match_catalogs`. Here's an 

86 example of how to call ``calcPhotRepeat`` directly given the Butler output 

87 repository generated by examples/runHscQuickTest.sh: 

88 >>> import lsst.daf.persistence as dafPersist 

89 >>> from lsst.afw.table import SourceCatalog, SchemaMapper, Field 

90 >>> from lsst.afw.table import MultiMatch, SourceRecord, GroupView 

91 >>> from from metric_pipeline_utils.phot_repeat import calcPhotRepeat 

92 >>> from lsst.validate.drp.util import discoverDataIds 

93 >>> import numpy as np 

94 >>> repo = 'HscQuick/output' 

95 >>> butler = dafPersist.Butler(repo) 

96 >>> dataset = 'src' 

97 >>> schema = butler.get(dataset + '_schema', immediate=True).schema 

98 >>> visitDataIds = discoverDataIds(repo) 

99 >>> mmatch = None 

100 >>> for vId in visitDataIds: 

101 ... cat = butler.get('src', vId) 

102 ... calib = butler.get('calexp_photoCalib', vId) 

103 ... cat = calib.calibrateCatalog(cat, ['modelfit_CModel']) 

104 ... if mmatch is None: 

105 ... mmatch = MultiMatch(cat.schema, 

106 ... dataIdFormat={'visit': np.int32, 'ccd': np.int32}, 

107 ... RecordClass=SourceRecord) 

108 ... mmatch.add(catalog=cat, dataId=vId) 

109 ... 

110 >>> matchCat = mmatch.finish() 

111 >>> allMatches = GroupView.build(matchCat) 

112 >>> magKey = allMatches.schema.find('slot_ModelFlux_mag').key 

113 >>> def matchFilter(cat): 

114 >>> if len(cat) < 2: 

115 >>> return False 

116 >>> return np.isfinite(cat.get(magKey)).all() 

117 >>> repeat = calcPhotRepeat(allMatches.where(matchFilter), magKey) 

118 """ 

119 rng = np.random.default_rng(randomSeed) 

120 mprSamples = [calcPhotRepeatSample(matches, magKey, rng=rng) 

121 for _ in range(numRandomShuffles)] 

122 

123 rms = np.array([mpr.rms for mpr in mprSamples]) * u.mmag 

124 iqr = np.array([mpr.iqr for mpr in mprSamples]) * u.mmag 

125 magDiff = np.array([mpr.magDiffs for mpr in mprSamples]) * u.mmag 

126 magMean = np.array([mpr.magMean for mpr in mprSamples]) * u.mag 

127 repeat = np.mean(iqr) 

128 return {'rms': rms, 'iqr': iqr, 'magDiff': magDiff, 'magMean': magMean, 'repeatability': repeat} 

129 

130 

131def calcPhotRepeatSample(matches, magKey, rng=None): 

132 """Compute one realization of repeatability by randomly sampling pairs of 

133 visits. 

134 Parameters 

135 ---------- 

136 matches : `lsst.afw.table.GroupView` 

137 `~lsst.afw.table.GroupView` of sources matched between visits, 

138 from MultiMatch, provided by 

139 `lsst.faro.utils.matcher.match_catalogs`. 

140 magKey : `lsst.afw.table` schema key 

141 Magnitude column key in the ``groupView``. 

142 E.g., ``magKey = allMatches.schema.find("base_PsfFlux_mag").key`` 

143 where ``allMatches`` is the result of 

144 `lsst.afw.table.MultiMatch.finish()`. 

145 rng : `numpy.random._generator.Generator` 

146 Input random number generator. 

147 Returns 

148 ------- 

149 metrics : `lsst.pipe.base.Struct` 

150 Metrics of pairs of sources matched between two visits. Fields are: 

151 - ``rms``: scalar RMS of differences of sources observed in this 

152 randomly sampled pair of visits. 

153 - ``iqr``: scalar inter-quartile range (IQR) of differences of sources 

154 observed in a randomly sampled pair of visits. 

155 - ``magDiffs`: array, shape ``(nMatches,)``, of magnitude differences 

156 (mmag) for observed sources across a randomly sampled pair of visits. 

157 - ``magMean``: array, shape ``(nMatches,)``, of mean magnitudes 

158 of sources observed across a randomly sampled pair of visits. 

159 See also 

160 -------- 

161 calcPhotRepeat : A wrapper that repeatedly calls this function to build 

162 the repeatability measurement. 

163 """ 

164 sampler = functools.partial(getRandomDiffRmsInMmags, rng=rng) 

165 magDiffs = matches.aggregate(sampler, field=magKey) 

166 magMean = matches.aggregate(np.mean, field=magKey) 

167 rms, iqr = computeWidths(magDiffs) 

168 return pipeBase.Struct(rms=rms, iqr=iqr, magDiffs=magDiffs, magMean=magMean,) 

169 

170 

171def computeWidths(array): 

172 """Compute the RMS and the scaled inter-quartile range of an array. 

173 Parameters 

174 ---------- 

175 array : `list` or `numpy.ndarray` 

176 Array. 

177 Returns 

178 ------- 

179 rms : `float` 

180 RMS 

181 iqr : `float` 

182 Scaled inter-quartile range (IQR, see *Notes*). 

183 Notes 

184 ----- 

185 We estimate the width of the histogram in two ways: 

186 - using a simple RMS, 

187 - using the interquartile range (IQR) 

188 The IQR is scaled by the IQR/RMS ratio for a Gaussian such that it 

189 if the array is Gaussian distributed, then the scaled IQR = RMS. 

190 """ 

191 # For scalars, math.sqrt is several times faster than numpy.sqrt. 

192 rmsSigma = math.sqrt(np.mean(array**2)) 

193 iqrSigma = np.subtract.reduce(np.percentile(array, [75, 25])) / (norm.ppf(0.75)*2) 

194 return rmsSigma, iqrSigma 

195 

196 

197def getRandomDiffRmsInMmags(array, rng=None): 

198 """Calculate the RMS difference in mmag between a random pairing of 

199 visits of a source. 

200 Parameters 

201 ---------- 

202 array : `list` or `numpy.ndarray` 

203 Magnitudes from which to select the pair [mag]. 

204 rng : `numpy.random._generator.Generator` 

205 Input random number generator. 

206 Returns 

207 ------- 

208 rmsMmags : `float` 

209 RMS difference in mmag from a random pair of visits. 

210 Notes 

211 ----- 

212 The LSST SRD recommends computing repeatability from a histogram of 

213 magnitude differences for the same source measured on two visits 

214 (using a median over the magDiffs to reject outliers). 

215 Because we have N>=2 measurements for each source, we select a random 

216 pair of visits for each source. We divide each difference by sqrt(2) 

217 to obtain the RMS about the (unknown) mean magnitude, 

218 instead of obtaining just the RMS difference. 

219 See Also 

220 -------- 

221 getRandomDiff : Get the difference between two randomly selected elements of an array. 

222 Examples 

223 -------- 

224 >>> mag = [24.2, 25.5] 

225 >>> rms = getRandomDiffRmsInMmags(mag) 

226 >>> print(rms) 

227 212.132034 

228 """ 

229 thousandDivSqrtTwo = 1000/math.sqrt(2) 

230 return thousandDivSqrtTwo * getRandomDiff(array, rng=rng) 

231 

232 

233def getRandomDiff(array, rng=None): 

234 """Get the difference between two randomly selected elements of an array. 

235 Parameters 

236 ---------- 

237 array : `list` or `numpy.ndarray` 

238 Input array. 

239 rng : `numpy.random._generator.Generator` 

240 Input random number generator. 

241 Returns 

242 ------- 

243 float or int 

244 Difference between two random elements of the array. 

245 """ 

246 if not rng: 

247 rng = np.random.default_rng() 

248 a, b = rng.choice(range(len(array)), 2, replace=False) 

249 return array[a] - array[b]