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