22 __all__ = [
"PcaPsfDeterminerConfig",
"PcaPsfDeterminerTask"]
36 from .psfDeterminer
import BasePsfDeterminerTask, psfDeterminerRegistry
37 from .psfCandidate
import PsfCandidateF
38 from .spatialModelPsf
import createKernelFromPsfCandidates, countPsfCandidates, \
39 fitSpatialKernelFromPsfCandidates, fitKernelParamsToImage
40 from .pcaPsf
import PcaPsf
45 """Return the number of PSF candidates to be rejected.
47 The number of candidates being rejected on each iteration gradually
48 increases, so that on the Nth of M iterations we reject N/M of the bad
53 numBadCandidates : `int`
54 Number of bad candidates under consideration.
57 The number of the current PSF iteration.
60 The total number of PSF iterations.
65 Number of candidates to reject.
67 return int(numBadCandidates*(numIter + 1)//totalIter + 0.5)
71 nonLinearSpatialFit = pexConfig.Field(
72 doc=
"Use non-linear fitter for spatial variation of Kernel",
76 nEigenComponents = pexConfig.Field(
77 doc=
"number of eigen components for PSF kernel creation",
81 spatialOrder = pexConfig.Field(
82 doc=
"specify spatial order for PSF kernel creation",
86 sizeCellX = pexConfig.Field(
87 doc=
"size of cell used to determine PSF (pixels, column direction)",
91 check=
lambda x: x >= 10,
93 sizeCellY = pexConfig.Field(
94 doc=
"size of cell used to determine PSF (pixels, row direction)",
96 default=sizeCellX.default,
98 check=
lambda x: x >= 10,
100 nStarPerCell = pexConfig.Field(
101 doc=
"number of stars per psf cell for PSF kernel creation",
105 borderWidth = pexConfig.Field(
106 doc=
"Number of pixels to ignore around the edge of PSF candidate postage stamps",
110 nStarPerCellSpatialFit = pexConfig.Field(
111 doc=
"number of stars per psf Cell for spatial fitting",
115 constantWeight = pexConfig.Field(
116 doc=
"Should each PSF candidate be given the same weight, independent of magnitude?",
120 nIterForPsf = pexConfig.Field(
121 doc=
"number of iterations of PSF candidate star list",
125 tolerance = pexConfig.Field(
126 doc=
"tolerance of spatial fitting",
130 lam = pexConfig.Field(
131 doc=
"floor for variance is lam*data",
135 reducedChi2ForPsfCandidates = pexConfig.Field(
136 doc=
"for psf candidate evaluation",
140 spatialReject = pexConfig.Field(
141 doc=
"Rejection threshold (stdev) for candidates based on spatial fit",
145 pixelThreshold = pexConfig.Field(
146 doc=
"Threshold (stdev) for rejecting extraneous pixels around candidate; applied if positive",
150 doRejectBlends = pexConfig.Field(
151 doc=
"Reject candidates that are blended?",
155 doMaskBlends = pexConfig.Field(
156 doc=
"Mask blends in image?",
163 """A measurePsfTask psf estimator.
165 ConfigClass = PcaPsfDeterminerConfig
167 def _fitPsf(self, exposure, psfCellSet, kernelSize, nEigenComponents):
168 PsfCandidateF.setPixelThreshold(self.config.pixelThreshold)
169 PsfCandidateF.setMaskBlends(self.config.doMaskBlends)
173 for nEigen
in range(nEigenComponents, 0, -1):
177 psfCellSet, exposure.getDimensions(), exposure.getXY0(), nEigen,
178 self.config.spatialOrder, kernelSize, self.config.nStarPerCell,
179 bool(self.config.constantWeight))
184 raise IndexError(
"No viable PSF candidates survive")
186 self.log.
warn(
"%s: reducing number of eigen components" % e.what())
191 size = kernelSize + 2*self.config.borderWidth
194 for val
in eigenValues]
198 kernel, psfCellSet, bool(self.config.nonLinearSpatialFit),
199 self.config.nStarPerCellSpatialFit, self.config.tolerance, self.config.lam)
203 return psf, eigenValues, nEigen, chi2
205 def determinePsf(self, exposure, psfCandidateList, metadata=None, flagKey=None):
206 """Determine a PCA PSF model for an exposure given a list of PSF candidates.
210 exposure : `lsst.afw.image.Exposure`
211 Exposure containing the psf candidates.
212 psfCandidateList : `list` of `lsst.meas.algorithms.PsfCandidate`
213 A sequence of PSF candidates typically obtained by detecting sources
214 and then running them through a star selector.
215 metadata : `lsst.daf.base import PropertyList` or `None`, optional
216 A home for interesting tidbits of information.
217 flagKey : `str`, optional
218 Schema key used to mark sources actually used in PSF determination.
222 psf : `lsst.meas.algorithms.PcaPsf`
224 psfCellSet : `lsst.afw.math.SpatialCellSet`
230 displayPsfCandidates =
lsstDebug.Info(__name__).displayPsfCandidates
232 displayPsfComponents =
lsstDebug.Info(__name__).displayPsfComponents
236 matchKernelAmplitudes =
lsstDebug.Info(__name__).matchKernelAmplitudes
238 keepMatplotlibPlots =
lsstDebug.Info(__name__).keepMatplotlibPlots
239 displayPsfSpatialModel =
lsstDebug.Info(__name__).displayPsfSpatialModel
246 afwDisplay.setDefaultMaskTransparency(75)
250 mi = exposure.getMaskedImage()
252 if len(psfCandidateList) == 0:
253 raise RuntimeError(
"No PSF candidates supplied.")
259 for i, psfCandidate
in enumerate(psfCandidateList):
260 if psfCandidate.getSource().getPsfFluxFlag():
264 psfCellSet.insertCandidate(psfCandidate)
265 except Exception
as e:
266 self.log.
debug(
"Skipping PSF candidate %d of %d: %s", i, len(psfCandidateList), e)
268 source = psfCandidate.getSource()
271 axes = afwEll.Axes(quad)
272 sizes.append(axes.getA())
274 raise RuntimeError(
"No usable PSF candidates supplied")
275 nEigenComponents = self.config.nEigenComponents
277 if self.config.kernelSize >= 15:
278 self.log.
warn(
"WARNING: NOT scaling kernelSize by stellar quadrupole moment "
279 "because config.kernelSize=%s >= 15; "
280 "using config.kernelSize as as the width, instead",
281 self.config.kernelSize)
282 actualKernelSize = int(self.config.kernelSize)
284 medSize = numpy.median(sizes)
285 actualKernelSize = 2*int(self.config.kernelSize*math.sqrt(medSize) + 0.5) + 1
286 if actualKernelSize < self.config.kernelSizeMin:
287 actualKernelSize = self.config.kernelSizeMin
288 if actualKernelSize > self.config.kernelSizeMax:
289 actualKernelSize = self.config.kernelSizeMax
292 print(
"Median size=%s" % (medSize,))
293 self.log.
trace(
"Kernel size=%s", actualKernelSize)
296 psfCandidateList[0].setHeight(actualKernelSize)
297 psfCandidateList[0].setWidth(actualKernelSize)
299 if self.config.doRejectBlends:
301 blendedCandidates = []
303 if len(cand.getSource().getFootprint().getPeaks()) > 1:
304 blendedCandidates.append((cell, cand))
307 print(
"Removing %d blended Psf candidates" % len(blendedCandidates))
308 for cell, cand
in blendedCandidates:
309 cell.removeCandidate(cand)
311 raise RuntimeError(
"All PSF candidates removed as blends")
315 disp = afwDisplay.Display(frame=0)
316 disp.mtv(exposure, title=
"psf determination")
317 utils.showPsfSpatialCells(exposure, psfCellSet, self.config.nStarPerCell, symb=
"o",
318 ctype=afwDisplay.CYAN, ctypeUnused=afwDisplay.YELLOW,
319 size=4, display=disp)
325 for iterNum
in range(self.config.nIterForPsf):
326 if display
and displayPsfCandidates:
329 for cell
in psfCellSet.getCellList():
330 for cand
in cell.begin(
not showBadCandidates):
332 im = cand.getMaskedImage()
334 chi2 = cand.getChi2()
338 stamps.append((im,
"%d%s" %
339 (utils.splitId(cand.getSource().getId(),
True)[
"objId"], chi2),
345 print(
"WARNING: No PSF candidates to show; try setting showBadCandidates=True")
347 mos = afwDisplay.utils.Mosaic()
348 for im, label, status
in stamps:
349 im =
type(im)(im,
True)
352 except NotImplementedError:
355 mos.append(im, label,
356 (afwDisplay.GREEN
if status == afwMath.SpatialCellCandidate.GOOD
else
357 afwDisplay.YELLOW
if status == afwMath.SpatialCellCandidate.UNKNOWN
else
360 disp8 = afwDisplay.Display(frame=8)
361 mos.makeMosaic(display=disp8, title=
"Psf Candidates")
370 psf, eigenValues, nEigenComponents, fitChi2 = \
371 self.
_fitPsf_fitPsf(exposure, psfCellSet, actualKernelSize, nEigenComponents)
376 for cell
in psfCellSet.getCellList():
378 for cand
in cell.begin(
False):
379 cand.setStatus(afwMath.SpatialCellCandidate.UNKNOWN)
380 rchi2 = cand.getChi2()
381 if not numpy.isfinite(rchi2)
or rchi2 <= 0:
383 awfulCandidates.append(cand)
385 self.log.
debug(
"chi^2=%s; id=%s",
386 cand.getChi2(), cand.getSource().getId())
387 for cand
in awfulCandidates:
389 print(
"Removing bad candidate: id=%d, chi^2=%f" %
390 (cand.getSource().getId(), cand.getChi2()))
391 cell.removeCandidate(cand)
396 badCandidates =
list()
397 for cell
in psfCellSet.getCellList():
398 for cand
in cell.begin(
False):
399 rchi2 = cand.getChi2()
401 if rchi2 > self.config.reducedChi2ForPsfCandidates:
402 badCandidates.append(cand)
404 badCandidates.sort(key=
lambda x: x.getChi2(), reverse=
True)
406 self.config.nIterForPsf)
407 for i, c
in zip(range(numBad), badCandidates):
413 print(
"Chi^2 clipping %-4d %.2g" % (c.getSource().getId(), chi2))
414 c.setStatus(afwMath.SpatialCellCandidate.BAD)
427 kernel = psf.getKernel()
428 noSpatialKernel = psf.getKernel()
429 for cell
in psfCellSet.getCellList():
430 for cand
in cell.begin(
False):
433 im = cand.getMaskedImage(kernel.getWidth(), kernel.getHeight())
441 for p, k
in zip(params, kernels):
444 predict = [kernel.getSpatialFunction(k)(candCenter.getX(), candCenter.getY())
for
445 k
in range(kernel.getNKernelParameters())]
447 residuals.append([a/amp - p
for a, p
in zip(params, predict)])
448 candidates.append(cand)
450 residuals = numpy.array(residuals)
452 for k
in range(kernel.getNKernelParameters()):
455 mean = residuals[:, k].mean()
456 rms = residuals[:, k].
std()
459 sr = numpy.sort(residuals[:, k])
460 mean = (sr[int(0.5*len(sr))]
if len(sr)%2
else
461 0.5*(sr[int(0.5*len(sr))] + sr[int(0.5*len(sr)) + 1]))
462 rms = 0.74*(sr[int(0.75*len(sr))] - sr[int(0.25*len(sr))])
465 mean = stats.getValue(afwMath.MEANCLIP)
466 rms = stats.getValue(afwMath.STDEVCLIP)
468 rms =
max(1.0e-4, rms)
471 print(
"Mean for component %d is %f" % (k, mean))
472 print(
"RMS for component %d is %f" % (k, rms))
473 badCandidates =
list()
474 for i, cand
in enumerate(candidates):
475 if numpy.fabs(residuals[i, k] - mean) > self.config.spatialReject*rms:
476 badCandidates.append(i)
478 badCandidates.sort(key=
lambda x: numpy.fabs(residuals[x, k] - mean), reverse=
True)
481 self.config.nIterForPsf)
483 for i, c
in zip(range(
min(len(badCandidates), numBad)), badCandidates):
486 print(
"Spatial clipping %d (%f,%f) based on %d: %f vs %f" %
487 (cand.getSource().getId(), cand.getXCenter(), cand.getYCenter(), k,
488 residuals[badCandidates[i], k], self.config.spatialReject*rms))
489 cand.setStatus(afwMath.SpatialCellCandidate.BAD)
494 if display
and displayIterations:
498 utils.showPsfSpatialCells(exposure, psfCellSet, self.config.nStarPerCell, showChi2=
True,
499 symb=
"o", size=8, display=disp, ctype=afwDisplay.YELLOW,
500 ctypeBad=afwDisplay.RED, ctypeUnused=afwDisplay.MAGENTA)
501 if self.config.nStarPerCellSpatialFit != self.config.nStarPerCell:
502 utils.showPsfSpatialCells(exposure, psfCellSet, self.config.nStarPerCellSpatialFit,
503 symb=
"o", size=10, display=disp,
504 ctype=afwDisplay.YELLOW, ctypeBad=afwDisplay.RED)
508 disp4 = afwDisplay.Display(frame=4)
509 utils.showPsfCandidates(exposure, psfCellSet, psf=psf, display=disp4,
510 normalize=normalizeResiduals,
511 showBadCandidates=showBadCandidates)
512 disp5 = afwDisplay.Display(frame=5)
513 utils.showPsfCandidates(exposure, psfCellSet, psf=psf, display=disp5,
514 normalize=normalizeResiduals,
515 showBadCandidates=showBadCandidates,
518 if not showBadCandidates:
519 showBadCandidates =
True
523 if displayPsfComponents:
524 disp6 = afwDisplay.Display(frame=6)
525 utils.showPsf(psf, eigenValues, display=disp6)
527 disp7 = afwDisplay.Display(frame=7)
528 utils.showPsfMosaic(exposure, psf, display=disp7, showFwhm=
True)
529 disp7.scale(
'linear', 0, 1)
530 if displayPsfSpatialModel:
531 utils.plotPsfSpatialModel(exposure, psf, psfCellSet, showBadCandidates=
True,
532 matchKernelAmplitudes=matchKernelAmplitudes,
533 keepPlots=keepMatplotlibPlots)
538 reply = input(
"Next iteration? [ynchpqQs] ").
strip()
542 reply = reply.split()
544 reply, args = reply[0], reply[1:]
548 if reply
in (
"",
"c",
"h",
"n",
"p",
"q",
"Q",
"s",
"y"):
552 print(
"c[ontinue without prompting] h[elp] n[o] p[db] q[uit displaying] "
553 "s[ave fileName] y[es]")
563 fileName = args.pop(0)
565 print(
"Please provide a filename")
568 print(
"Saving to %s" % fileName)
569 utils.saveSpatialCellSet(psfCellSet, fileName=fileName)
573 print(
"Unrecognised response: %s" % reply, file=sys.stderr)
579 psf, eigenValues, nEigenComponents, fitChi2 = \
580 self.
_fitPsf_fitPsf(exposure, psfCellSet, actualKernelSize, nEigenComponents)
585 if display
and reply !=
"n":
586 disp = afwDisplay.Display(frame=0)
588 utils.showPsfSpatialCells(exposure, psfCellSet, self.config.nStarPerCell, showChi2=
True,
589 symb=
"o", ctype=afwDisplay.YELLOW, ctypeBad=afwDisplay.RED,
590 size=8, display=disp)
591 if self.config.nStarPerCellSpatialFit != self.config.nStarPerCell:
592 utils.showPsfSpatialCells(exposure, psfCellSet, self.config.nStarPerCellSpatialFit,
593 symb=
"o", ctype=afwDisplay.YELLOW, ctypeBad=afwDisplay.RED,
594 size=10, display=disp)
596 disp4 = afwDisplay.Display(frame=4)
597 utils.showPsfCandidates(exposure, psfCellSet, psf=psf, display=disp4,
598 normalize=normalizeResiduals,
599 showBadCandidates=showBadCandidates)
601 if displayPsfComponents:
602 disp6 = afwDisplay.Display(frame=6)
603 utils.showPsf(psf, eigenValues, display=disp6)
606 disp7 = afwDisplay.Display(frame=7)
607 utils.showPsfMosaic(exposure, psf, display=disp7, showFwhm=
True)
608 disp7.scale(
"linear", 0, 1)
609 if displayPsfSpatialModel:
610 utils.plotPsfSpatialModel(exposure, psf, psfCellSet, showBadCandidates=
True,
611 matchKernelAmplitudes=matchKernelAmplitudes,
612 keepPlots=keepMatplotlibPlots)
624 for cell
in psfCellSet.getCellList():
625 for cand
in cell.begin(
False):
628 for cand
in cell.begin(
True):
629 src = cand.getSource()
630 if flagKey
is not None:
631 src.set(flagKey,
True)
639 if metadata
is not None:
640 metadata.set(
"spatialFitChi2", fitChi2)
641 metadata.set(
"numGoodStars", numGoodStars)
642 metadata.set(
"numAvailStars", numAvailStars)
643 metadata.set(
"avgX", avgX)
644 metadata.set(
"avgY", avgY)
648 return psf, psfCellSet
652 """Generator for Psf candidates.
654 This allows two 'for' loops to be reduced to one.
658 psfCellSet : `lsst.afw.math.SpatialCellSet`
659 SpatialCellSet of PSF candidates.
660 ignoreBad : `bool`, optional
661 Ignore candidates flagged as BAD?
665 cell : `lsst.afw.math.SpatialCell`
667 cand : `lsst.meas.algorithms.PsfCandidate`
670 for cell
in psfCellSet.getCellList():
671 for cand
in cell.begin(ignoreBad):
675 psfDeterminerRegistry.register(
"pca", PcaPsfDeterminerTask)
An ellipse core with quadrupole moments as parameters.
A collection of SpatialCells covering an entire image.
def determinePsf(self, exposure, psfCandidateList, metadata=None, flagKey=None)
def _fitPsf(self, exposure, psfCellSet, kernelSize, nEigenComponents)
Base class for PSF determiners.
Reports attempts to exceed implementation-defined length limits for some classes.
daf::base::PropertyList * list
Statistics makeStatistics(lsst::afw::image::Image< Pixel > const &img, lsst::afw::image::Mask< image::MaskPixel > const &msk, int const flags, StatisticsControl const &sctrl=StatisticsControl())
Handle a watered-down front-end to the constructor (no variance)
def numCandidatesToReject(numBadCandidates, numIter, totalIter)
def candidatesIter(psfCellSet, ignoreBad=True)
std::pair< std::shared_ptr< afw::math::LinearCombinationKernel >, std::vector< double > > createKernelFromPsfCandidates(afw::math::SpatialCellSet const &psfCells, geom::Extent2I const &dims, geom::Point2I const &xy0, int const nEigenComponents, int const spatialOrder, int const ksize, int const nStarPerCell=-1, bool const constantWeight=true, int const border=3)
Return a Kernel pointer and a list of eigenvalues resulting from analysing the provided SpatialCellSe...
int countPsfCandidates(afw::math::SpatialCellSet const &psfCells, int const nStarPerCell=-1)
Count the number of candidates in use.
std::pair< std::vector< double >, afw::math::KernelList > fitKernelParamsToImage(afw::math::LinearCombinationKernel const &kernel, Image const &image, geom::Point2D const &pos)
Fit a LinearCombinationKernel to an Image, allowing the coefficients of the components to vary.
std::pair< bool, double > fitSpatialKernelFromPsfCandidates(afw::math::Kernel *kernel, afw::math::SpatialCellSet const &psfCells, int const nStarPerCell=-1, double const tolerance=1e-5, double const lambda=0.0)
Fit spatial kernel using full-nonlinear optimization to estimate candidate amplitudes.