lsst.pipe.tasks  15.0-19-gcae3f1ab+2
photoCal.py
Go to the documentation of this file.
1 #
2 # LSST Data Management System
3 # Copyright 2008-2016 AURA/LSST.
4 #
5 # This product includes software developed by the
6 # LSST Project (http://www.lsst.org/).
7 #
8 # This program is free software: you can redistribute it and/or modify
9 # it under the terms of the GNU General Public License as published by
10 # the Free Software Foundation, either version 3 of the License, or
11 # (at your option) any later version.
12 #
13 # This program is distributed in the hope that it will be useful,
14 # but WITHOUT ANY WARRANTY; without even the implied warranty of
15 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 # GNU General Public License for more details.
17 #
18 # You should have received a copy of the LSST License Statement and
19 # the GNU General Public License along with this program. If not,
20 # see <https://www.lsstcorp.org/LegalNotices/>.
21 #
22 # @package lsst.pipe.tasks.
23 import math
24 import sys
25 
26 import numpy as np
27 
28 import lsst.pex.config as pexConf
29 import lsst.pipe.base as pipeBase
30 from lsst.afw.image import abMagFromFlux, abMagErrFromFluxErr, Calib
31 from lsst.meas.astrom import DirectMatchTask, DirectMatchConfigWithoutLoader
32 import lsst.afw.display.ds9 as ds9
33 from lsst.meas.algorithms import getRefFluxField, ReserveSourcesTask
34 from .colorterms import ColortermLibrary
35 
36 __all__ = ["PhotoCalTask", "PhotoCalConfig"]
37 
38 
39 class PhotoCalConfig(pexConf.Config):
40  """Config for PhotoCal"""
41  match = pexConf.ConfigField("Match to reference catalog",
42  DirectMatchConfigWithoutLoader)
43  reserve = pexConf.ConfigurableField(target=ReserveSourcesTask, doc="Reserve sources from fitting")
44  fluxField = pexConf.Field(
45  dtype=str,
46  default="slot_CalibFlux_flux",
47  doc=("Name of the source flux field to use. The associated flag field\n"
48  "('<name>_flags') will be implicitly included in badFlags."),
49  )
50  applyColorTerms = pexConf.Field(
51  dtype=bool,
52  default=None,
53  doc=("Apply photometric color terms to reference stars? One of:\n"
54  "None: apply if colorterms and photoCatName are not None;\n"
55  " fail if color term data is not available for the specified ref catalog and filter.\n"
56  "True: always apply colorterms; fail if color term data is not available for the\n"
57  " specified reference catalog and filter.\n"
58  "False: do not apply."),
59  optional=True,
60  )
61  sigmaMax = pexConf.Field(
62  dtype=float,
63  default=0.25,
64  doc="maximum sigma to use when clipping",
65  optional=True,
66  )
67  nSigma = pexConf.Field(
68  dtype=float,
69  default=3.0,
70  doc="clip at nSigma",
71  )
72  useMedian = pexConf.Field(
73  dtype=bool,
74  default=True,
75  doc="use median instead of mean to compute zeropoint",
76  )
77  nIter = pexConf.Field(
78  dtype=int,
79  default=20,
80  doc="number of iterations",
81  )
82  colorterms = pexConf.ConfigField(
83  dtype=ColortermLibrary,
84  doc="Library of photometric reference catalog name: color term dict",
85  )
86  photoCatName = pexConf.Field(
87  dtype=str,
88  optional=True,
89  doc=("Name of photometric reference catalog; used to select a color term dict in colorterms."
90  " see also applyColorTerms"),
91  )
92  magErrFloor = pexConf.RangeField(
93  dtype=float,
94  default=0.0,
95  doc="Additional magnitude uncertainty to be added in quadrature with measurement errors.",
96  min=0.0,
97  )
98 
99  def validate(self):
100  pexConf.Config.validate(self)
101  if self.applyColorTerms and self.photoCatName is None:
102  raise RuntimeError("applyColorTerms=True requires photoCatName is non-None")
103  if self.applyColorTerms and len(self.colorterms.data) == 0:
104  raise RuntimeError("applyColorTerms=True requires colorterms be provided")
105 
106  def setDefaults(self):
107  pexConf.Config.setDefaults(self)
108  self.match.sourceSelection.doFlags = True
109  self.match.sourceSelection.flags.bad = [
110  "base_PixelFlags_flag_edge",
111  "base_PixelFlags_flag_interpolated",
112  "base_PixelFlags_flag_saturated",
113  ]
114  self.match.sourceSelection.doUnresolved = True
115 
116 
117 
123 
124 class PhotoCalTask(pipeBase.Task):
125  r"""!
126 @anchor PhotoCalTask_
127 
128 @brief Calculate the zero point of an exposure given a lsst.afw.table.ReferenceMatchVector.
129 
130 @section pipe_tasks_photocal_Contents Contents
131 
132  - @ref pipe_tasks_photocal_Purpose
133  - @ref pipe_tasks_photocal_Initialize
134  - @ref pipe_tasks_photocal_IO
135  - @ref pipe_tasks_photocal_Config
136  - @ref pipe_tasks_photocal_Debug
137  - @ref pipe_tasks_photocal_Example
138 
139 @section pipe_tasks_photocal_Purpose Description
140 
141 @copybrief PhotoCalTask
142 
143 Calculate an Exposure's zero-point given a set of flux measurements of stars matched to an input catalogue.
144 The type of flux to use is specified by PhotoCalConfig.fluxField.
145 
146 The algorithm clips outliers iteratively, with parameters set in the configuration.
147 
148 @note This task can adds fields to the schema, so any code calling this task must ensure that
149 these columns are indeed present in the input match list; see @ref pipe_tasks_photocal_Example
150 
151 @section pipe_tasks_photocal_Initialize Task initialisation
152 
153 @copydoc \_\_init\_\_
154 
155 @section pipe_tasks_photocal_IO Inputs/Outputs to the run method
156 
157 @copydoc run
158 
159 @section pipe_tasks_photocal_Config Configuration parameters
160 
161 See @ref PhotoCalConfig
162 
163 @section pipe_tasks_photocal_Debug Debug variables
164 
165 The @link lsst.pipe.base.cmdLineTask.CmdLineTask command line task@endlink interface supports a
166 flag @c -d to import @b debug.py from your @c PYTHONPATH; see @ref baseDebug for more about @b debug.py files.
167 
168 The available variables in PhotoCalTask are:
169 <DL>
170  <DT> @c display
171  <DD> If True enable other debug outputs
172  <DT> @c displaySources
173  <DD> If True, display the exposure on ds9's frame 1 and overlay the source catalogue.
174  <DL>
175  <DT> red o
176  <DD> Reserved objects
177  <DT> green o
178  <DD> Objects used in the photometric calibration
179  </DL>
180  <DT> @c scatterPlot
181  <DD> Make a scatter plot of flux v. reference magnitude as a function of reference magnitude.
182  - good objects in blue
183  - rejected objects in red
184  (if @c scatterPlot is 2 or more, prompt to continue after each iteration)
185 </DL>
186 
187 @section pipe_tasks_photocal_Example A complete example of using PhotoCalTask
188 
189 This code is in @link examples/photoCalTask.py@endlink, and can be run as @em e.g.
190 @code
191 examples/photoCalTask.py
192 @endcode
193 @dontinclude photoCalTask.py
194 
195 Import the tasks (there are some other standard imports; read the file for details)
196 @skipline from lsst.pipe.tasks.astrometry
197 @skipline measPhotocal
198 
199 We need to create both our tasks before processing any data as the task constructors
200 can add extra columns to the schema which we get from the input catalogue, @c scrCat:
201 @skipline getSchema
202 
203 Astrometry first:
204 @skip AstrometryTask.ConfigClass
205 @until aTask
206 (that @c filterMap line is because our test code doesn't use a filter that the reference catalogue recognises,
207 so we tell it to use the @c r band)
208 
209 Then photometry:
210 @skip measPhotocal
211 @until pTask
212 
213 If the schema has indeed changed we need to add the new columns to the source table
214 (yes; this should be easier!)
215 @skip srcCat
216 @until srcCat = cat
217 
218 We're now ready to process the data (we could loop over multiple exposures/catalogues using the same
219 task objects):
220 @skip matches
221 @until result
222 
223 We can then unpack and use the results:
224 @skip calib
225 @until np.log
226 
227 <HR>
228 To investigate the @ref pipe_tasks_photocal_Debug, put something like
229 @code{.py}
230  import lsstDebug
231  def DebugInfo(name):
232  di = lsstDebug.getInfo(name) # N.b. lsstDebug.Info(name) would call us recursively
233  if name.endswith(".PhotoCal"):
234  di.display = 1
235 
236  return di
237 
238  lsstDebug.Info = DebugInfo
239 @endcode
240 into your debug.py file and run photoCalTask.py with the @c --debug flag.
241  """
242  ConfigClass = PhotoCalConfig
243  _DefaultName = "photoCal"
244 
245  def __init__(self, refObjLoader, schema=None, **kwds):
246  """!Create the photometric calibration task. See PhotoCalTask.init for documentation
247  """
248  pipeBase.Task.__init__(self, **kwds)
249  self.scatterPlot = None
250  self.fig = None
251  if schema is not None:
252  self.usedKey = schema.addField("calib_photometry_used", type="Flag",
253  doc="set if source was used in photometric calibration")
254  else:
255  self.usedKey = None
256  self.match = DirectMatchTask(self.config.match, refObjLoader=refObjLoader,
257  name="match", parentTask=self)
258  self.makeSubtask("reserve", columnName="calib_photometry", schema=schema,
259  doc="set if source was reserved from photometric calibration")
260 
261  def getSourceKeys(self, schema):
262  """!Return a struct containing the source catalog keys for fields used by PhotoCalTask.
263 
264  Returned fields include:
265  - flux
266  - fluxErr
267  """
268  flux = schema.find(self.config.fluxField).key
269  fluxErr = schema.find(self.config.fluxField + "Sigma").key
270  return pipeBase.Struct(flux=flux, fluxErr=fluxErr)
271 
272  @pipeBase.timeMethod
273  def extractMagArrays(self, matches, filterName, sourceKeys):
274  """!Extract magnitude and magnitude error arrays from the given matches.
275 
276  @param[in] matches Reference/source matches, a @link lsst::afw::table::ReferenceMatchVector@endlink
277  @param[in] filterName Name of filter being calibrated
278  @param[in] sourceKeys Struct of source catalog keys, as returned by getSourceKeys()
279 
280  @return Struct containing srcMag, refMag, srcMagErr, refMagErr, and magErr numpy arrays
281  where magErr is an error in the magnitude; the error in srcMag - refMag
282  If nonzero, config.magErrFloor will be added to magErr *only* (not srcMagErr or refMagErr), as
283  magErr is what is later used to determine the zero point.
284  Struct also contains refFluxFieldList: a list of field names of the reference catalog used for fluxes
285  (1 or 2 strings)
286  @note These magnitude arrays are the @em inputs to the photometric calibration, some may have been
287  discarded by clipping while estimating the calibration (https://jira.lsstcorp.org/browse/DM-813)
288  """
289  srcFluxArr = np.array([m.second.get(sourceKeys.flux) for m in matches])
290  srcFluxErrArr = np.array([m.second.get(sourceKeys.fluxErr) for m in matches])
291  if not np.all(np.isfinite(srcFluxErrArr)):
292  # this is an unpleasant hack; see DM-2308 requesting a better solution
293  self.log.warn("Source catalog does not have flux uncertainties; using sqrt(flux).")
294  srcFluxErrArr = np.sqrt(srcFluxArr)
295 
296  # convert source flux from DN to an estimate of Jy
297  JanskysPerABFlux = 3631.0
298  srcFluxArr = srcFluxArr * JanskysPerABFlux
299  srcFluxErrArr = srcFluxErrArr * JanskysPerABFlux
300 
301  if not matches:
302  raise RuntimeError("No reference stars are available")
303  refSchema = matches[0].first.schema
304 
305  applyColorTerms = self.config.applyColorTerms
306  applyCTReason = "config.applyColorTerms is %s" % (self.config.applyColorTerms,)
307  if self.config.applyColorTerms is None:
308  # apply color terms if color term data is available and photoCatName specified
309  ctDataAvail = len(self.config.colorterms.data) > 0
310  photoCatSpecified = self.config.photoCatName is not None
311  applyCTReason += " and data %s available" % ("is" if ctDataAvail else "is not")
312  applyCTReason += " and photoRefCat %s provided" % ("is" if photoCatSpecified else "is not")
313  applyColorTerms = ctDataAvail and photoCatSpecified
314 
315  if applyColorTerms:
316  self.log.info("Applying color terms for filterName=%r, config.photoCatName=%s because %s",
317  filterName, self.config.photoCatName, applyCTReason)
318  ct = self.config.colorterms.getColorterm(
319  filterName=filterName, photoCatName=self.config.photoCatName, doRaise=True)
320  else:
321  self.log.info("Not applying color terms because %s", applyCTReason)
322  ct = None
323 
324  if ct: # we have a color term to worry about
325  fluxFieldList = [getRefFluxField(refSchema, filt) for filt in (ct.primary, ct.secondary)]
326  missingFluxFieldList = []
327  for fluxField in fluxFieldList:
328  try:
329  refSchema.find(fluxField).key
330  except KeyError:
331  missingFluxFieldList.append(fluxField)
332 
333  if missingFluxFieldList:
334  self.log.warn("Source catalog does not have fluxes for %s; ignoring color terms",
335  " ".join(missingFluxFieldList))
336  ct = None
337 
338  if not ct:
339  fluxFieldList = [getRefFluxField(refSchema, filterName)]
340 
341  refFluxArrList = [] # list of ref arrays, one per flux field
342  refFluxErrArrList = [] # list of ref flux arrays, one per flux field
343  for fluxField in fluxFieldList:
344  fluxKey = refSchema.find(fluxField).key
345  refFluxArr = np.array([m.first.get(fluxKey) for m in matches])
346  try:
347  fluxErrKey = refSchema.find(fluxField + "Sigma").key
348  refFluxErrArr = np.array([m.first.get(fluxErrKey) for m in matches])
349  except KeyError:
350  # Reference catalogue may not have flux uncertainties; HACK
351  self.log.warn("Reference catalog does not have flux uncertainties for %s; using sqrt(flux).",
352  fluxField)
353  refFluxErrArr = np.sqrt(refFluxArr)
354 
355  refFluxArrList.append(refFluxArr)
356  refFluxErrArrList.append(refFluxErrArr)
357 
358  if ct: # we have a color term to worry about
359  refMagArr1 = np.array([abMagFromFlux(rf1) for rf1 in refFluxArrList[0]]) # primary
360  refMagArr2 = np.array([abMagFromFlux(rf2) for rf2 in refFluxArrList[1]]) # secondary
361 
362  refMagArr = ct.transformMags(refMagArr1, refMagArr2)
363  refFluxErrArr = ct.propagateFluxErrors(refFluxErrArrList[0], refFluxErrArrList[1])
364  else:
365  refMagArr = np.array([abMagFromFlux(rf) for rf in refFluxArrList[0]])
366 
367  srcMagArr = np.array([abMagFromFlux(sf) for sf in srcFluxArr])
368 
369  # Fitting with error bars in both axes is hard
370  # for now ignore reference flux error, but ticket DM-2308 is a request for a better solution
371  magErrArr = np.array([abMagErrFromFluxErr(fe, sf) for fe, sf in zip(srcFluxErrArr, srcFluxArr)])
372  if self.config.magErrFloor != 0.0:
373  magErrArr = (magErrArr**2 + self.config.magErrFloor**2)**0.5
374 
375  srcMagErrArr = np.array([abMagErrFromFluxErr(sfe, sf) for sfe, sf in zip(srcFluxErrArr, srcFluxArr)])
376  refMagErrArr = np.array([abMagErrFromFluxErr(rfe, rf) for rfe, rf in zip(refFluxErrArr, refFluxArr)])
377 
378  good = np.isfinite(srcMagArr) & np.isfinite(refMagArr)
379 
380  return pipeBase.Struct(
381  srcMag=srcMagArr[good],
382  refMag=refMagArr[good],
383  magErr=magErrArr[good],
384  srcMagErr=srcMagErrArr[good],
385  refMagErr=refMagErrArr[good],
386  refFluxFieldList=fluxFieldList,
387  )
388 
389  @pipeBase.timeMethod
390  def run(self, exposure, sourceCat, expId=0):
391  """!Do photometric calibration - select matches to use and (possibly iteratively) compute
392  the zero point.
393 
394  @param[in] exposure Exposure upon which the sources in the matches were detected.
395  @param[in] sourceCat A catalog of sources to use in the calibration
396  (@em i.e. a list of lsst.afw.table.Match with
397  @c first being of type lsst.afw.table.SimpleRecord and @c second type lsst.afw.table.SourceRecord ---
398  the reference object and matched object respectively).
399  (will not be modified except to set the outputField if requested.).
400 
401  @return Struct of:
402  - calib ------- @link lsst::afw::image::Calib@endlink object containing the zero point
403  - arrays ------ Magnitude arrays returned be PhotoCalTask.extractMagArrays
404  - matches ----- Final ReferenceMatchVector, as returned by PhotoCalTask.selectMatches.
405  - zp ---------- Photometric zero point (mag)
406  - sigma ------- Standard deviation of fit of photometric zero point (mag)
407  - ngood ------- Number of sources used to fit photometric zero point
408 
409  The exposure is only used to provide the name of the filter being calibrated (it may also be
410  used to generate debugging plots).
411 
412  The reference objects:
413  - Must include a field @c photometric; True for objects which should be considered as
414  photometric standards
415  - Must include a field @c flux; the flux used to impose a magnitude limit and also to calibrate
416  the data to (unless a color term is specified, in which case ColorTerm.primary is used;
417  See https://jira.lsstcorp.org/browse/DM-933)
418  - May include a field @c stargal; if present, True means that the object is a star
419  - May include a field @c var; if present, True means that the object is variable
420 
421  The measured sources:
422  - Must include PhotoCalConfig.fluxField; the flux measurement to be used for calibration
423 
424  @throws RuntimeError with the following strings:
425 
426  <DL>
427  <DT> No matches to use for photocal
428  <DD> No matches are available (perhaps no sources/references were selected by the matcher).
429  <DT> No reference stars are available
430  <DD> No matches are available from which to extract magnitudes.
431  </DL>
432  """
433  import lsstDebug
434 
435  display = lsstDebug.Info(__name__).display
436  displaySources = display and lsstDebug.Info(__name__).displaySources
437  self.scatterPlot = display and lsstDebug.Info(__name__).scatterPlot
438 
439  if self.scatterPlot:
440  from matplotlib import pyplot
441  try:
442  self.fig.clf()
443  except Exception:
444  self.fig = pyplot.figure()
445 
446  filterName = exposure.getFilter().getName()
447 
448  # Match sources
449  matchResults = self.match.run(sourceCat, filterName)
450  matches = matchResults.matches
451  reserveResults = self.reserve.run([mm.second for mm in matches], expId=expId)
452  if displaySources:
453  self.displaySources(exposure, matches, reserveResults.reserved)
454  if reserveResults.reserved.sum() > 0:
455  matches = [mm for mm, use in zip(matches, reserveResults.use) if use]
456  if len(matches) == 0:
457  raise RuntimeError("No matches to use for photocal")
458  if self.usedKey is not None:
459  for mm in matches:
460  mm.second.set(self.usedKey, True)
461 
462  # Prepare for fitting
463  sourceKeys = self.getSourceKeys(matches[0].second.schema)
464  arrays = self.extractMagArrays(matches=matches, filterName=filterName, sourceKeys=sourceKeys)
465 
466  # Fit for zeropoint
467  r = self.getZeroPoint(arrays.srcMag, arrays.refMag, arrays.magErr)
468  self.log.info("Magnitude zero point: %f +/- %f from %d stars", r.zp, r.sigma, r.ngood)
469 
470  # Prepare the results
471  flux0 = 10**(0.4*r.zp) # Flux of mag=0 star
472  flux0err = 0.4*math.log(10)*flux0*r.sigma # Error in flux0
473  calib = Calib()
474  calib.setFluxMag0(flux0, flux0err)
475 
476  return pipeBase.Struct(
477  calib=calib,
478  arrays=arrays,
479  matches=matches,
480  zp=r.zp,
481  sigma=r.sigma,
482  ngood=r.ngood,
483  )
484 
485  def displaySources(self, exposure, matches, reserved, frame=1):
486  """Display sources we'll use for photocal
487 
488  Sources that will be actually used will be green.
489  Sources reserved from the fit will be red.
490 
491  Parameters
492  ----------
493  exposure : `lsst.afw.image.ExposureF`
494  Exposure to display.
495  matches : `list` of `lsst.afw.table.RefMatch`
496  Matches used for photocal.
497  reserved : `numpy.ndarray` of type `bool`
498  Boolean array indicating sources that are reserved.
499  frame : `int`
500  Frame number for display.
501  """
502  ds9.mtv(exposure, frame=frame, title="photocal")
503  with ds9.Buffering():
504  for mm, rr in zip(matches, reserved):
505  x, y = mm.second.getCentroid()
506  ctype = ds9.RED if rr else ds9.GREEN
507  ds9.dot("o", x, y, size=4, frame=frame, ctype=ctype)
508 
509  def getZeroPoint(self, src, ref, srcErr=None, zp0=None):
510  """!Flux calibration code, returning (ZeroPoint, Distribution Width, Number of stars)
511 
512  We perform nIter iterations of a simple sigma-clipping algorithm with a couple of twists:
513  1. We use the median/interquartile range to estimate the position to clip around, and the
514  "sigma" to use.
515  2. We never allow sigma to go _above_ a critical value sigmaMax --- if we do, a sufficiently
516  large estimate will prevent the clipping from ever taking effect.
517  3. Rather than start with the median we start with a crude mode. This means that a set of magnitude
518  residuals with a tight core and asymmetrical outliers will start in the core. We use the width of
519  this core to set our maximum sigma (see 2.)
520 
521  @return Struct of:
522  - zp ---------- Photometric zero point (mag)
523  - sigma ------- Standard deviation of fit of zero point (mag)
524  - ngood ------- Number of sources used to fit zero point
525  """
526  sigmaMax = self.config.sigmaMax
527 
528  dmag = ref - src
529 
530  indArr = np.argsort(dmag)
531  dmag = dmag[indArr]
532 
533  if srcErr is not None:
534  dmagErr = srcErr[indArr]
535  else:
536  dmagErr = np.ones(len(dmag))
537 
538  # need to remove nan elements to avoid errors in stats calculation with numpy
539  ind_noNan = np.array([i for i in range(len(dmag))
540  if (not np.isnan(dmag[i]) and not np.isnan(dmagErr[i]))])
541  dmag = dmag[ind_noNan]
542  dmagErr = dmagErr[ind_noNan]
543 
544  IQ_TO_STDEV = 0.741301109252802 # 1 sigma in units of interquartile (assume Gaussian)
545 
546  npt = len(dmag)
547  ngood = npt
548  good = None # set at end of first iteration
549  for i in range(self.config.nIter):
550  if i > 0:
551  npt = sum(good)
552 
553  center = None
554  if i == 0:
555  #
556  # Start by finding the mode
557  #
558  nhist = 20
559  try:
560  hist, edges = np.histogram(dmag, nhist, new=True)
561  except TypeError:
562  hist, edges = np.histogram(dmag, nhist) # they removed new=True around numpy 1.5
563  imode = np.arange(nhist)[np.where(hist == hist.max())]
564 
565  if imode[-1] - imode[0] + 1 == len(imode): # Multiple modes, but all contiguous
566  if zp0:
567  center = zp0
568  else:
569  center = 0.5*(edges[imode[0]] + edges[imode[-1] + 1])
570 
571  peak = sum(hist[imode])/len(imode) # peak height
572 
573  # Estimate FWHM of mode
574  j = imode[0]
575  while j >= 0 and hist[j] > 0.5*peak:
576  j -= 1
577  j = max(j, 0)
578  q1 = dmag[sum(hist[range(j)])]
579 
580  j = imode[-1]
581  while j < nhist and hist[j] > 0.5*peak:
582  j += 1
583  j = min(j, nhist - 1)
584  j = min(sum(hist[range(j)]), npt - 1)
585  q3 = dmag[j]
586 
587  if q1 == q3:
588  q1 = dmag[int(0.25*npt)]
589  q3 = dmag[int(0.75*npt)]
590 
591  sig = (q3 - q1)/2.3 # estimate of standard deviation (based on FWHM; 2.358 for Gaussian)
592 
593  if sigmaMax is None:
594  sigmaMax = 2*sig # upper bound on st. dev. for clipping. multiplier is a heuristic
595 
596  self.log.debug("Photo calibration histogram: center = %.2f, sig = %.2f", center, sig)
597 
598  else:
599  if sigmaMax is None:
600  sigmaMax = dmag[-1] - dmag[0]
601 
602  center = np.median(dmag)
603  q1 = dmag[int(0.25*npt)]
604  q3 = dmag[int(0.75*npt)]
605  sig = (q3 - q1)/2.3 # estimate of standard deviation (based on FWHM; 2.358 for Gaussian)
606 
607  if center is None: # usually equivalent to (i > 0)
608  gdmag = dmag[good]
609  if self.config.useMedian:
610  center = np.median(gdmag)
611  else:
612  gdmagErr = dmagErr[good]
613  center = np.average(gdmag, weights=gdmagErr)
614 
615  q3 = gdmag[min(int(0.75*npt + 0.5), npt - 1)]
616  q1 = gdmag[min(int(0.25*npt + 0.5), npt - 1)]
617 
618  sig = IQ_TO_STDEV*(q3 - q1) # estimate of standard deviation
619 
620  good = abs(dmag - center) < self.config.nSigma*min(sig, sigmaMax) # don't clip too softly
621 
622  # =-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
623  if self.scatterPlot:
624  try:
625  self.fig.clf()
626 
627  axes = self.fig.add_axes((0.1, 0.1, 0.85, 0.80))
628 
629  axes.plot(ref[good], dmag[good] - center, "b+")
630  axes.errorbar(ref[good], dmag[good] - center, yerr=dmagErr[good],
631  linestyle='', color='b')
632 
633  bad = np.logical_not(good)
634  if len(ref[bad]) > 0:
635  axes.plot(ref[bad], dmag[bad] - center, "r+")
636  axes.errorbar(ref[bad], dmag[bad] - center, yerr=dmagErr[bad],
637  linestyle='', color='r')
638 
639  axes.plot((-100, 100), (0, 0), "g-")
640  for x in (-1, 1):
641  axes.plot((-100, 100), x*0.05*np.ones(2), "g--")
642 
643  axes.set_ylim(-1.1, 1.1)
644  axes.set_xlim(24, 13)
645  axes.set_xlabel("Reference")
646  axes.set_ylabel("Reference - Instrumental")
647 
648  self.fig.show()
649 
650  if self.scatterPlot > 1:
651  reply = None
652  while i == 0 or reply != "c":
653  try:
654  reply = input("Next iteration? [ynhpc] ")
655  except EOFError:
656  reply = "n"
657 
658  if reply == "h":
659  print("Options: c[ontinue] h[elp] n[o] p[db] y[es]", file=sys.stderr)
660  continue
661 
662  if reply in ("", "c", "n", "p", "y"):
663  break
664  else:
665  print("Unrecognised response: %s" % reply, file=sys.stderr)
666 
667  if reply == "n":
668  break
669  elif reply == "p":
670  import pdb
671  pdb.set_trace()
672  except Exception as e:
673  print("Error plotting in PhotoCal.getZeroPoint: %s" % e, file=sys.stderr)
674 
675  # =-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
676 
677  old_ngood = ngood
678  ngood = sum(good)
679  if ngood == 0:
680  msg = "PhotoCal.getZeroPoint: no good stars remain"
681 
682  if i == 0: # failed the first time round -- probably all fell in one bin
683  center = np.average(dmag, weights=dmagErr)
684  msg += " on first iteration; using average of all calibration stars"
685 
686  self.log.warn(msg)
687 
688  return pipeBase.Struct(
689  zp=center,
690  sigma=sig,
691  ngood=len(dmag))
692  elif ngood == old_ngood:
693  break
694 
695  if False:
696  ref = ref[good]
697  dmag = dmag[good]
698  dmagErr = dmagErr[good]
699 
700  dmag = dmag[good]
701  dmagErr = dmagErr[good]
702  zp, weightSum = np.average(dmag, weights=1/dmagErr**2, returned=True)
703  sigma = np.sqrt(1.0/weightSum)
704  return pipeBase.Struct(
705  zp=zp,
706  sigma=sigma,
707  ngood=len(dmag),
708  )
def __init__(self, refObjLoader, schema=None, kwds)
Create the photometric calibration task.
Definition: photoCal.py:245
def run(self, exposure, sourceCat, expId=0)
Do photometric calibration - select matches to use and (possibly iteratively) compute the zero point...
Definition: photoCal.py:390
def getSourceKeys(self, schema)
Return a struct containing the source catalog keys for fields used by PhotoCalTask.
Definition: photoCal.py:261
def extractMagArrays(self, matches, filterName, sourceKeys)
Extract magnitude and magnitude error arrays from the given matches.
Definition: photoCal.py:273
Calculate the zero point of an exposure given a lsst.afw.table.ReferenceMatchVector.
Definition: photoCal.py:124
def getZeroPoint(self, src, ref, srcErr=None, zp0=None)
Flux calibration code, returning (ZeroPoint, Distribution Width, Number of stars) ...
Definition: photoCal.py:509
def displaySources(self, exposure, matches, reserved, frame=1)
Definition: photoCal.py:485