Coverage for python/lsst/sims/photUtils/SignalToNoise.py : 15%

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 numpy
2from .Sed import Sed
3from . import LSSTdefaults
6__all__ = ["FWHMeff2FWHMgeom", "FWHMgeom2FWHMeff",
7 "calcNeff", "calcInstrNoiseSq", "calcTotalNonSourceNoiseSq", "calcSNR_sed",
8 "calcM5", "calcSkyCountsPerPixelForM5", "calcGamma", "calcSNR_m5",
9 "calcAstrometricError", "magErrorFromSNR", "calcMagError_m5", "calcMagError_sed"]
12def FWHMeff2FWHMgeom(FWHMeff):
13 """
14 Convert FWHMeff to FWHMgeom.
15 This conversion was calculated by Bo Xin and Zeljko Ivezic
16 (and will be in an update on the LSE-40 and overview papers).
18 Parameters
19 ----------
20 FWHMeff: float
21 the single-gaussian equivalent FWHM value, appropriate for calcNeff, in arcseconds
23 Returns
24 -------
25 float
26 FWHM geom, the geometric FWHM value as measured from a typical PSF profile in arcseconds.
27 """
28 FWHMgeom = 0.822*FWHMeff + 0.052
29 return FWHMgeom
32def FWHMgeom2FWHMeff(FWHMgeom):
33 """
34 Convert FWHMgeom to FWHMeff.
35 This conversion was calculated by Bo Xin and Zeljko Ivezic
36 (and will be in an update on the LSE-40 and overview papers).
38 Parameters
39 ----------
40 FWHMgeom: float
41 The geometric FWHM value, as measured from a typical PSF profile, in arcseconds.
43 Returns
44 -------
45 float
46 FWHM effective, the single-gaussian equivalent FWHM value, appropriate for calcNeff, in arcseconds.
47 """
48 FWHMeff = (FWHMgeom - 0.052)/0.822
49 return FWHMeff
52def calcNeff(FWHMeff, platescale):
53 """
54 Calculate the effective number of pixels in a single gaussian PSF.
55 This equation comes from LSE-40, equation 27.
56 https://docushare.lsstcorp.org/docushare/dsweb/ImageStoreViewer/LSE-40
58 Parameters
59 ----------
60 FWHMeff: float
61 The width of a single-gaussian that produces correct Neff for typical PSF profile.
62 platescale: float
63 The platescale in arcseconds per pixel (0.2 for LSST)
65 Returns
66 -------
67 float
68 The effective number of pixels contained in the PSF
70 The FWHMeff is a way to represent the equivalent seeing value, if the
71 atmosphere could be simply represented as a single gaussian (instead of a more
72 complicated von Karman profile for the atmosphere, convolved properly with the
73 telescope hardware additional blurring of 0.4").
74 A translation from the geometric FWHM to the FWHMeff is provided in FWHMgeom2FWHMeff.
75 """
76 return 2.266*(FWHMeff/platescale)**2
79def calcInstrNoiseSq(photParams):
80 """
81 Combine all of the noise due to intrumentation into one value
83 @param [in] photParams is an instantiation of the
84 PhotometricParameters class that carries details about the
85 photometric response of the telescope.
87 @param [out] The noise due to all of these sources added in quadrature
88 in ADU counts
89 """
91 # instrumental squared noise in electrons
92 instNoiseSq = photParams.nexp*photParams.readnoise**2 + \
93 photParams.darkcurrent*photParams.exptime*photParams.nexp + \
94 photParams.nexp*photParams.othernoise**2
96 # convert to ADU counts
97 instNoiseSq = instNoiseSq/(photParams.gain*photParams.gain)
99 return instNoiseSq
102def calcTotalNonSourceNoiseSq(skySed, hardwarebandpass, photParams, FWHMeff):
103 """
104 Calculate the noise due to things that are not the source being observed
105 (i.e. intrumentation and sky background)
107 @param [in] skySed -- an instantiation of the Sed class representing the sky
108 (normalized so that skySed.calcMag() gives the sky brightness in magnitudes
109 per square arcsecond)
111 @param [in] hardwarebandpass -- an instantiation of the Bandpass class representing
112 just the instrumentation throughputs
114 @param [in] photParams is an instantiation of the
115 PhotometricParameters class that carries details about the
116 photometric response of the telescope.
118 @param [in] FWHMeff in arcseconds
120 @param [out] total non-source noise squared (in ADU counts)
121 (this is simga^2_tot * neff in equation 41 of the SNR document
122 https://docushare.lsstcorp.org/docushare/dsweb/ImageStoreViewer/LSE-40 )
123 """
125 # Calculate the effective number of pixels for double-Gaussian PSF
126 neff = calcNeff(FWHMeff, photParams.platescale)
128 # Calculate the counts from the sky.
129 # We multiply by two factors of the platescale because we expect the
130 # skySed to be normalized such that calcADU gives counts per
131 # square arc second, and we need to convert to counts per pixel.
133 skycounts = skySed.calcADU(hardwarebandpass, photParams=photParams) \
134 * photParams.platescale * photParams.platescale
136 # Calculate the square of the noise due to instrumental effects.
137 # Include the readout noise as many times as there are exposures
139 noise_instr_sq = calcInstrNoiseSq(photParams=photParams)
141 # Calculate the square of the noise due to sky background poisson noise
142 noise_sky_sq = skycounts/photParams.gain
144 # Discount error in sky measurement for now
145 noise_skymeasurement_sq = 0
147 total_noise_sq = neff*(noise_sky_sq + noise_instr_sq + noise_skymeasurement_sq)
149 return total_noise_sq
152def calcSkyCountsPerPixelForM5(m5target, totalBandpass, photParams, FWHMeff=None):
153 """
154 Calculate the number of sky counts per pixel expected for a given
155 value of the 5-sigma limiting magnitude (m5)
157 The 5-sigma limiting magnitude (m5) for an observation is
158 determined by a combination of the telescope and camera parameters
159 (such as diameter of the mirrors and the readnoise) together with the
160 sky background.
162 @param [in] the desired value of m5
164 @param [in] totalBandpass is an instantiation of the Bandpass class
165 representing the total throughput of the telescope (instrumentation
166 plus atmosphere)
168 @param [in] photParams is an instantiation of the
169 PhotometricParameters class that carries details about the
170 photometric response of the telescope.
172 @param [in] FWHMeff in arcseconds
174 @param [out] returns the expected number of sky counts per pixel
175 """
177 if FWHMeff is None:
178 FWHMeff = LSSTdefaults().FWHMeff('r')
180 # instantiate a flat SED
181 flatSed = Sed()
182 flatSed.setFlatSED()
184 # normalize the SED so that it has a magnitude equal to the desired m5
185 fNorm = flatSed.calcFluxNorm(m5target, totalBandpass)
186 flatSed.multiplyFluxNorm(fNorm)
187 sourceCounts = flatSed.calcADU(totalBandpass, photParams=photParams)
189 # calculate the effective number of pixels for a double-Gaussian PSF
190 neff = calcNeff(FWHMeff, photParams.platescale)
192 # calculate the square of the noise due to the instrument
193 noise_instr_sq = calcInstrNoiseSq(photParams=photParams)
195 # now solve equation 41 of the SNR document for the neff * sigma_total^2 term
196 # given snr=5 and counts as calculated above
197 # SNR document can be found at
198 # https://docushare.lsstcorp.org/docushare/dsweb/ImageStoreViewer/LSE-40
200 nSigmaSq = (sourceCounts*sourceCounts)/25.0 - sourceCounts/photParams.gain
202 skyNoiseTarget = nSigmaSq/neff - noise_instr_sq
203 skyCountsTarget = skyNoiseTarget*photParams.gain
205 # TODO:
206 # This method should throw an error if skyCountsTarget is negative
207 # unfortunately, that currently happens for default values of
208 # m5 as taken from arXiv:0805.2366, table 2. Adding the error
209 # should probably wait for a later issue in which we hash out what
210 # the units are for all of the parameters stored in PhotometricDefaults.
212 return skyCountsTarget
215def calcM5(skysed, totalBandpass, hardware, photParams, FWHMeff=None):
216 """
217 Calculate the AB magnitude of a 5-sigma above sky background source.
219 The 5-sigma limiting magnitude (m5) for an observation is determined by
220 a combination of the telescope and camera parameters (such as diameter
221 of the mirrors and the readnoise) together with the sky background. This
222 method (calcM5) calculates the expected m5 value for an observation given
223 a sky background Sed and hardware parameters.
225 @param [in] skysed is an instantiation of the Sed class representing
226 sky emission, normalized so that skysed.calcMag gives the sky brightness
227 in magnitudes per square arcsecond.
229 @param [in] totalBandpass is an instantiation of the Bandpass class
230 representing the total throughput of the telescope (instrumentation
231 plus atmosphere)
233 @param [in] hardware is an instantiation of the Bandpass class representing
234 the throughput due solely to instrumentation.
236 @param [in] photParams is an instantiation of the
237 PhotometricParameters class that carries details about the
238 photometric response of the telescope.
240 @param [in] FWHMeff in arcseconds
242 @param [out] returns the value of m5 for the given bandpass and sky SED
243 """
244 # This comes from equation 45 of the SNR document (v1.2, May 2010)
245 # https://docushare.lsstcorp.org/docushare/dsweb/ImageStoreViewer/LSE-40
247 if FWHMeff is None:
248 FWHMeff = LSSTdefaults().FWHMeff('r')
250 # create a flat fnu source
251 flatsource = Sed()
252 flatsource.setFlatSED()
253 snr = 5.0
254 v_n = calcTotalNonSourceNoiseSq(skysed, hardware, photParams, FWHMeff)
256 counts_5sigma = (snr**2)/2.0/photParams.gain + \
257 numpy.sqrt((snr**4)/4.0/photParams.gain + (snr**2)*v_n)
259 # renormalize flatsource so that it has the required counts to be a 5-sigma detection
260 # given the specified background
261 counts_flat = flatsource.calcADU(totalBandpass, photParams=photParams)
262 flatsource.multiplyFluxNorm(counts_5sigma/counts_flat)
264 # Calculate the AB magnitude of this source.
265 mag_5sigma = flatsource.calcMag(totalBandpass)
266 return mag_5sigma
269def magErrorFromSNR(snr):
270 """
271 convert flux signal to noise ratio to an error in magnitude
273 @param [in] snr is the signal to noise ratio in flux
275 @param [out] the resulting error in magnitude
276 """
278 #see www.ucolick.org/~bolte/AY257/s_n.pdf section 3.1
279 return 2.5*numpy.log10(1.0+1.0/snr)
282def calcGamma(bandpass, m5, photParams):
284 """
285 Calculate the gamma parameter used for determining photometric
286 signal to noise in equation 5 of the LSST overview paper
287 (arXiv:0805.2366)
289 @param [in] bandpass is an instantiation of the Bandpass class
290 representing the bandpass for which you desire to calculate the
291 gamma parameter
293 @param [in] m5 is the magnitude at which a 5-sigma detection occurs
294 in this Bandpass
296 @param [in] photParams is an instantiation of the
297 PhotometricParameters class that carries details about the
298 photometric response of the telescope.
300 @param [out] gamma
301 """
302 # This is based on the LSST SNR document (v1.2, May 2010)
303 # https://docushare.lsstcorp.org/docushare/dsweb/ImageStoreViewer/LSE-40
304 # as well as equations 4-6 of the overview paper (arXiv:0805.2366)
306 # instantiate a flat SED
307 flatSed = Sed()
308 flatSed.setFlatSED()
310 # normalize the SED so that it has a magnitude equal to the desired m5
311 fNorm = flatSed.calcFluxNorm(m5, bandpass)
312 flatSed.multiplyFluxNorm(fNorm)
313 counts = flatSed.calcADU(bandpass, photParams=photParams)
315 # The expression for gamma below comes from:
316 #
317 # 1) Take the approximation N^2 = N0^2 + alpha S from footnote 88 in the overview paper
318 # where N is the noise in flux of a source, N0 is the noise in flux due to sky brightness
319 # and instrumentation, S is the number of counts registered from the source and alpha
320 # is some constant
321 #
322 # 2) Divide by S^2 and demand that N/S = 0.2 for a source detected at m5. Solve
323 # the resulting equation for alpha in terms of N0 and S5 (the number of counts from
324 # a source at m5)
325 #
326 # 3) Substitute this expression for alpha back into the equation for (N/S)^2
327 # for a general source. Re-factor the equation so that it looks like equation
328 # 5 of the overview paper (note that x = S5/S). This should give you gamma = (N0/S5)^2
329 #
330 # 4) Solve equation 41 of the SNR document for the neff * sigma_total^2 term
331 # given snr=5 and counts as calculated above. Note that neff * sigma_total^2
332 # is N0^2 in the equation above
333 #
334 # This should give you
336 gamma = 0.04 - 1.0/(counts*photParams.gain)
338 return gamma
341def calcSNR_m5(magnitude, bandpass, m5, photParams, gamma=None):
342 """
343 Calculate signal to noise in flux using the model from equation (5) of arXiv:0805.2366
345 @param [in] magnitude of the sources whose signal to noise you are calculating
346 (can be a numpy array)
348 @param [in] bandpass (an instantiation of the class Bandpass) in which the magnitude
349 was calculated
351 @param [in] m5 is the 5-sigma limiting magnitude for the bandpass
353 @param [in] photParams is an instantiation of the
354 PhotometricParameters class that carries details about the
355 photometric response of the telescope.
357 @param [in] gamma (optional) is the gamma parameter from equation(5) of
358 arXiv:0805.2366. If not provided, this method will calculate it.
360 @param [out] snr is the signal to noise ratio corresponding to
361 the input magnitude.
363 @param [out] gamma is the calculated gamma parameter for the
364 bandpass used here (in case the user wants to call this method again).
366 Note: You can also pass in a numpy array of magnitudes calculated
367 in the same bandpass with the same m5 and get a numpy array of SNR out.
368 """
370 if gamma is None:
371 gamma = calcGamma(bandpass, m5, photParams=photParams)
373 dummySed = Sed()
374 m5Flux = dummySed.fluxFromMag(m5)
375 sourceFlux = dummySed.fluxFromMag(magnitude)
377 fluxRatio = m5Flux/sourceFlux
379 noise = numpy.sqrt((0.04-gamma)*fluxRatio+gamma*fluxRatio*fluxRatio)
381 return 1.0/noise, gamma
384def calcMagError_m5(magnitude, bandpass, m5, photParams, gamma=None):
385 """
386 Calculate magnitude error using the model from equation (5) of arXiv:0805.2366
388 @param [in] magnitude of the source whose error you want
389 to calculate (can be a numpy array)
391 @param [in] bandpass (an instantiation of the Bandpass class) in question
393 @param [in] m5 is the 5-sigma limiting magnitude in that bandpass
395 @param [in] photParams is an instantiation of the
396 PhotometricParameters class that carries details about the
397 photometric response of the telescope.
399 @param [in] gamma (optional) is the gamma parameter from equation(5) of
400 arXiv:0805.2366. If not provided, this method will calculate it.
402 @param [out] the error associated with the magnitude
404 @param [out] gamma is the calculated gamma parameter for the
405 bandpass used here (in case the user wants to call this method again).
407 Note: you can also pass in a numpy of array of magnitudes calculated in
408 the same Bandpass with the same m5 and get a numpy array of errors out.
409 """
411 snr, gamma = calcSNR_m5(magnitude, bandpass, m5, photParams, gamma=gamma)
413 if photParams.sigmaSys is not None:
414 return numpy.sqrt(numpy.power(magErrorFromSNR(snr),2) + numpy.power(photParams.sigmaSys,2)), gamma
415 else:
416 return magErrorFromSNR(snr), gamma
419def calcSNR_sed(sourceSed, totalbandpass, skysed, hardwarebandpass,
420 photParams, FWHMeff, verbose=False):
421 """
422 Calculate the signal to noise ratio for a source, given the bandpass(es) and sky SED.
424 For a given source, sky sed, total bandpass and hardware bandpass, as well as
425 FWHMeff / exptime, calculates the SNR with optimal PSF extraction
426 assuming a double-gaussian PSF.
428 @param [in] sourceSed is an instantiation of the Sed class containing the SED of
429 the object whose signal to noise ratio is being calculated
431 @param [in] totalbandpass is an instantiation of the Bandpass class
432 representing the total throughput (system + atmosphere)
434 @param [in] skysed is an instantiation of the Sed class representing
435 the sky emission per square arcsecond.
437 @param [in] hardwarebandpass is an instantiation of the Bandpass class
438 representing just the throughput of the system hardware.
440 @param [in] photParams is an instantiation of the
441 PhotometricParameters class that carries details about the
442 photometric response of the telescope.
444 @param [in] FWHMeff in arcseconds
446 @param [in] verbose is a boolean
448 @param [out] signal to noise ratio
449 """
451 # Calculate the counts from the source.
452 sourcecounts = sourceSed.calcADU(totalbandpass, photParams=photParams)
454 # Calculate the (square of the) noise due to signal poisson noise.
455 noise_source_sq = sourcecounts/photParams.gain
457 non_source_noise_sq = calcTotalNonSourceNoiseSq(skysed, hardwarebandpass, photParams, FWHMeff)
459 # Calculate total noise
460 noise = numpy.sqrt(noise_source_sq + non_source_noise_sq)
461 # Calculate the signal to noise ratio.
462 snr = sourcecounts / noise
463 if verbose:
464 skycounts = skysed.calcADU(hardwarebandpass, photParams) * (photParams.platescale**2)
465 noise_sky_sq = skycounts/photParams.gain
466 neff = calcNeff(FWHMeff, photParams.platescale)
467 noise_instr_sq = calcInstrNoiseSq(photParams)
469 print("For Nexp %.1f of time %.1f: " % (photParams.nexp, photParams.exptime))
470 print("Counts from source: %.2f Counts from sky: %.2f" %(sourcecounts, skycounts))
471 print("FWHMeff: %.2f('') Neff pixels: %.3f(pix)" %(FWHMeff, neff))
472 print("Noise from sky: %.2f Noise from instrument: %.2f" \
473 %(numpy.sqrt(noise_sky_sq), numpy.sqrt(noise_instr_sq)))
474 print("Noise from source: %.2f" %(numpy.sqrt(noise_source_sq)))
475 print(" Total Signal: %.2f Total Noise: %.2f SNR: %.2f" %(sourcecounts, noise, snr))
476 # Return the signal to noise value.
477 return snr
480def calcMagError_sed(sourceSed, totalbandpass, skysed, hardwarebandpass,
481 photParams, FWHMeff, verbose=False):
482 """
483 Calculate the magnitudeError for a source, given the bandpass(es) and sky SED.
485 For a given source, sky sed, total bandpass and hardware bandpass, as well as
486 FWHMeff / exptime, calculates the SNR with optimal PSF extraction
487 assuming a double-gaussian PSF.
489 @param [in] sourceSed is an instantiation of the Sed class containing the SED of
490 the object whose signal to noise ratio is being calculated
492 @param [in] totalbandpass is an instantiation of the Bandpass class
493 representing the total throughput (system + atmosphere)
495 @param [in] skysed is an instantiation of the Sed class representing
496 the sky emission per square arcsecond.
498 @param [in] hardwarebandpass is an instantiation of the Bandpass class
499 representing just the throughput of the system hardware.
501 @param [in] photParams is an instantiation of the
502 PhotometricParameters class that carries details about the
503 photometric response of the telescope.
505 @param [in] FWHMeff in arcseconds
507 @param [in] verbose is a boolean
509 @param [out] magnitude error
510 """
512 snr = calcSNR_sed(sourceSed, totalbandpass, skysed, hardwarebandpass,
513 photParams, FWHMeff, verbose=verbose)
515 if photParams.sigmaSys is not None:
516 return numpy.sqrt(numpy.power(magErrorFromSNR(snr),2) + numpy.power(photParams.sigmaSys,2))
517 else:
518 return magErrorFromSNR(snr)
521def calcAstrometricError(mag, m5, fwhmGeom=0.7, nvisit=1, systematicFloor=10):
522 """
523 Calculate an expected astrometric error.
524 Can be used to estimate this for general catalog purposes (use typical FWHM and nvisit=Number of visit).
525 Or can be used for a single visit, use actual FWHM and nvisit=1.
527 Parameters
528 ----------
529 mag: float
530 Magnitude of the source
531 m5: float
532 Point source five sigma limiting magnitude of the image (or typical depth per image).
533 fwhmGeom: float, opt
534 The geometric (physical) FWHM of the image, in arcseconds. Default 0.7.
535 nvisit: int, opt
536 The number of visits/measurement. Default 1.
537 If this is >1, the random error contribution is reduced by sqrt(nvisits).
538 systematicFloor: float, opt
539 The systematic noise floor for the astrometric measurements, in mas. Default 10mas.
541 Returns
542 -------
543 float
544 Astrometric error for a given SNR, in mas.
545 """
546 # The astrometric error can be applied to parallax or proper motion (for nvisit>1).
547 # If applying to proper motion, should also divide by the # of years of the survey.
548 # This is also referenced in the astroph/0805.2366 paper.
549 # D. Monet suggests sqrt(Nvisit/2) for first 3 years, sqrt(N) for longer, in reduction of error
550 # because of the astrometric measurement method, the systematic and random error are both reduced.
551 # Zeljko says 'be conservative', so removing this reduction for now.
552 rgamma = 0.039
553 xval = numpy.power(10, 0.4*(mag-m5))
554 # The average FWHMeff is 0.7" (or 700 mas), but user can specify. Convert to mas.
555 seeing = fwhmGeom * 1000.
556 error_rand = seeing * numpy.sqrt((0.04-rgamma)*xval + rgamma*xval*xval)
557 error_rand = error_rand / numpy.sqrt(nvisit)
558 # The systematic error floor in astrometry (mas).
559 error_sys = systematicFloor
560 astrom_error = numpy.sqrt(error_sys * error_sys + error_rand*error_rand)
561 return astrom_error