32 from .
import utils
as diUtils
33 from .
import diffimLib
36 """!Configuration for detecting sources on images for building a PSF-matching kernel
38 Configuration for turning detected lsst.afw.detection.FootPrints into an acceptable
39 (unmasked, high signal-to-noise, not too large or not too small) list of
40 lsst.ip.diffim.KernelSources that are used to build the Psf-matching kernel"""
42 detThreshold = pexConfig.Field(
44 doc =
"Value of footprint detection threshold",
46 check =
lambda x : x >= 3.0
48 detThresholdType = pexConfig.ChoiceField(
50 doc =
"Type of detection threshold",
51 default =
"pixel_stdev",
53 "value" :
"Use counts as the detection threshold type",
54 "stdev" :
"Use standard deviation of image plane",
55 "variance" :
"Use variance of image plane",
56 "pixel_stdev" :
"Use stdev derived from variance plane"
59 detOnTemplate = pexConfig.Field(
61 doc =
"""If true run detection on the template (image to convolve);
62 if false run detection on the science image""",
65 badMaskPlanes = pexConfig.ListField(
67 doc =
"""Mask planes that lead to an invalid detection.
68 Options: NO_DATA EDGE SAT BAD CR INTRP""",
69 default = (
"NO_DATA",
"EDGE",
"SAT")
71 fpNpixMin = pexConfig.Field(
73 doc =
"Minimum number of pixels in an acceptable Footprint",
75 check =
lambda x : x >= 5
77 fpNpixMax = pexConfig.Field(
79 doc =
"""Maximum number of pixels in an acceptable Footprint;
80 too big and the subsequent convolutions become unwieldy""",
82 check =
lambda x : x <= 500
84 fpGrowKernelScaling = pexConfig.Field(
86 doc =
"""If config.scaleByFwhm, grow the footprint based on
87 the final kernelSize. Each footprint will be
88 2*fpGrowKernelScaling*kernelSize x
89 2*fpGrowKernelScaling*kernelSize. With the value
90 of 1.0, the remaining pixels in each KernelCandiate
91 after convolution by the basis functions will be
92 eqaul to the kernel size iteslf.""",
94 check =
lambda x : x >= 1.0
96 fpGrowPix = pexConfig.Field(
98 doc =
"""Growing radius (in pixels) for each raw detection
99 footprint. The smaller the faster; however the
100 kernel sum does not converge if the stamp is too
101 small; and the kernel is not constrained at all if
102 the stamp is the size of the kernel. The grown stamp
103 is 2 * fpGrowPix pixels larger in each dimension.
104 This is overridden by fpGrowKernelScaling if scaleByFwhm""",
106 check =
lambda x : x >= 10
108 scaleByFwhm = pexConfig.Field(
110 doc =
"Scale fpGrowPix by input Fwhm?",
116 """!Base configuration for Psf-matching
118 The base configuration of the Psf-matching kernel, and of the warping, detection,
119 and background modeling subTasks."""
121 warpingConfig = pexConfig.ConfigField(
"Config for warping exposures to a common alignment",
122 afwMath.warper.WarperConfig)
123 detectionConfig = pexConfig.ConfigField(
"Controlling the detection of sources for kernel building",
125 afwBackgroundConfig = pexConfig.ConfigField(
"Controlling the Afw background fitting",
128 useAfwBackground = pexConfig.Field(
130 doc =
"Use afw background subtraction instead of ip_diffim",
133 fitForBackground = pexConfig.Field(
135 doc =
"Include terms (including kernel cross terms) for background in ip_diffim",
138 kernelBasisSet = pexConfig.ChoiceField(
140 doc =
"Type of basis set for PSF matching kernel.",
141 default =
"alard-lupton",
143 "alard-lupton" :
"""Alard-Lupton sum-of-gaussians basis set,
144 * The first term has no spatial variation
145 * The kernel sum is conserved
146 * You may want to turn off 'usePcaForSpatialKernel'""",
147 "delta-function" :
"""Delta-function kernel basis set,
148 * You may enable the option useRegularization
149 * You should seriously consider usePcaForSpatialKernel, which will also
150 enable kernel sum conservation for the delta function kernels"""
153 kernelSize = pexConfig.Field(
155 doc =
"""Number of rows/columns in the convolution kernel; should be odd-valued.
156 Modified by kernelSizeFwhmScaling if scaleByFwhm = true""",
159 scaleByFwhm = pexConfig.Field(
161 doc =
"Scale kernelSize, alardGaussians by input Fwhm",
164 kernelSizeFwhmScaling = pexConfig.Field(
166 doc =
"""How much to scale the kernel size based on the largest AL Sigma""",
168 check =
lambda x : x >= 1.0
170 kernelSizeMin = pexConfig.Field(
172 doc =
"""Minimum Kernel Size""",
175 kernelSizeMax = pexConfig.Field(
177 doc =
"""Maximum Kernel Size""",
180 spatialModelType = pexConfig.ChoiceField(
182 doc =
"Type of spatial functions for kernel and background",
183 default =
"chebyshev1",
185 "chebyshev1" :
"Chebyshev polynomial of the first kind",
186 "polynomial" :
"Standard x,y polynomial",
189 spatialKernelOrder = pexConfig.Field(
191 doc =
"Spatial order of convolution kernel variation",
193 check =
lambda x : x >= 0
195 spatialBgOrder = pexConfig.Field(
197 doc =
"Spatial order of differential background variation",
199 check =
lambda x : x >= 0
201 sizeCellX = pexConfig.Field(
203 doc =
"Size (rows) in pixels of each SpatialCell for spatial modeling",
205 check =
lambda x : x >= 32
207 sizeCellY = pexConfig.Field(
209 doc =
"Size (columns) in pixels of each SpatialCell for spatial modeling",
211 check =
lambda x : x >= 32
213 nStarPerCell = pexConfig.Field(
215 doc =
"Number of KernelCandidates in each SpatialCell to use in the spatial fitting",
217 check =
lambda x : x >= 1
219 maxSpatialIterations = pexConfig.Field(
221 doc =
"Maximum number of iterations for rejecting bad KernelCandidates in spatial fitting",
223 check =
lambda x : x >= 1
and x <= 5
225 usePcaForSpatialKernel = pexConfig.Field(
227 doc =
"""Use Pca to reduce the dimensionality of the kernel basis sets.
228 This is particularly useful for delta-function kernels.
229 Functionally, after all Cells have their raw kernels determined, we run
230 a Pca on these Kernels, re-fit the Cells using the eigenKernels and then
231 fit those for spatial variation using the same technique as for Alard-Lupton kernels.
232 If this option is used, the first term will have no spatial variation and the
233 kernel sum will be conserved.""",
236 subtractMeanForPca = pexConfig.Field(
238 doc =
"Subtract off the mean feature before doing the Pca",
241 numPrincipalComponents = pexConfig.Field(
243 doc =
"""Number of principal components to use for Pca basis, including the
244 mean kernel if requested.""",
246 check =
lambda x : x >= 3
248 singleKernelClipping = pexConfig.Field(
250 doc =
"Do sigma clipping on each raw kernel candidate",
253 kernelSumClipping = pexConfig.Field(
255 doc =
"Do sigma clipping on the ensemble of kernel sums",
258 spatialKernelClipping = pexConfig.Field(
260 doc =
"Do sigma clipping after building the spatial model",
263 checkConditionNumber = pexConfig.Field(
265 doc =
"""Test for maximum condition number when inverting a kernel matrix.
266 Anything above maxConditionNumber is not used and the candidate is set as BAD.
267 Also used to truncate inverse matrix in estimateBiasedRisk. However,
268 if you are doing any deconvolution you will want to turn this off, or use
269 a large maxConditionNumber""",
272 badMaskPlanes = pexConfig.ListField(
274 doc =
"""Mask planes to ignore when calculating diffim statistics
275 Options: NO_DATA EDGE SAT BAD CR INTRP""",
276 default = (
"NO_DATA",
"EDGE",
"SAT")
278 candidateResidualMeanMax = pexConfig.Field(
280 doc =
"""Rejects KernelCandidates yielding bad difference image quality.
281 Used by BuildSingleKernelVisitor, AssessSpatialKernelVisitor.
282 Represents average over pixels of (image/sqrt(variance)).""",
284 check =
lambda x : x >= 0.0
286 candidateResidualStdMax = pexConfig.Field(
288 doc =
"""Rejects KernelCandidates yielding bad difference image quality.
289 Used by BuildSingleKernelVisitor, AssessSpatialKernelVisitor.
290 Represents stddev over pixels of (image/sqrt(variance)).""",
292 check =
lambda x : x >= 0.0
294 useCoreStats = pexConfig.Field(
296 doc =
"""Use the core of the footprint for the quality statistics, instead of the entire footprint.
297 WARNING: if there is deconvolution we probably will need to turn this off""",
300 candidateCoreRadius = pexConfig.Field(
302 doc =
"""Radius for calculation of stats in 'core' of KernelCandidate diffim.
303 Total number of pixels used will be (2*radius)**2.
304 This is used both for 'core' diffim quality as well as ranking of
305 KernelCandidates by their total flux in this core""",
307 check =
lambda x : x >= 1
309 maxKsumSigma = pexConfig.Field(
311 doc =
"""Maximum allowed sigma for outliers from kernel sum distribution.
312 Used to reject variable objects from the kernel model""",
314 check =
lambda x : x >= 0.0
316 maxConditionNumber = pexConfig.Field(
318 doc =
"Maximum condition number for a well conditioned matrix",
320 check =
lambda x : x >= 0.0
322 conditionNumberType = pexConfig.ChoiceField(
324 doc =
"Use singular values (SVD) or eigen values (EIGENVALUE) to determine condition number",
325 default =
"EIGENVALUE",
327 "SVD" :
"Use singular values",
328 "EIGENVALUE" :
"Use eigen values (faster)",
331 maxSpatialConditionNumber = pexConfig.Field(
333 doc =
"Maximum condition number for a well conditioned spatial matrix",
335 check =
lambda x : x >= 0.0
337 iterateSingleKernel = pexConfig.Field(
339 doc =
"""Remake KernelCandidate using better variance estimate after first pass?
340 Primarily useful when convolving a single-depth image, otherwise not necessary.""",
343 constantVarianceWeighting = pexConfig.Field(
345 doc =
"""Use constant variance weighting in single kernel fitting?
346 In some cases this is better for bright star residuals.""",
349 calculateKernelUncertainty = pexConfig.Field(
351 doc =
"""Calculate kernel and background uncertainties for each kernel candidate?
352 This comes from the inverse of the covariance matrix.
353 Warning: regularization can cause problems for this step.""",
356 useBicForKernelBasis = pexConfig.Field(
358 doc =
"""Use Bayesian Information Criterion to select the number of bases going into the kernel""",
364 """!The parameters specific to the "Alard-Lupton" (sum-of-Gaussian) Psf-matching basis"""
367 PsfMatchConfig.setDefaults(self)
371 alardNGauss = pexConfig.Field(
373 doc =
"Number of Gaussians in alard-lupton basis",
375 check =
lambda x : x >= 1
377 alardDegGauss = pexConfig.ListField(
379 doc =
"Polynomial order of spatial modification of Gaussians. Must in number equal alardNGauss",
382 alardSigGauss = pexConfig.ListField(
384 doc =
"""Sigma in pixels of Gaussians (FWHM = 2.35 sigma). Must in number equal alardNGauss""",
385 default = (0.7, 1.5, 3.0),
387 alardGaussBeta = pexConfig.Field(
389 doc =
"""Default scale factor between Gaussian sigmas """,
391 check =
lambda x: x >= 0.0,
393 alardMinSig = pexConfig.Field(
395 doc =
"""Minimum Sigma (pixels) for Gaussians""",
397 check =
lambda x : x >= 0.25
399 alardDegGaussDeconv = pexConfig.Field(
401 doc =
"""Degree of spatial modification of ALL gaussians in AL basis during deconvolution""",
403 check =
lambda x : x >= 1
405 alardMinSigDeconv = pexConfig.Field(
407 doc =
"""Minimum Sigma (pixels) for Gaussians during deconvolution;
408 make smaller than alardMinSig as this is only indirectly used""",
410 check =
lambda x : x >= 0.25
412 alardNGaussDeconv = pexConfig.Field(
414 doc =
"Number of Gaussians in AL basis during deconvolution",
416 check =
lambda x : x >= 1
422 """!The parameters specific to the delta-function (one basis per-pixel) Psf-matching basis"""
425 PsfMatchConfig.setDefaults(self)
432 useRegularization = pexConfig.Field(
434 doc =
"Use regularization to smooth the delta function kernels",
437 regularizationType = pexConfig.ChoiceField(
439 doc =
"Type of regularization.",
440 default =
"centralDifference",
442 "centralDifference":
"Penalize second derivative using 2-D stencil of central finite difference",
443 "forwardDifference":
"Penalize first, second, third derivatives using forward finite differeces"
446 centralRegularizationStencil = pexConfig.ChoiceField(
448 doc =
"Type of stencil to approximate central derivative (for centralDifference only)",
451 5 :
"5-point stencil including only adjacent-in-x,y elements",
452 9 :
"9-point stencil including diagonal elements"
455 forwardRegularizationOrders = pexConfig.ListField(
457 doc =
"Array showing which order derivatives to penalize (for forwardDifference only)",
459 itemCheck =
lambda x: (x > 0)
and (x < 4)
461 regularizationBorderPenalty = pexConfig.Field(
463 doc =
"Value of the penalty for kernel border pixels",
465 check =
lambda x : x >= 0.0
467 lambdaType = pexConfig.ChoiceField(
469 doc =
"How to choose the value of the regularization strength",
470 default =
"absolute",
472 "absolute" :
"Use lambdaValue as the value of regularization strength",
473 "relative" :
"Use lambdaValue as fraction of the default regularization strength (N.R. 18.5.8)",
474 "minimizeBiasedRisk" :
"Minimize biased risk estimate",
475 "minimizeUnbiasedRisk" :
"Minimize unbiased risk estimate",
478 lambdaValue = pexConfig.Field(
480 doc =
"Value used for absolute determinations of regularization strength",
483 lambdaScaling = pexConfig.Field(
485 doc =
"Fraction of the default lambda strength (N.R. 18.5.8) to use. 1e-4 or 1e-5",
488 lambdaStepType = pexConfig.ChoiceField(
490 doc =
"""If a scan through lambda is needed (minimizeBiasedRisk, minimizeUnbiasedRisk),
491 use log or linear steps""",
494 "log" :
"Step in log intervals; e.g. lambdaMin, lambdaMax, lambdaStep = -1.0, 2.0, 0.1",
495 "linear" :
"Step in linear intervals; e.g. lambdaMin, lambdaMax, lambdaStep = 0.1, 100, 0.1",
498 lambdaMin = pexConfig.Field(
500 doc =
"""If scan through lambda needed (minimizeBiasedRisk, minimizeUnbiasedRisk),
501 start at this value. If lambdaStepType = log:linear, suggest -1:0.1""",
504 lambdaMax = pexConfig.Field(
506 doc =
"""If scan through lambda needed (minimizeBiasedRisk, minimizeUnbiasedRisk),
507 stop at this value. If lambdaStepType = log:linear, suggest 2:100""",
510 lambdaStep = pexConfig.Field(
512 doc =
"""If scan through lambda needed (minimizeBiasedRisk, minimizeUnbiasedRisk),
513 step in these increments. If lambdaStepType = log:linear, suggest 0.1:0.1""",
528 \anchor PsfMatchTask_
530 \brief Base class for Psf Matching; should not be called directly
532 \section ip_diffim_psfmatch_Contents Contents
534 - \ref ip_diffim_psfmatch_Purpose
535 - \ref ip_diffim_psfmatch_Initialize
536 - \ref ip_diffim_psfmatch_IO
537 - \ref ip_diffim_psfmatch_Config
538 - \ref ip_diffim_psfmatch_Metadata
539 - \ref ip_diffim_psfmatch_Debug
540 - \ref ip_diffim_psfmatch_Example
542 #-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
544 \section ip_diffim_psfmatch_Purpose Description
546 PsfMatchTask is a base class that implements the core functionality for matching the
547 Psfs of two images using a spatially varying Psf-matching lsst.afw.math.LinearCombinationKernel.
548 The Task requires the user to provide an instance of an lsst.afw.math.SpatialCellSet,
549 filled with lsst.ip.diffim.KernelCandidate instances, and an lsst.afw.math.KernelList
550 of basis shapes that will be used for the decomposition. If requested, the Task
551 also performs background matching and returns the differential background model as an
552 lsst.afw.math.Kernel.SpatialFunction.
554 #-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
556 \section ip_diffim_psfmatch_Initialize Task initialization
558 \copydoc \_\_init\_\_
560 #-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
562 \section ip_diffim_psfmatch_IO Invoking the Task
564 As a base class, this Task is not directly invoked. However, run() methods that are
565 implemented on derived classes will make use of the core _solve() functionality,
566 which defines a sequence of lsst.afw.math.CandidateVisitor classes that iterate
567 through the KernelCandidates, first building up a per-candidate solution and then
568 building up a spatial model from the ensemble of candidates. Sigma clipping is
569 performed using the mean and standard deviation of all kernel sums (to reject
570 variable objects), on the per-candidate substamp diffim residuals
571 (to indicate a bad choice of kernel basis shapes for that particular object),
572 and on the substamp diffim residuals using the spatial kernel fit (to indicate a bad
573 choice of spatial kernel order, or poor constraints on the spatial model). The
574 _diagnostic() method logs information on the quality of the spatial fit, and also
575 modifies the Task metadata.
577 #-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
579 \section ip_diffim_psfmatch_Config Configuration parameters
581 See \ref PsfMatchConfig, \ref PsfMatchConfigAL, \ref PsfMatchConfigDF, and \ref DetectionConfig.
583 #-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
585 \section ip_diffim_psfmatch_Metadata Quantities set in Metadata
589 <DT> spatialConditionNum <DD> Condition number of the spatial kernel fit;
590 via \link lsst.ip.diffim.PsfMatchTask._diagnostic PsfMatchTask._diagnostic \endlink </DD> </DT>
591 <DT> spatialKernelSum <DD> Kernel sum (10^{-0.4 * Δ zeropoint}) of the spatial Psf-matching kernel;
592 via \link lsst.ip.diffim.PsfMatchTask._diagnostic PsfMatchTask._diagnostic \endlink </DD> </DT>
594 <DT> ALBasisNGauss <DD> If using sum-of-Gaussian basis, the number of gaussians used;
595 via \link lsst.ip.diffim.makeKernelBasisList.generateAlardLuptonBasisList
596 generateAlardLuptonBasisList\endlink </DD> </DT>
597 <DT> ALBasisDegGauss <DD> If using sum-of-Gaussian basis, the degree of spatial variation of the Gaussians;
598 via \link lsst.ip.diffim.makeKernelBasisList.generateAlardLuptonBasisList
599 generateAlardLuptonBasisList\endlink </DD> </DT>
600 <DT> ALBasisSigGauss <DD> If using sum-of-Gaussian basis, the widths (sigma) of the Gaussians;
601 via \link lsst.ip.diffim.makeKernelBasisList.generateAlardLuptonBasisList
602 generateAlardLuptonBasisList\endlink </DD> </DT>
603 <DT> ALKernelSize <DD> If using sum-of-Gaussian basis, the kernel size;
604 via \link lsst.ip.diffim.makeKernelBasisList.generateAlardLuptonBasisList
605 generateAlardLuptonBasisList\endlink </DD> </DT>
607 <DT> NFalsePositivesTotal <DD> Total number of diaSources;
608 via \link lsst.ip.diffim.KernelCandidateQa.aggregate KernelCandidateQa.aggregate\endlink </DD> </DT>
609 <DT> NFalsePositivesRefAssociated <DD> Number of diaSources that associate with the reference catalog;
610 via \link lsst.ip.diffim.KernelCandidateQa.aggregate KernelCandidateQa.aggregate\endlink </DD> </DT>
611 <DT> NFalsePositivesRefAssociated <DD> Number of diaSources that associate with the source catalog;
612 via \link lsst.ip.diffim.KernelCandidateQa.aggregate KernelCandidateQa.aggregate\endlink </DD> </DT>
613 <DT> NFalsePositivesUnassociated <DD> Number of diaSources that are orphans;
614 via \link lsst.ip.diffim.KernelCandidateQa.aggregate KernelCandidateQa.aggregate\endlink </DD> </DT>
615 <DT> metric_MEAN <DD> Mean value of substamp diffim quality metrics across all KernelCandidates,
616 for both the per-candidate (LOCAL) and SPATIAL residuals;
617 via \link lsst.ip.diffim.KernelCandidateQa.aggregate KernelCandidateQa.aggregate\endlink </DD> </DT>
618 <DT> metric_MEDIAN <DD> Median value of substamp diffim quality metrics across all KernelCandidates,
619 for both the per-candidate (LOCAL) and SPATIAL residuals;
620 via \link lsst.ip.diffim.KernelCandidateQa.aggregate KernelCandidateQa.aggregate\endlink </DD> </DT>
621 <DT> metric_STDEV <DD> Standard deviation of substamp diffim quality metrics across all KernelCandidates,
622 for both the per-candidate (LOCAL) and SPATIAL residuals;
623 via \link lsst.ip.diffim.KernelCandidateQa.aggregate KernelCandidateQa.aggregate\endlink </DD> </DT>
627 #-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
629 \section ip_diffim_psfmatch_Debug Debug variables
632 The \link lsst.pipe.base.cmdLineTask.CmdLineTask command line task\endlink interface supports a
633 flag \c -d/--debug to import \b debug.py from your \c PYTHONPATH. The relevant contents of debug.py
634 for this Task include:
640 di = lsstDebug.getInfo(name)
641 if name == "lsst.ip.diffim.psfMatch":
642 di.display = True # enable debug output
643 di.maskTransparency = 80 # ds9 mask transparency
644 di.displayCandidates = True # show all the candidates and residuals
645 di.displayKernelBasis = False # show kernel basis functions
646 di.displayKernelMosaic = True # show kernel realized across the image
647 di.plotKernelSpatialModel = False # show coefficients of spatial model
648 di.showBadCandidates = True # show the bad candidates (red) along with good (green)
650 lsstDebug.Info = DebugInfo
655 Note that if you want addional logging info, you may add to your scripts:
657 import lsst.pex.logging as pexLog
658 pexLog.Trace_setVerbosity('lsst.ip.diffim', 5)
661 #-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
663 \section ip_diffim_psfmatch_Example Example code
665 As a base class, there is no example code for PsfMatchTask.
666 However, see \link lsst.ip.diffim.imagePsfMatch.ImagePsfMatchTask ImagePsfMatchTask\endlink,
667 \link lsst.ip.diffim.snapPsfMatch.SnapPsfMatchTask SnapPsfMatchTask\endlink, and
668 \link lsst.ip.diffim.modelPsfMatch.ModelPsfMatchTask ModelPsfMatchTask\endlink.
670 #-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
673 ConfigClass = PsfMatchConfig
674 _DefaultName =
"psfMatch"
677 """!Create the psf-matching Task
679 \param *args arguments to be passed to lsst.pipe.base.task.Task.__init__
680 \param **kwargs keyword arguments to be passed to lsst.pipe.base.task.Task.__init__
682 The initialization sets the Psf-matching kernel configuration using the value of
683 self.config.kernel.active. If the kernel is requested with regularization to moderate
684 the bias/variance tradeoff, currently only used when a delta function kernel basis
685 is provided, it creates a regularization matrix stored as member variable
688 pipeBase.Task.__init__(self, *args, **kwargs)
692 if 'useRegularization' in self.kConfig.keys():
698 self.
hMat = diffimLib.makeRegularizationMatrix(pexConfig.makePolicy(self.
kConfig))
700 def _diagnostic(self, kernelCellSet, spatialSolution, spatialKernel, spatialBg):
701 """!Provide logging diagnostics on quality of spatial kernel fit
703 @param kernelCellSet: Cellset that contains the KernelCandidates used in the fitting
704 @param spatialSolution: KernelSolution of best-fit
705 @param spatialKernel: Best-fit spatial Kernel model
706 @param spatialBg: Best-fit spatial background model
710 kImage = afwImage.ImageD(spatialKernel.getDimensions())
711 kSum = spatialKernel.computeImage(kImage,
False)
712 self.log.info(
"Final spatial kernel sum %.3f" % (kSum))
715 conditionNum = spatialSolution.getConditionNumber(
716 getattr(diffimLib.KernelSolution, self.kConfig.conditionNumberType))
717 self.log.info(
"Spatial model condition number %.3e" % (conditionNum))
719 if conditionNum < 0.0:
720 self.log.warn(
"Condition number is negative (%.3e)" % (conditionNum))
721 if conditionNum > self.kConfig.maxSpatialConditionNumber:
722 self.log.warn(
"Spatial solution exceeds max condition number (%.3e > %.3e)" % (
723 conditionNum, self.kConfig.maxSpatialConditionNumber))
725 self.metadata.set(
"spatialConditionNum", conditionNum)
726 self.metadata.set(
"spatialKernelSum", kSum)
729 nBasisKernels = spatialKernel.getNBasisKernels()
730 nKernelTerms = spatialKernel.getNSpatialParameters()
731 if nKernelTerms == 0:
735 nBgTerms = spatialBg.getNParameters()
737 if spatialBg.getParameters()[0] == 0.0:
743 for cell
in kernelCellSet.getCellList():
744 for cand
in cell.begin(
False):
745 cand = diffimLib.cast_KernelCandidateF(cand)
747 if cand.getStatus() == afwMath.SpatialCellCandidate.GOOD:
749 if cand.getStatus() == afwMath.SpatialCellCandidate.BAD:
752 self.log.info(
"Doing stats of kernel candidates used in the spatial fit.")
756 self.log.warn(
"Many more candidates rejected than accepted; %d total, %d rejected, %d used" % (
759 self.log.info(
"%d candidates total, %d rejected, %d used" % (nTot, nBad, nGood))
762 if nGood < nKernelTerms:
763 self.log.warn(
"Spatial kernel model underconstrained; %d candidates, %d terms, %d bases" % (
764 nGood, nKernelTerms, nBasisKernels))
765 self.log.warn(
"Consider lowering the spatial order")
766 elif nGood <= 2*nKernelTerms:
767 self.log.warn(
"Spatial kernel model poorly constrained; %d candidates, %d terms, %d bases" % (
768 nGood, nKernelTerms, nBasisKernels))
769 self.log.warn(
"Consider lowering the spatial order")
771 self.log.info(
"Spatial kernel model well constrained; %d candidates, %d terms, %d bases" % (
772 nGood, nKernelTerms, nBasisKernels))
775 self.log.warn(
"Spatial background model underconstrained; %d candidates, %d terms" % (
777 self.log.warn(
"Consider lowering the spatial order")
778 elif nGood <= 2*nBgTerms:
779 self.log.warn(
"Spatial background model poorly constrained; %d candidates, %d terms" % (
781 self.log.warn(
"Consider lowering the spatial order")
783 self.log.info(
"Spatial background model appears well constrained; %d candidates, %d terms" % (
787 """!Provide visualization of the inputs and ouputs to the Psf-matching code
789 @param kernelCellSet: the SpatialCellSet used in determining the matching kernel and background
790 @param spatialKernel: spatially varying Psf-matching kernel
791 @param spatialBackground: spatially varying background-matching function
797 displayKernelMosaic =
lsstDebug.Info(__name__).displayKernelMosaic
798 plotKernelSpatialModel =
lsstDebug.Info(__name__).plotKernelSpatialModel
801 if not maskTransparency:
803 ds9.setMaskTransparency(maskTransparency)
805 if displayCandidates:
806 diUtils.showKernelCandidates(kernelCellSet, kernel=spatialKernel, background=spatialBackground,
807 frame=lsstDebug.frame,
808 showBadCandidates=showBadCandidates)
810 diUtils.showKernelCandidates(kernelCellSet, kernel=spatialKernel, background=spatialBackground,
811 frame=lsstDebug.frame,
812 showBadCandidates=showBadCandidates,
815 diUtils.showKernelCandidates(kernelCellSet, kernel=spatialKernel, background=spatialBackground,
816 frame=lsstDebug.frame,
817 showBadCandidates=showBadCandidates,
821 if displayKernelBasis:
822 diUtils.showKernelBasis(spatialKernel, frame=lsstDebug.frame)
825 if displayKernelMosaic:
826 diUtils.showKernelMosaic(kernelCellSet.getBBox(), spatialKernel, frame=lsstDebug.frame)
829 if plotKernelSpatialModel:
830 diUtils.plotKernelSpatialModel(spatialKernel, kernelCellSet, showBadCandidates=showBadCandidates)
834 """!Create Principal Component basis
836 If a principal component analysis is requested, typically when using a delta function basis,
837 perform the PCA here and return a new basis list containing the new principal components.
839 @param kernelCellSet: a SpatialCellSet containing KernelCandidates, from which components are derived
840 @param nStarPerCell: the number of stars per cell to visit when doing the PCA
841 @param policy: input policy controlling the single kernel visitor
844 - nRejectedPca: number of KernelCandidates rejected during PCA loop
845 - spatialBasisList: basis list containing the principal shapes as Kernels
848 nComponents = self.kConfig.numPrincipalComponents
849 imagePca = diffimLib.KernelPcaD()
850 importStarVisitor = diffimLib.KernelPcaVisitorF(imagePca)
851 kernelCellSet.visitCandidates(importStarVisitor, nStarPerCell)
852 if self.kConfig.subtractMeanForPca:
853 importStarVisitor.subtractMean()
856 eigenValues = imagePca.getEigenValues()
857 pcaBasisList = importStarVisitor.getEigenKernels()
859 eSum = num.sum(eigenValues)
861 raise RuntimeError(
"Eigenvalues sum to zero")
862 for j
in range(len(eigenValues)):
864 "Eigenvalue %d : %f (%f)" % (j, eigenValues[j], eigenValues[j]/eSum))
866 nToUse = min(nComponents, len(eigenValues))
868 for j
in range(nToUse):
870 kimage = afwImage.ImageD(pcaBasisList[j].getDimensions())
871 pcaBasisList[j].computeImage(kimage,
False)
872 if not (
True in num.isnan(kimage.getArray())):
873 trimBasisList.push_back(pcaBasisList[j])
876 spatialBasisList = diffimLib.renormalizeKernelList(trimBasisList)
879 singlekvPca = diffimLib.BuildSingleKernelVisitorF(spatialBasisList, policy)
880 singlekvPca.setSkipBuilt(
False)
881 kernelCellSet.visitCandidates(singlekvPca, nStarPerCell)
882 singlekvPca.setSkipBuilt(
True)
883 nRejectedPca = singlekvPca.getNRejected()
885 return nRejectedPca, spatialBasisList
889 """!Fill a SpatialCellSet with KernelCandidates for the Psf-matching process;
890 override in derived classes"""
894 def _solve(self, kernelCellSet, basisList, returnOnExcept=False):
895 """!Solve for the PSF matching kernel
897 @param kernelCellSet: a SpatialCellSet to use in determining the matching kernel
898 (typically as provided by _buildCellSet)
899 @param basisList: list of Kernels to be used in the decomposition of the spatially varying kernel
900 (typically as provided by makeKernelBasisList)
901 @param returnOnExcept: if True then return (None, None) if an error occurs, else raise the exception
904 - psfMatchingKernel: PSF matching kernel
905 - backgroundModel: differential background model
907 Raise Exception if unable to determine PSF matching kernel and returnOnExcept False
913 maxSpatialIterations = self.kConfig.maxSpatialIterations
914 nStarPerCell = self.kConfig.nStarPerCell
915 usePcaForSpatialKernel = self.kConfig.usePcaForSpatialKernel
918 policy = pexConfig.makePolicy(self.
kConfig)
920 singlekv = diffimLib.BuildSingleKernelVisitorF(basisList, policy, self.
hMat)
922 singlekv = diffimLib.BuildSingleKernelVisitorF(basisList, policy)
925 ksv = diffimLib.KernelSumVisitorF(policy)
932 while (thisIteration < maxSpatialIterations):
936 while (nRejectedSkf != 0):
938 "Building single kernels...")
939 kernelCellSet.visitCandidates(singlekv, nStarPerCell)
940 nRejectedSkf = singlekv.getNRejected()
942 "Iteration %d, rejected %d candidates due to initial kernel fit" %
943 (thisIteration, nRejectedSkf))
947 ksv.setMode(diffimLib.KernelSumVisitorF.AGGREGATE)
948 kernelCellSet.visitCandidates(ksv, nStarPerCell)
949 ksv.processKsumDistribution()
950 ksv.setMode(diffimLib.KernelSumVisitorF.REJECT)
951 kernelCellSet.visitCandidates(ksv, nStarPerCell)
953 nRejectedKsum = ksv.getNRejected()
955 "Iteration %d, rejected %d candidates due to kernel sum" %
956 (thisIteration, nRejectedKsum))
960 if nRejectedKsum > 0:
969 if (usePcaForSpatialKernel):
970 pexLog.Trace(self.log.getName()+
"._solve", 1,
"Building Pca basis")
972 nRejectedPca, spatialBasisList = self.
_createPcaBasis(kernelCellSet, nStarPerCell, policy)
974 "Iteration %d, rejected %d candidates due to Pca kernel fit" % (
975 thisIteration, nRejectedPca))
986 if (nRejectedPca > 0):
990 spatialBasisList = basisList
993 regionBBox = kernelCellSet.getBBox()
994 spatialkv = diffimLib.BuildSpatialKernelVisitorF(spatialBasisList, regionBBox, policy)
995 kernelCellSet.visitCandidates(spatialkv, nStarPerCell)
996 spatialkv.solveLinearEquation()
998 "Spatial kernel built with %d candidates" % (spatialkv.getNCandidates()))
999 spatialKernel, spatialBackground = spatialkv.getSolutionPair()
1002 assesskv = diffimLib.AssessSpatialKernelVisitorF(spatialKernel, spatialBackground, policy)
1003 kernelCellSet.visitCandidates(assesskv, nStarPerCell)
1004 nRejectedSpatial = assesskv.getNRejected()
1005 nGoodSpatial = assesskv.getNGood()
1007 "Iteration %d, rejected %d candidates due to spatial kernel fit" % (
1008 thisIteration, nRejectedSpatial))
1010 "%d candidates used in fit" % (nGoodSpatial))
1013 if nGoodSpatial == 0
and nRejectedSpatial == 0:
1014 raise RuntimeError(
"No kernel candidates for spatial fit")
1016 if nRejectedSpatial == 0:
1024 if (nRejectedSpatial > 0)
and (thisIteration == maxSpatialIterations):
1025 pexLog.Trace(self.log.getName()+
"._solve", 2,
"Final spatial fit")
1026 if (usePcaForSpatialKernel):
1027 nRejectedPca, spatialBasisList = self.
_createPcaBasis(kernelCellSet, nStarPerCell, policy)
1028 regionBBox = kernelCellSet.getBBox()
1029 spatialkv = diffimLib.BuildSpatialKernelVisitorF(spatialBasisList, regionBBox, policy)
1030 kernelCellSet.visitCandidates(spatialkv, nStarPerCell)
1031 spatialkv.solveLinearEquation()
1033 "Spatial kernel built with %d candidates" % (spatialkv.getNCandidates()))
1034 spatialKernel, spatialBackground = spatialkv.getSolutionPair()
1036 spatialSolution = spatialkv.getKernelSolution()
1039 pexLog.Trace(self.log.getName()+
"._solve", 1,
"ERROR: Unable to calculate psf matching kernel")
1040 pexLog.Trace(self.log.getName()+
"._solve", 2, e.args[0].what())
1042 except Exception, e:
1043 pexLog.Trace(self.log.getName()+
"._solve", 1,
"ERROR: Unable to calculate psf matching kernel")
1044 pexLog.Trace(self.log.getName()+
"._solve", 2, e.args[0])
1049 "Total time to compute the spatial kernel : %.2f s" % (t1 - t0))
1052 self.
_displayDebug(kernelCellSet, spatialKernel, spatialBackground)
1054 self.
_diagnostic(kernelCellSet, spatialSolution, spatialKernel, spatialBackground)
1056 return spatialSolution, spatialKernel, spatialBackground
1058 PsfMatch=PsfMatchTask
Configuration for detecting sources on images for building a PSF-matching kernel. ...
The parameters specific to the "Alard-Lupton" (sum-of-Gaussian) Psf-matching basis.
Base configuration for Psf-matching.
The parameters specific to the delta-function (one basis per-pixel) Psf-matching basis.
Base class for Psf Matching; should not be called directly.
def _displayDebug
Provide visualization of the inputs and ouputs to the Psf-matching code.
def _diagnostic
Provide logging diagnostics on quality of spatial kernel fit.
limited backward compatibility to the DC2 run-time trace facilities
def __init__
Create the psf-matching Task.
def _buildCellSet
Fill a SpatialCellSet with KernelCandidates for the Psf-matching process; override in derived classes...
def _createPcaBasis
Create Principal Component basis.
std::vector< boost::shared_ptr< Kernel > > KernelList
tuple usePcaForSpatialKernel
def _solve
Solve for the PSF matching kernel.
tuple useBicForKernelBasis