lsst.pipe.tasks g0099ee1360+6048f86b6d
photoCal.py
Go to the documentation of this file.
1# This file is part of pipe_tasks.
2#
3# Developed for the LSST Data Management System.
4# This product includes software developed by the LSST Project
5# (https://www.lsst.org).
6# See the COPYRIGHT file at the top-level directory of this distribution
7# for details of code ownership.
8#
9# This program is free software: you can redistribute it and/or modify
10# it under the terms of the GNU General Public License as published by
11# the Free Software Foundation, either version 3 of the License, or
12# (at your option) any later version.
13#
14# This program is distributed in the hope that it will be useful,
15# but WITHOUT ANY WARRANTY; without even the implied warranty of
16# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17# GNU General Public License for more details.
18#
19# You should have received a copy of the GNU General Public License
20# along with this program. If not, see <https://www.gnu.org/licenses/>.
21#
22# @package lsst.pipe.tasks.
23import math
24import sys
25
26import numpy as np
27import astropy.units as u
28
29import lsst.pex.config as pexConf
30import lsst.pipe.base as pipeBase
31from lsst.afw.image import abMagErrFromFluxErr, makePhotoCalibFromCalibZeroPoint
32import lsst.afw.table as afwTable
33from lsst.meas.astrom import DirectMatchTask, DirectMatchConfigWithoutLoader
34import lsst.afw.display as afwDisplay
35from lsst.meas.algorithms import getRefFluxField, ReserveSourcesTask
36from .colorterms import ColortermLibrary
37
38__all__ = ["PhotoCalTask", "PhotoCalConfig"]
39
40
41class 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.applyColorTermsapplyColorTerms and self.photoCatNamephotoCatName is None:
104 raise RuntimeError("applyColorTerms=True requires photoCatName is non-None")
105 if self.applyColorTermsapplyColorTerms and len(self.colortermscolorterms.data) == 0:
106 raise RuntimeError("applyColorTerms=True requires colorterms be provided")
107
108 def setDefaults(self):
109 pexConf.Config.setDefaults(self)
110 self.matchmatch.sourceSelection.doFlags = True
111 self.matchmatch.sourceSelection.flags.bad = [
112 "base_PixelFlags_flag_edge",
113 "base_PixelFlags_flag_interpolated",
114 "base_PixelFlags_flag_saturated",
115 ]
116 self.matchmatch.sourceSelection.doUnresolved = True
117
118
119
125
126class 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
145Calculate an Exposure's zero-point given a set of flux measurements of stars matched to an input catalogue.
146The type of flux to use is specified by PhotoCalConfig.fluxField.
147
148The 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
151these 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
163See @ref PhotoCalConfig
164
165@section pipe_tasks_photocal_Debug Debug variables
166
167The @link lsst.pipe.base.cmdLineTask.CmdLineTask command line task@endlink interface supports a
168flag @c -d to import @b debug.py from your @c PYTHONPATH; see @ref baseDebug for more about @b debug.py files.
169
170The 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 the display'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
191This code is in @link examples/photoCalTask.py@endlink, and can be run as @em e.g.
192@code
193examples/photoCalTask.py
194@endcode
195@dontinclude photoCalTask.py
196
197Import the tasks (there are some other standard imports; read the file for details)
198@skipline from lsst.pipe.tasks.astrometry
199@skipline measPhotocal
200
201We need to create both our tasks before processing any data as the task constructors
202can add extra columns to the schema which we get from the input catalogue, @c scrCat:
203@skipline getSchema
204
205Astrometry 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,
209so we tell it to use the @c r band)
210
211Then photometry:
212@skip measPhotocal
213@until pTask
214
215If 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
220We're now ready to process the data (we could loop over multiple exposures/catalogues using the same
221task objects):
222@skip matches
223@until result
224
225We can then unpack and use the results:
226@skip calib
227@until np.log
228
229<HR>
230To 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
242into 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.scatterPlotscatterPlot = None
252 self.figfig = None
253 if schema is not None:
254 self.usedKeyusedKey = schema.addField("calib_photometry_used", type="Flag",
255 doc="set if source was used in photometric calibration")
256 else:
257 self.usedKeyusedKey = None
258 self.matchmatch = 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, filterLabel, 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] filterLabel Label 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.warning("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 nJy
310 referenceFlux = (0*u.ABmag).to_value(u.nJy)
311 srcInstFluxArr = srcInstFluxArr * referenceFlux
312 srcInstFluxErrArr = srcInstFluxErrArr * referenceFlux
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 filter=%r, config.photoCatName=%s because %s",
330 filterLabel.physicalLabel, self.config.photoCatName, applyCTReason)
331 colorterm = self.config.colorterms.getColorterm(filterLabel.physicalLabel,
332 self.config.photoCatName,
333 doRaise=True)
334 refCat = afwTable.SimpleCatalog(matches[0].first.schema)
335
336 # extract the matched refCat as a Catalog for the colorterm code
337 refCat.reserve(len(matches))
338 for x in matches:
339 record = refCat.addNew()
340 record.assign(x.first)
341
342 refMagArr, refMagErrArr = colorterm.getCorrectedMagnitudes(refCat)
343 fluxFieldList = [getRefFluxField(refSchema, filt) for filt in (colorterm.primary,
344 colorterm.secondary)]
345 else:
346 # no colorterms to apply
347 self.log.info("Not applying color terms because %s", applyCTReason)
348 colorterm = None
349
350 fluxFieldList = [getRefFluxField(refSchema, filterLabel.bandLabel)]
351 fluxField = getRefFluxField(refSchema, filterLabel.bandLabel)
352 fluxKey = refSchema.find(fluxField).key
353 refFluxArr = np.array([m.first.get(fluxKey) for m in matches])
354
355 try:
356 fluxErrKey = refSchema.find(fluxField + "Err").key
357 refFluxErrArr = np.array([m.first.get(fluxErrKey) for m in matches])
358 except KeyError:
359 # Reference catalogue may not have flux uncertainties; HACK DM-2308
360 self.log.warning("Reference catalog does not have flux uncertainties for %s;"
361 " using sqrt(flux).", fluxField)
362 refFluxErrArr = np.sqrt(refFluxArr)
363
364 refMagArr = u.Quantity(refFluxArr, u.nJy).to_value(u.ABmag)
365 # HACK convert to Jy until we have a replacement for this (DM-16903)
366 refMagErrArr = abMagErrFromFluxErr(refFluxErrArr*1e-9, refFluxArr*1e-9)
367
368 # compute the source catalog magnitudes and errors
369 srcMagArr = u.Quantity(srcInstFluxArr, u.nJy).to_value(u.ABmag)
370 # Fitting with error bars in both axes is hard
371 # for now ignore reference flux error, but ticket DM-2308 is a request for a better solution
372 # HACK convert to Jy until we have a replacement for this (DM-16903)
373 magErrArr = abMagErrFromFluxErr(srcInstFluxErrArr*1e-9, srcInstFluxArr*1e-9)
374 if self.config.magErrFloor != 0.0:
375 magErrArr = (magErrArr**2 + self.config.magErrFloor**2)**0.5
376
377 srcMagErrArr = abMagErrFromFluxErr(srcInstFluxErrArr*1e-9, srcInstFluxArr*1e-9)
378
379 good = np.isfinite(srcMagArr) & np.isfinite(refMagArr)
380
381 return pipeBase.Struct(
382 srcMag=srcMagArr[good],
383 refMag=refMagArr[good],
384 magErr=magErrArr[good],
385 srcMagErr=srcMagErrArr[good],
386 refMagErr=refMagErrArr[good],
387 refFluxFieldList=fluxFieldList,
388 )
389
390 @pipeBase.timeMethod
391 def run(self, exposure, sourceCat, expId=0):
392 """!Do photometric calibration - select matches to use and (possibly iteratively) compute
393 the zero point.
394
395 @param[in] exposure Exposure upon which the sources in the matches were detected.
396 @param[in] sourceCat A catalog of sources to use in the calibration
397 (@em i.e. a list of lsst.afw.table.Match with
398 @c first being of type lsst.afw.table.SimpleRecord and @c second type lsst.afw.table.SourceRecord ---
399 the reference object and matched object respectively).
400 (will not be modified except to set the outputField if requested.).
401
402 @return Struct of:
403 - photoCalib -- @link lsst::afw::image::PhotoCalib@endlink object containing the calibration
404 - arrays ------ Magnitude arrays returned be PhotoCalTask.extractMagArrays
405 - matches ----- Final ReferenceMatchVector, as returned by PhotoCalTask.selectMatches.
406 - zp ---------- Photometric zero point (mag)
407 - sigma ------- Standard deviation of fit of photometric zero point (mag)
408 - ngood ------- Number of sources used to fit photometric zero point
409
410 The exposure is only used to provide the name of the filter being calibrated (it may also be
411 used to generate debugging plots).
412
413 The reference objects:
414 - Must include a field @c photometric; True for objects which should be considered as
415 photometric standards
416 - Must include a field @c flux; the flux used to impose a magnitude limit and also to calibrate
417 the data to (unless a color term is specified, in which case ColorTerm.primary is used;
418 See https://jira.lsstcorp.org/browse/DM-933)
419 - May include a field @c stargal; if present, True means that the object is a star
420 - May include a field @c var; if present, True means that the object is variable
421
422 The measured sources:
423 - Must include PhotoCalConfig.fluxField; the flux measurement to be used for calibration
424
425 @throws RuntimeError with the following strings:
426
427 <DL>
428 <DT> No matches to use for photocal
429 <DD> No matches are available (perhaps no sources/references were selected by the matcher).
430 <DT> No reference stars are available
431 <DD> No matches are available from which to extract magnitudes.
432 </DL>
433 """
434 import lsstDebug
435
436 display = lsstDebug.Info(__name__).display
437 displaySources = display and lsstDebug.Info(__name__).displaySources
438 self.scatterPlotscatterPlot = display and lsstDebug.Info(__name__).scatterPlot
439
440 if self.scatterPlotscatterPlot:
441 from matplotlib import pyplot
442 try:
443 self.figfig.clf()
444 except Exception:
445 self.figfig = pyplot.figure()
446
447 filterLabel = exposure.getFilterLabel()
448
449 # Match sources
450 matchResults = self.matchmatch.run(sourceCat, filterLabel.bandLabel)
451 matches = matchResults.matches
452
453 reserveResults = self.reserve.run([mm.second for mm in matches], expId=expId)
454 if displaySources:
455 self.displaySourcesdisplaySources(exposure, matches, reserveResults.reserved)
456 if reserveResults.reserved.sum() > 0:
457 matches = [mm for mm, use in zip(matches, reserveResults.use) if use]
458 if len(matches) == 0:
459 raise RuntimeError("No matches to use for photocal")
460 if self.usedKeyusedKey is not None:
461 for mm in matches:
462 mm.second.set(self.usedKeyusedKey, True)
463
464 # Prepare for fitting
465 sourceKeys = self.getSourceKeysgetSourceKeys(matches[0].second.schema)
466 arrays = self.extractMagArraysextractMagArrays(matches, filterLabel, sourceKeys)
467
468 # Fit for zeropoint
469 r = self.getZeroPointgetZeroPoint(arrays.srcMag, arrays.refMag, arrays.magErr)
470 self.log.info("Magnitude zero point: %f +/- %f from %d stars", r.zp, r.sigma, r.ngood)
471
472 # Prepare the results
473 flux0 = 10**(0.4*r.zp) # Flux of mag=0 star
474 flux0err = 0.4*math.log(10)*flux0*r.sigma # Error in flux0
475 photoCalib = makePhotoCalibFromCalibZeroPoint(flux0, flux0err)
476
477 return pipeBase.Struct(
478 photoCalib=photoCalib,
479 arrays=arrays,
480 matches=matches,
481 zp=r.zp,
482 sigma=r.sigma,
483 ngood=r.ngood,
484 )
485
486 def displaySources(self, exposure, matches, reserved, frame=1):
487 """Display sources we'll use for photocal
488
489 Sources that will be actually used will be green.
490 Sources reserved from the fit will be red.
491
492 Parameters
493 ----------
494 exposure : `lsst.afw.image.ExposureF`
495 Exposure to display.
496 matches : `list` of `lsst.afw.table.RefMatch`
497 Matches used for photocal.
498 reserved : `numpy.ndarray` of type `bool`
499 Boolean array indicating sources that are reserved.
500 frame : `int`
501 Frame number for display.
502 """
503 disp = afwDisplay.getDisplay(frame=frame)
504 disp.mtv(exposure, title="photocal")
505 with disp.Buffering():
506 for mm, rr in zip(matches, reserved):
507 x, y = mm.second.getCentroid()
508 ctype = afwDisplay.RED if rr else afwDisplay.GREEN
509 disp.dot("o", x, y, size=4, ctype=ctype)
510
511 def getZeroPoint(self, src, ref, srcErr=None, zp0=None):
512 """!Flux calibration code, returning (ZeroPoint, Distribution Width, Number of stars)
513
514 We perform nIter iterations of a simple sigma-clipping algorithm with a couple of twists:
515 1. We use the median/interquartile range to estimate the position to clip around, and the
516 "sigma" to use.
517 2. We never allow sigma to go _above_ a critical value sigmaMax --- if we do, a sufficiently
518 large estimate will prevent the clipping from ever taking effect.
519 3. Rather than start with the median we start with a crude mode. This means that a set of magnitude
520 residuals with a tight core and asymmetrical outliers will start in the core. We use the width of
521 this core to set our maximum sigma (see 2.)
522
523 @return Struct of:
524 - zp ---------- Photometric zero point (mag)
525 - sigma ------- Standard deviation of fit of zero point (mag)
526 - ngood ------- Number of sources used to fit zero point
527 """
528 sigmaMax = self.config.sigmaMax
529
530 dmag = ref - src
531
532 indArr = np.argsort(dmag)
533 dmag = dmag[indArr]
534
535 if srcErr is not None:
536 dmagErr = srcErr[indArr]
537 else:
538 dmagErr = np.ones(len(dmag))
539
540 # need to remove nan elements to avoid errors in stats calculation with numpy
541 ind_noNan = np.array([i for i in range(len(dmag))
542 if (not np.isnan(dmag[i]) and not np.isnan(dmagErr[i]))])
543 dmag = dmag[ind_noNan]
544 dmagErr = dmagErr[ind_noNan]
545
546 IQ_TO_STDEV = 0.741301109252802 # 1 sigma in units of interquartile (assume Gaussian)
547
548 npt = len(dmag)
549 ngood = npt
550 good = None # set at end of first iteration
551 for i in range(self.config.nIter):
552 if i > 0:
553 npt = sum(good)
554
555 center = None
556 if i == 0:
557 #
558 # Start by finding the mode
559 #
560 nhist = 20
561 try:
562 hist, edges = np.histogram(dmag, nhist, new=True)
563 except TypeError:
564 hist, edges = np.histogram(dmag, nhist) # they removed new=True around numpy 1.5
565 imode = np.arange(nhist)[np.where(hist == hist.max())]
566
567 if imode[-1] - imode[0] + 1 == len(imode): # Multiple modes, but all contiguous
568 if zp0:
569 center = zp0
570 else:
571 center = 0.5*(edges[imode[0]] + edges[imode[-1] + 1])
572
573 peak = sum(hist[imode])/len(imode) # peak height
574
575 # Estimate FWHM of mode
576 j = imode[0]
577 while j >= 0 and hist[j] > 0.5*peak:
578 j -= 1
579 j = max(j, 0)
580 q1 = dmag[sum(hist[range(j)])]
581
582 j = imode[-1]
583 while j < nhist and hist[j] > 0.5*peak:
584 j += 1
585 j = min(j, nhist - 1)
586 j = min(sum(hist[range(j)]), npt - 1)
587 q3 = dmag[j]
588
589 if q1 == q3:
590 q1 = dmag[int(0.25*npt)]
591 q3 = dmag[int(0.75*npt)]
592
593 sig = (q3 - q1)/2.3 # estimate of standard deviation (based on FWHM; 2.358 for Gaussian)
594
595 if sigmaMax is None:
596 sigmaMax = 2*sig # upper bound on st. dev. for clipping. multiplier is a heuristic
597
598 self.log.debug("Photo calibration histogram: center = %.2f, sig = %.2f", center, sig)
599
600 else:
601 if sigmaMax is None:
602 sigmaMax = dmag[-1] - dmag[0]
603
604 center = np.median(dmag)
605 q1 = dmag[int(0.25*npt)]
606 q3 = dmag[int(0.75*npt)]
607 sig = (q3 - q1)/2.3 # estimate of standard deviation (based on FWHM; 2.358 for Gaussian)
608
609 if center is None: # usually equivalent to (i > 0)
610 gdmag = dmag[good]
611 if self.config.useMedian:
612 center = np.median(gdmag)
613 else:
614 gdmagErr = dmagErr[good]
615 center = np.average(gdmag, weights=gdmagErr)
616
617 q3 = gdmag[min(int(0.75*npt + 0.5), npt - 1)]
618 q1 = gdmag[min(int(0.25*npt + 0.5), npt - 1)]
619
620 sig = IQ_TO_STDEV*(q3 - q1) # estimate of standard deviation
621
622 good = abs(dmag - center) < self.config.nSigma*min(sig, sigmaMax) # don't clip too softly
623
624 # =-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
625 if self.scatterPlotscatterPlot:
626 try:
627 self.figfig.clf()
628
629 axes = self.figfig.add_axes((0.1, 0.1, 0.85, 0.80))
630
631 axes.plot(ref[good], dmag[good] - center, "b+")
632 axes.errorbar(ref[good], dmag[good] - center, yerr=dmagErr[good],
633 linestyle='', color='b')
634
635 bad = np.logical_not(good)
636 if len(ref[bad]) > 0:
637 axes.plot(ref[bad], dmag[bad] - center, "r+")
638 axes.errorbar(ref[bad], dmag[bad] - center, yerr=dmagErr[bad],
639 linestyle='', color='r')
640
641 axes.plot((-100, 100), (0, 0), "g-")
642 for x in (-1, 1):
643 axes.plot((-100, 100), x*0.05*np.ones(2), "g--")
644
645 axes.set_ylim(-1.1, 1.1)
646 axes.set_xlim(24, 13)
647 axes.set_xlabel("Reference")
648 axes.set_ylabel("Reference - Instrumental")
649
650 self.figfig.show()
651
652 if self.scatterPlotscatterPlot > 1:
653 reply = None
654 while i == 0 or reply != "c":
655 try:
656 reply = input("Next iteration? [ynhpc] ")
657 except EOFError:
658 reply = "n"
659
660 if reply == "h":
661 print("Options: c[ontinue] h[elp] n[o] p[db] y[es]", file=sys.stderr)
662 continue
663
664 if reply in ("", "c", "n", "p", "y"):
665 break
666 else:
667 print("Unrecognised response: %s" % reply, file=sys.stderr)
668
669 if reply == "n":
670 break
671 elif reply == "p":
672 import pdb
673 pdb.set_trace()
674 except Exception as e:
675 print("Error plotting in PhotoCal.getZeroPoint: %s" % e, file=sys.stderr)
676
677 # =-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
678
679 old_ngood = ngood
680 ngood = sum(good)
681 if ngood == 0:
682 msg = "PhotoCal.getZeroPoint: no good stars remain"
683
684 if i == 0: # failed the first time round -- probably all fell in one bin
685 center = np.average(dmag, weights=dmagErr)
686 msg += " on first iteration; using average of all calibration stars"
687
688 self.log.warning(msg)
689
690 return pipeBase.Struct(
691 zp=center,
692 sigma=sig,
693 ngood=len(dmag))
694 elif ngood == old_ngood:
695 break
696
697 if False:
698 ref = ref[good]
699 dmag = dmag[good]
700 dmagErr = dmagErr[good]
701
702 dmag = dmag[good]
703 dmagErr = dmagErr[good]
704 zp, weightSum = np.average(dmag, weights=1/dmagErr**2, returned=True)
705 sigma = np.sqrt(1.0/weightSum)
706 return pipeBase.Struct(
707 zp=zp,
708 sigma=sigma,
709 ngood=len(dmag),
710 )
def getZeroPoint(self, src, ref, srcErr=None, zp0=None)
Flux calibration code, returning (ZeroPoint, Distribution Width, Number of stars)
Definition: photoCal.py:511
def extractMagArrays(self, matches, filterLabel, sourceKeys)
Extract magnitude and magnitude error arrays from the given matches.
Definition: photoCal.py:286
def __init__(self, refObjLoader, schema=None, **kwds)
Create the photometric calibration task.
Definition: photoCal.py:247
def displaySources(self, exposure, matches, reserved, frame=1)
Definition: photoCal.py:486
def run(self, exposure, sourceCat, expId=0)
Do photometric calibration - select matches to use and (possibly iteratively) compute the zero point.
Definition: photoCal.py:391