lsst.ip.diffim  17.0.1-6-ga5f3de7+6
diffimTools.py
Go to the documentation of this file.
1 # This file is part of ip_diffim.
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 __all__ = ["backgroundSubtract", "writeKernelCellSet", "sourceToFootprintList", "NbasisEvaluator"]
23 
24 # python
25 import time
26 import os
27 from collections import Counter
28 import numpy as np
29 
30 # all the c++ level classes and routines
31 from . import diffimLib
32 
33 # all the other LSST packages
34 import lsst.afw.geom as afwGeom
35 import lsst.afw.image as afwImage
36 import lsst.afw.table as afwTable
37 import lsst.afw.detection as afwDetect
38 import lsst.afw.math.mathLib as afwMath
39 from lsst.log import Log
40 import lsst.pex.config as pexConfig
41 from .makeKernelBasisList import makeKernelBasisList
42 
43 # Helper functions for ipDiffim; mostly viewing of results and writing
44 # debugging info to disk.
45 
46 
49 
50 
51 def makeFlatNoiseImage(mi, seedStat=afwMath.MAX):
52  img = mi.getImage()
53  seed = int(10.*afwMath.makeStatistics(mi.getImage(), seedStat).getValue() + 1)
54  rdm = afwMath.Random(afwMath.Random.MT19937, seed)
55  rdmImage = img.Factory(img.getDimensions())
56  afwMath.randomGaussianImage(rdmImage, rdm)
57  return rdmImage
58 
59 
61  """Return a Poisson noise image based on im
62 
63  Parameters
64  ----------
65  im : `lsst.afw.image.Image`
66  image; the output image has the same dimensions and shape
67  and its expectation value is the value of ``im`` at each pixel
68 
69  Returns
70  -------
71  noiseIm : `lsst.afw.image.Image`
72  TODO: DM-17458
73 
74  Notes
75  -----
76  - Warning: This uses an undocumented numpy API (the documented API
77  uses a single float expectation value instead of an array).
78 
79  - Uses numpy.random; you may wish to call numpy.random.seed first.
80  """
81  import numpy.random as rand
82  imArr = im.getArray()
83  noiseIm = im.Factory(im.getBBox())
84  noiseArr = noiseIm.getArray()
85 
86  with np.errstate(invalid='ignore'):
87  intNoiseArr = rand.poisson(imArr)
88 
89  noiseArr[:, :] = intNoiseArr.astype(noiseArr.dtype)
90  return noiseIm
91 
92 
97 
98 
99 def fakeCoeffs():
100  kCoeffs = ((1.0, 0.0, 0.0),
101  (0.005, -0.000001, 0.000001),
102  (0.005, 0.000004, 0.000004),
103  (-0.001, -0.000030, 0.000030),
104  (-0.001, 0.000015, 0.000015),
105  (-0.005, -0.000050, 0.000050))
106  return kCoeffs
107 
108 
109 def makeFakeKernelSet(sizeCell=128, nCell=3,
110  deltaFunctionCounts=1.e4, tGaussianWidth=1.0,
111  addNoise=True, bgValue=100., display=False):
112  """TODO: DM-17458
113 
114  Parameters
115  ----------
116  sizeCell : `int`, optional
117  TODO: DM-17458
118  nCell : `int`, optional
119  TODO: DM-17458
120  deltaFunctionCounts : `float`, optional
121  TODO: DM-17458
122  tGaussianWidth : `float`, optional
123  TODO: DM-17458
124  addNoise : `bool`, optional
125  TODO: DM-17458
126  bgValue : `float`, optional
127  TODO: DM-17458
128  display : `bool`, optional
129  TODO: DM-17458
130 
131  Returns
132  -------
133  TODO: DM-17458
134  TODO: DM-17458
135  """
136  from . import imagePsfMatch
137  configFake = imagePsfMatch.ImagePsfMatchConfig()
138  configFake.kernel.name = "AL"
139  subconfigFake = configFake.kernel.active
140  subconfigFake.alardNGauss = 1
141  subconfigFake.alardSigGauss = [2.5, ]
142  subconfigFake.alardDegGauss = [2, ]
143  subconfigFake.sizeCellX = sizeCell
144  subconfigFake.sizeCellY = sizeCell
145  subconfigFake.spatialKernelOrder = 1
146  subconfigFake.spatialModelType = "polynomial"
147  subconfigFake.singleKernelClipping = False # variance is a hack
148  subconfigFake.spatialKernelClipping = False # variance is a hack
149  if bgValue > 0.0:
150  subconfigFake.fitForBackground = True
151 
152  policyFake = pexConfig.makePolicy(subconfigFake)
153 
154  basisList = makeKernelBasisList(subconfigFake)
155  kSize = subconfigFake.kernelSize
156 
157  # This sets the final extent of each convolved delta function
158  gaussKernelWidth = sizeCell//2
159 
160  # This sets the scale over which pixels are correlated in the
161  # spatial convolution; should be at least as big as the kernel you
162  # are trying to fit for
163  spatialKernelWidth = kSize
164 
165  # Number of bad pixels due to convolutions
166  border = (gaussKernelWidth + spatialKernelWidth)//2
167 
168  # Make a fake image with a matrix of delta functions
169  totalSize = nCell*sizeCell + 2*border
170  tim = afwImage.ImageF(afwGeom.Extent2I(totalSize, totalSize))
171  for x in range(nCell):
172  for y in range(nCell):
173  tim[x*sizeCell + sizeCell//2 + border - 1,
174  y*sizeCell + sizeCell//2 + border - 1,
175  afwImage.LOCAL] = deltaFunctionCounts
176 
177  # Turn this into stars with a narrow width; conserve counts
178  gaussFunction = afwMath.GaussianFunction2D(tGaussianWidth, tGaussianWidth)
179  gaussKernel = afwMath.AnalyticKernel(gaussKernelWidth, gaussKernelWidth, gaussFunction)
180  cim = afwImage.ImageF(tim.getDimensions())
181  afwMath.convolve(cim, tim, gaussKernel, True)
182  tim = cim
183 
184  # Trim off border pixels
185  bbox = gaussKernel.shrinkBBox(tim.getBBox(afwImage.LOCAL))
186  tim = afwImage.ImageF(tim, bbox, afwImage.LOCAL)
187 
188  # Now make a science image which is this convolved with some
189  # spatial function. Use input basis list.
190  polyFunc = afwMath.PolynomialFunction2D(1)
191  kCoeffs = fakeCoeffs()
192  nToUse = min(len(kCoeffs), len(basisList))
193 
194  # Make the full convolved science image
195  sKernel = afwMath.LinearCombinationKernel(basisList[:nToUse], polyFunc)
196  sKernel.setSpatialParameters(kCoeffs[:nToUse])
197  sim = afwImage.ImageF(tim.getDimensions())
198  afwMath.convolve(sim, tim, sKernel, True)
199 
200  # Get the good subregion
201  bbox = sKernel.shrinkBBox(sim.getBBox(afwImage.LOCAL))
202 
203  # Add background
204  sim += bgValue
205 
206  # Watch out for negative values
207  tim += 2*np.abs(np.min(tim.getArray()))
208 
209  # Add noise?
210  if addNoise:
211  sim = makePoissonNoiseImage(sim)
212  tim = makePoissonNoiseImage(tim)
213 
214  # And turn into MaskedImages
215  sim = afwImage.ImageF(sim, bbox, afwImage.LOCAL)
216  svar = afwImage.ImageF(sim, True)
217  smask = afwImage.Mask(sim.getDimensions())
218  smask.set(0x0)
219  sMi = afwImage.MaskedImageF(sim, smask, svar)
220 
221  tim = afwImage.ImageF(tim, bbox, afwImage.LOCAL)
222  tvar = afwImage.ImageF(tim, True)
223  tmask = afwImage.Mask(tim.getDimensions())
224  tmask.set(0x0)
225  tMi = afwImage.MaskedImageF(tim, tmask, tvar)
226 
227  if display:
228  import lsst.afw.display as afwDisplay
229  afwDisplay.Display(frame=1).mtv(tMi)
230  afwDisplay.Display(frame=2).mtv(sMi)
231 
232  # Finally, make a kernelSet from these 2 images
233  kernelCellSet = afwMath.SpatialCellSet(afwGeom.Box2I(afwGeom.Point2I(0, 0),
234  afwGeom.Extent2I(sizeCell*nCell,
235  sizeCell*nCell)),
236  sizeCell,
237  sizeCell)
238  stampHalfWidth = 2*kSize
239  for x in range(nCell):
240  for y in range(nCell):
241  xCoord = x*sizeCell + sizeCell//2
242  yCoord = y*sizeCell + sizeCell//2
243  p0 = afwGeom.Point2I(xCoord - stampHalfWidth,
244  yCoord - stampHalfWidth)
245  p1 = afwGeom.Point2I(xCoord + stampHalfWidth,
246  yCoord + stampHalfWidth)
247  bbox = afwGeom.Box2I(p0, p1)
248  tsi = afwImage.MaskedImageF(tMi, bbox, origin=afwImage.LOCAL)
249  ssi = afwImage.MaskedImageF(sMi, bbox, origin=afwImage.LOCAL)
250 
251  kc = diffimLib.makeKernelCandidate(xCoord, yCoord, tsi, ssi, policyFake)
252  kernelCellSet.insertCandidate(kc)
253 
254  tMi.setXY0(0, 0)
255  sMi.setXY0(0, 0)
256  return tMi, sMi, sKernel, kernelCellSet, configFake
257 
258 
259 
262 
263 def backgroundSubtract(config, maskedImages):
264  """Subtract the background from masked images.
265 
266  Parameters
267  ----------
268  config : TODO: DM-17458
269  TODO: DM-17458
270  maskedImages : `list` of `lsst.afw.image.MaskedImage`
271  TODO: DM-17458
272 
273  Returns
274  -------
275  TODO: DM-17458
276  TODO: DM-17458
277  """
278  backgrounds = []
279  t0 = time.time()
280  algorithm = config.algorithm
281  binsize = config.binSize
282  undersample = config.undersampleStyle
283  bctrl = afwMath.BackgroundControl(algorithm)
284  bctrl.setUndersampleStyle(undersample)
285  for maskedImage in maskedImages:
286  bctrl.setNxSample(maskedImage.getWidth()//binsize + 1)
287  bctrl.setNySample(maskedImage.getHeight()//binsize + 1)
288  image = maskedImage.getImage()
289  backobj = afwMath.makeBackground(image, bctrl)
290 
291  image -= backobj.getImageF()
292  backgrounds.append(backobj.getImageF())
293  del backobj
294 
295  t1 = time.time()
296  logger = Log.getLogger("ip.diffim.backgroundSubtract")
297  logger.debug("Total time for background subtraction : %.2f s", (t1 - t0))
298  return backgrounds
299 
300 
303 
304 
305 def writeKernelCellSet(kernelCellSet, psfMatchingKernel, backgroundModel, outdir):
306  """TODO: DM-17458
307 
308  Parameters
309  ----------
310  kernelCellSet : TODO: DM-17458
311  TODO: DM-17458
312  psfMatchingKernel : TODO: DM-17458
313  TODO: DM-17458
314  backgroundModel : TODO: DM-17458
315  TODO: DM-17458
316  outdir : TODO: DM-17458
317  TODO: DM-17458
318  """
319  if not os.path.isdir(outdir):
320  os.makedirs(outdir)
321 
322  for cell in kernelCellSet.getCellList():
323  for cand in cell.begin(False): # False = include bad candidates
324  if cand.getStatus() == afwMath.SpatialCellCandidate.GOOD:
325  xCand = int(cand.getXCenter())
326  yCand = int(cand.getYCenter())
327  idCand = cand.getId()
328  diffIm = cand.getDifferenceImage(diffimLib.KernelCandidateF.ORIG)
329  kernel = cand.getKernelImage(diffimLib.KernelCandidateF.ORIG)
330  diffIm.writeFits(os.path.join(outdir, 'diffim_c%d_x%d_y%d.fits' % (idCand, xCand, yCand)))
331  kernel.writeFits(os.path.join(outdir, 'kernel_c%d_x%d_y%d.fits' % (idCand, xCand, yCand)))
332 
333  # Diffim from spatial model
334  ski = afwImage.ImageD(kernel.getDimensions())
335  psfMatchingKernel.computeImage(ski, False, xCand, yCand)
336  sk = afwMath.FixedKernel(ski)
337  sbg = backgroundModel(xCand, yCand)
338  sdmi = cand.getDifferenceImage(sk, sbg)
339  sdmi.writeFits(os.path.join(outdir, 'sdiffim_c%d_x%d_y%d.fits' % (idCand, xCand, yCand)))
340 
341 
344 
345 
346 def sourceToFootprintList(candidateInList, templateExposure, scienceExposure, kernelSize, config, log):
347  """Convert a list of sources for the PSF-matching Kernel to Footprints.
348 
349  Parameters
350  ----------
351  candidateInList : TODO: DM-17458
352  Input list of Sources
353  templateExposure : TODO: DM-17458
354  Template image, to be checked for Mask bits in Source Footprint
355  scienceExposure : TODO: DM-17458
356  Science image, to be checked for Mask bits in Source Footprint
357  kernelSize : TODO: DM-17458
358  TODO: DM-17458
359  config : TODO: DM-17458
360  Config that defines the Mask planes that indicate an invalid Source and Bbox grow radius
361  log : TODO: DM-17458
362  Log for output
363 
364  Returns
365  -------
366  candidateOutList : `list`
367  a list of dicts having a "source" and "footprint" field, to be used for Psf-matching
368 
369  Raises
370  ------
371  RuntimeError
372  TODO: DM-17458
373 
374  Notes
375  -----
376  Takes an input list of Sources that were selected to constrain
377  the Psf-matching Kernel and turns them into a List of Footprints,
378  which are used to seed a set of KernelCandidates. The function
379  checks both the template and science image for masked pixels,
380  rejecting the Source if certain Mask bits (defined in config) are
381  set within the Footprint.
382  """
383 
384  candidateOutList = []
385  fsb = diffimLib.FindSetBitsU()
386  badBitMask = 0
387  for mp in config.badMaskPlanes:
388  badBitMask |= afwImage.Mask.getPlaneBitMask(mp)
389  bbox = scienceExposure.getBBox()
390 
391  # Size to grow Sources
392  if config.scaleByFwhm:
393  fpGrowPix = int(config.fpGrowKernelScaling*kernelSize + 0.5)
394  else:
395  fpGrowPix = config.fpGrowPix
396  log.info("Growing %d kernel candidate stars by %d pixels", len(candidateInList), fpGrowPix)
397 
398  for kernelCandidate in candidateInList:
399  if not type(kernelCandidate) == afwTable.SourceRecord:
400  raise RuntimeError("Candiate not of type afwTable.SourceRecord")
401  bm1 = 0
402  bm2 = 0
403  center = afwGeom.Point2I(scienceExposure.getWcs().skyToPixel(kernelCandidate.getCoord()))
404  if center[0] < bbox.getMinX() or center[0] > bbox.getMaxX():
405  continue
406  if center[1] < bbox.getMinY() or center[1] > bbox.getMaxY():
407  continue
408 
409  xmin = center[0] - fpGrowPix
410  xmax = center[0] + fpGrowPix
411  ymin = center[1] - fpGrowPix
412  ymax = center[1] + fpGrowPix
413 
414  # Keep object centered
415  if (xmin - bbox.getMinX()) < 0:
416  xmax += (xmin - bbox.getMinX())
417  xmin -= (xmin - bbox.getMinX())
418  if (ymin - bbox.getMinY()) < 0:
419  ymax += (ymin - bbox.getMinY())
420  ymin -= (ymin - bbox.getMinY())
421  if (bbox.getMaxX() - xmax) < 0:
422  xmin -= (bbox.getMaxX() - xmax)
423  xmax += (bbox.getMaxX() - xmax)
424  if (bbox.getMaxY() - ymax) < 0:
425  ymin -= (bbox.getMaxY() - ymax)
426  ymax += (bbox.getMaxY() - ymax)
427  if xmin > xmax or ymin > ymax:
428  continue
429 
430  kbbox = afwGeom.Box2I(afwGeom.Point2I(xmin, ymin), afwGeom.Point2I(xmax, ymax))
431  try:
432  fsb.apply(afwImage.MaskedImageF(templateExposure.getMaskedImage(), kbbox, deep=False).getMask())
433  bm1 = fsb.getBits()
434  fsb.apply(afwImage.MaskedImageF(scienceExposure.getMaskedImage(), kbbox, deep=False).getMask())
435  bm2 = fsb.getBits()
436  except Exception:
437  pass
438  else:
439  if not((bm1 & badBitMask) or (bm2 & badBitMask)):
440  candidateOutList.append({'source': kernelCandidate,
441  'footprint': afwDetect.Footprint(afwGeom.SpanSet(kbbox))})
442  log.info("Selected %d / %d sources for KernelCandidacy", len(candidateOutList), len(candidateInList))
443  return candidateOutList
444 
445 
446 def sourceTableToCandidateList(sourceTable, templateExposure, scienceExposure, kConfig, dConfig, log,
447  basisList, doBuild=False):
448  """Convert a list of Sources into KernelCandidates.
449 
450  The KernelCandidates are used for fitting the Psf-matching kernel.
451 
452  Parameters
453  ----------
454  sourceTable : TODO: DM-17458
455  TODO: DM-17458
456  templateExposure : TODO: DM-17458
457  TODO: DM-17458
458  scienceExposure : TODO: DM-17458
459  TODO: DM-17458
460  kConfig : TODO: DM-17458
461  TODO: DM-17458
462  dConfig : TODO: DM-17458
463  TODO: DM-17458
464  log : TODO: DM-17458
465  TODO: DM-17458
466  basisList : TODO: DM-17458
467  TODO: DM-17458
468  doBuild : `bool`, optional
469  TODO: DM-17458
470 
471  Returns
472  -------
473  TODO: DM-17458
474  TODO: DM-17458
475  """
476  kernelSize = basisList[0].getWidth()
477  footprintList = sourceToFootprintList(list(sourceTable), templateExposure, scienceExposure,
478  kernelSize, dConfig, log)
479  candList = []
480 
481  if doBuild and not basisList:
482  doBuild = False
483  else:
484  policy = pexConfig.makePolicy(kConfig)
485  visitor = diffimLib.BuildSingleKernelVisitorF(basisList, policy)
486 
487  policy = pexConfig.makePolicy(kConfig)
488  for cand in footprintList:
489  bbox = cand['footprint'].getBBox()
490  tmi = afwImage.MaskedImageF(templateExposure.getMaskedImage(), bbox)
491  smi = afwImage.MaskedImageF(scienceExposure.getMaskedImage(), bbox)
492  kCand = diffimLib.makeKernelCandidate(cand['source'], tmi, smi, policy)
493  if doBuild:
494  visitor.processCandidate(kCand)
495  kCand.setStatus(afwMath.SpatialCellCandidate.UNKNOWN)
496  candList.append(kCand)
497  return candList
498 
499 
500 
503 
504 
505 class NbasisEvaluator(object):
506  """A functor to evaluate the Bayesian Information Criterion for the number of basis sets
507  going into the kernel fitting"""
508 
509  def __init__(self, psfMatchConfig, psfFwhmPixTc, psfFwhmPixTnc):
510  self.psfMatchConfig = psfMatchConfig
511  self.psfFwhmPixTc = psfFwhmPixTc
512  self.psfFwhmPixTnc = psfFwhmPixTnc
513  if not self.psfMatchConfig.kernelBasisSet == "alard-lupton":
514  raise RuntimeError("BIC only implemnted for AL (alard lupton) basis")
515 
516  def __call__(self, kernelCellSet, log):
517  d1, d2, d3 = self.psfMatchConfig.alardDegGauss
518  bicArray = {}
519  for d1i in range(1, d1 + 1):
520  for d2i in range(1, d2 + 1):
521  for d3i in range(1, d3 + 1):
522  dList = [d1i, d2i, d3i]
523  bicConfig = type(self.psfMatchConfig)(self.psfMatchConfig, alardDegGauss=dList)
524  kList = makeKernelBasisList(bicConfig, self.psfFwhmPixTc, self.psfFwhmPixTnc)
525  k = len(kList)
526  visitor = diffimLib.BuildSingleKernelVisitorF(kList, pexConfig.makePolicy(bicConfig))
527  visitor.setSkipBuilt(False)
528  kernelCellSet.visitCandidates(visitor, bicConfig.nStarPerCell)
529 
530  for cell in kernelCellSet.getCellList():
531  for cand in cell.begin(False): # False = include bad candidates
532  if cand.getStatus() != afwMath.SpatialCellCandidate.GOOD:
533  continue
534  diffIm = cand.getDifferenceImage(diffimLib.KernelCandidateF.RECENT)
535  bbox = cand.getKernel(diffimLib.KernelCandidateF.RECENT).shrinkBBox(
536  diffIm.getBBox(afwImage.LOCAL))
537  diffIm = type(diffIm)(diffIm, bbox, True)
538  chi2 = diffIm.getImage().getArray()**2/diffIm.getVariance().getArray()
539  n = chi2.shape[0]*chi2.shape[1]
540  bic = np.sum(chi2) + k*np.log(n)
541  if cand.getId() not in bicArray:
542  bicArray[cand.getId()] = {}
543  bicArray[cand.getId()][(d1i, d2i, d3i)] = bic
544 
545  bestConfigs = []
546  for candId in bicArray:
547  cconfig, cvals = list(bicArray[candId].keys()), list(bicArray[candId].values())
548  idx = np.argsort(cvals)
549  bestConfig = cconfig[idx[0]]
550  bestConfigs.append(bestConfig)
551 
552  counter = Counter(bestConfigs).most_common(3)
553  log.info("B.I.C. prefers basis complexity %s %d times; %s %d times; %s %d times",
554  counter[0][0], counter[0][1],
555  counter[1][0], counter[1][1],
556  counter[2][0], counter[2][1])
557  return counter[0][0], counter[1][0], counter[2][0]
def makeFlatNoiseImage(mi, seedStat=afwMath.MAX)
Add noise.
Definition: diffimTools.py:51
def makeKernelBasisList(config, targetFwhmPix=None, referenceFwhmPix=None, basisDegGauss=None, metadata=None)
def fakeCoeffs()
Make fake images for testing; one is a delta function (or narrow gaussian) and the other is a convolu...
Definition: diffimTools.py:99
std::shared_ptr< Background > makeBackground(ImageT const &img, BackgroundControl const &bgCtrl)
def backgroundSubtract(config, maskedImages)
Background subtraction for ip_diffim.
Definition: diffimTools.py:263
def makeFakeKernelSet(sizeCell=128, nCell=3, deltaFunctionCounts=1.e4, tGaussianWidth=1.0, addNoise=True, bgValue=100., display=False)
Definition: diffimTools.py:111
Statistics makeStatistics(lsst::afw::math::MaskedVector< EntryT > const &mv, std::vector< WeightPixel > const &vweights, int const flags, StatisticsControl const &sctrl=StatisticsControl())
def writeKernelCellSet(kernelCellSet, psfMatchingKernel, backgroundModel, outdir)
More coarse debugging.
Definition: diffimTools.py:305
void randomGaussianImage(ImageT *image, Random &rand)
def sourceTableToCandidateList(sourceTable, templateExposure, scienceExposure, kConfig, dConfig, log, basisList, doBuild=False)
Definition: diffimTools.py:447
def __init__(self, psfMatchConfig, psfFwhmPixTc, psfFwhmPixTnc)
Definition: diffimTools.py:509
void convolve(OutImageT &convolvedImage, InImageT const &inImage, KernelT const &kernel, bool doNormalize, bool doCopyEdge=false)
def __call__(self, kernelCellSet, log)
Definition: diffimTools.py:516
def sourceToFootprintList(candidateInList, templateExposure, scienceExposure, kernelSize, config, log)
Converting types.
Definition: diffimTools.py:346