lsst.ip.diffim  14.0-6-gf4bc96c+3
diffimTools.py
Go to the documentation of this file.
1 from builtins import range
2 from builtins import object
3 #
4 # LSST Data Management System
5 # Copyright 2008-2016 LSST Corporation.
6 #
7 # This product includes software developed by the
8 # LSST Project (http://www.lsst.org/).
9 #
10 # This program is free software: you can redistribute it and/or modify
11 # it under the terms of the GNU General Public License as published by
12 # the Free Software Foundation, either version 3 of the License, or
13 # (at your option) any later version.
14 #
15 # This program is distributed in the hope that it will be useful,
16 # but WITHOUT ANY WARRANTY; without even the implied warranty of
17 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 # GNU General Public License for more details.
19 #
20 # You should have received a copy of the LSST License Statement and
21 # the GNU General Public License along with this program. If not,
22 # see <http://www.lsstcorp.org/LegalNotices/>.
23 #
24 
25 # python
26 import time
27 import os
28 from collections import Counter
29 import numpy as np
30 
31 # all the c++ level classes and routines
32 from . import diffimLib
33 
34 # all the other LSST packages
35 import lsst.afw.geom as afwGeom
36 import lsst.afw.image as afwImage
37 import lsst.afw.table as afwTable
38 import lsst.afw.detection as afwDetect
39 import lsst.afw.math.mathLib as afwMath
40 from lsst.log import Log
41 import lsst.pex.config as pexConfig
42 from .makeKernelBasisList import makeKernelBasisList
43 
44 # Helper functions for ipDiffim; mostly viewing of results and writing
45 # debugging info to disk.
46 
47 
50 
51 
52 def makeFlatNoiseImage(mi, seedStat=afwMath.MAX):
53  img = mi.getImage()
54  seed = int(10. * afwMath.makeStatistics(mi.getImage(), seedStat).getValue() + 1)
55  rdm = afwMath.Random(afwMath.Random.MT19937, seed)
56  rdmImage = img.Factory(img.getDimensions())
57  afwMath.randomGaussianImage(rdmImage, rdm)
58  return rdmImage
59 
60 
62  """Return a Poisson noise image based on im
63 
64  Uses numpy.random; you may wish to call numpy.random.seed first.
65 
66  @warning This uses an undocumented numpy API (the documented API
67  uses a single float expectation value instead of an array).
68 
69  @param[in] im image; the output image has the same dimensions and shape
70  and its expectation value is the value of im at each pixel
71  """
72  import numpy.random as rand
73  imArr = im.getArray()
74  noiseIm = im.Factory(im.getBBox())
75  noiseArr = noiseIm.getArray()
76 
77  intNoiseArr = rand.poisson(imArr)
78  noiseArr[:, :] = intNoiseArr.astype(noiseArr.dtype)
79  return noiseIm
80 
81 
86 
87 
88 def fakeCoeffs():
89  kCoeffs = ((1.0, 0.0, 0.0),
90  (0.005, -0.000001, 0.000001),
91  (0.005, 0.000004, 0.000004),
92  (-0.001, -0.000030, 0.000030),
93  (-0.001, 0.000015, 0.000015),
94  (-0.005, -0.000050, 0.000050))
95  return kCoeffs
96 
97 
98 def makeFakeKernelSet(sizeCell=128, nCell=3,
99  deltaFunctionCounts=1.e4, tGaussianWidth=1.0,
100  addNoise=True, bgValue=100., display=False):
101 
102  from . import imagePsfMatch
103  configFake = imagePsfMatch.ImagePsfMatchConfig()
104  configFake.kernel.name = "AL"
105  subconfigFake = configFake.kernel.active
106  subconfigFake.alardNGauss = 1
107  subconfigFake.alardSigGauss = [2.5, ]
108  subconfigFake.alardDegGauss = [2, ]
109  subconfigFake.sizeCellX = sizeCell
110  subconfigFake.sizeCellY = sizeCell
111  subconfigFake.spatialKernelOrder = 1
112  subconfigFake.spatialModelType = "polynomial"
113  subconfigFake.singleKernelClipping = False # variance is a hack
114  subconfigFake.spatialKernelClipping = False # variance is a hack
115  if bgValue > 0.0:
116  subconfigFake.fitForBackground = True
117 
118  policyFake = pexConfig.makePolicy(subconfigFake)
119 
120  basisList = makeKernelBasisList(subconfigFake)
121  kSize = subconfigFake.kernelSize
122 
123  # This sets the final extent of each convolved delta function
124  gaussKernelWidth = sizeCell//2
125 
126  # This sets the scale over which pixels are correlated in the
127  # spatial convolution; should be at least as big as the kernel you
128  # are trying to fit for
129  spatialKernelWidth = kSize
130 
131  # Number of bad pixels due to convolutions
132  border = (gaussKernelWidth + spatialKernelWidth)//2
133 
134  # Make a fake image with a matrix of delta functions
135  totalSize = nCell * sizeCell + 2*border
136  tim = afwImage.ImageF(afwGeom.Extent2I(totalSize, totalSize))
137  for x in range(nCell):
138  for y in range(nCell):
139  tim.set(x*sizeCell + sizeCell//2 + border - 1,
140  y*sizeCell + sizeCell//2 + border - 1,
141  deltaFunctionCounts)
142 
143  # Turn this into stars with a narrow width; conserve counts
144  gaussFunction = afwMath.GaussianFunction2D(tGaussianWidth, tGaussianWidth)
145  gaussKernel = afwMath.AnalyticKernel(gaussKernelWidth, gaussKernelWidth, gaussFunction)
146  cim = afwImage.ImageF(tim.getDimensions())
147  afwMath.convolve(cim, tim, gaussKernel, True)
148  tim = cim
149 
150  # Trim off border pixels
151  bbox = gaussKernel.shrinkBBox(tim.getBBox(afwImage.LOCAL))
152  tim = afwImage.ImageF(tim, bbox, afwImage.LOCAL)
153 
154  # Now make a science image which is this convolved with some
155  # spatial function. Use input basis list.
156  polyFunc = afwMath.PolynomialFunction2D(1)
157  kCoeffs = fakeCoeffs()
158  nToUse = min(len(kCoeffs), len(basisList))
159 
160  # Make the full convolved science image
161  sKernel = afwMath.LinearCombinationKernel(basisList[:nToUse], polyFunc)
162  sKernel.setSpatialParameters(kCoeffs[:nToUse])
163  sim = afwImage.ImageF(tim.getDimensions())
164  afwMath.convolve(sim, tim, sKernel, True)
165 
166  # Get the good subregion
167  bbox = sKernel.shrinkBBox(sim.getBBox(afwImage.LOCAL))
168 
169  # Add background
170  sim += bgValue
171 
172  # Watch out for negative values
173  tim += 2 * np.abs(np.min(tim.getArray()))
174 
175  # Add noise?
176  if addNoise:
177  sim = makePoissonNoiseImage(sim)
178  tim = makePoissonNoiseImage(tim)
179 
180  # And turn into MaskedImages
181  sim = afwImage.ImageF(sim, bbox, afwImage.LOCAL)
182  svar = afwImage.ImageF(sim, True)
183  smask = afwImage.Mask(sim.getDimensions())
184  smask.set(0x0)
185  sMi = afwImage.MaskedImageF(sim, smask, svar)
186 
187  tim = afwImage.ImageF(tim, bbox, afwImage.LOCAL)
188  tvar = afwImage.ImageF(tim, True)
189  tmask = afwImage.Mask(tim.getDimensions())
190  tmask.set(0x0)
191  tMi = afwImage.MaskedImageF(tim, tmask, tvar)
192 
193  if display:
194  import lsst.afw.display.ds9 as ds9
195  ds9.mtv(tMi, frame=1)
196  ds9.mtv(sMi, frame=2)
197 
198  # Finally, make a kernelSet from these 2 images
200  afwGeom.Extent2I(sizeCell * nCell,
201  sizeCell * nCell)),
202  sizeCell,
203  sizeCell)
204  stampHalfWidth = 2 * kSize
205  for x in range(nCell):
206  for y in range(nCell):
207  xCoord = x * sizeCell + sizeCell // 2
208  yCoord = y * sizeCell + sizeCell // 2
209  p0 = afwGeom.Point2I(xCoord - stampHalfWidth,
210  yCoord - stampHalfWidth)
211  p1 = afwGeom.Point2I(xCoord + stampHalfWidth,
212  yCoord + stampHalfWidth)
213  bbox = afwGeom.Box2I(p0, p1)
214  tsi = afwImage.MaskedImageF(tMi, bbox, origin=afwImage.LOCAL)
215  ssi = afwImage.MaskedImageF(sMi, bbox, origin=afwImage.LOCAL)
216 
217  kc = diffimLib.makeKernelCandidate(xCoord, yCoord, tsi, ssi, policyFake)
218  kernelCellSet.insertCandidate(kc)
219 
220  tMi.setXY0(0, 0)
221  sMi.setXY0(0, 0)
222  return tMi, sMi, sKernel, kernelCellSet, configFake
223 
224 
225 
228 
229 def backgroundSubtract(config, maskedImages):
230  backgrounds = []
231  t0 = time.time()
232  algorithm = config.algorithm
233  binsize = config.binSize
234  undersample = config.undersampleStyle
235  bctrl = afwMath.BackgroundControl(algorithm)
236  bctrl.setUndersampleStyle(undersample)
237  for maskedImage in maskedImages:
238  bctrl.setNxSample(maskedImage.getWidth()//binsize + 1)
239  bctrl.setNySample(maskedImage.getHeight()//binsize + 1)
240  image = maskedImage.getImage()
241  backobj = afwMath.makeBackground(image, bctrl)
242 
243  image -= backobj.getImageF()
244  backgrounds.append(backobj.getImageF())
245  del backobj
246 
247  t1 = time.time()
248  logger = Log.getLogger("ip.diffim.backgroundSubtract")
249  logger.debug("Total time for background subtraction : %.2f s", (t1-t0))
250  return backgrounds
251 
252 
255 
256 
257 def writeKernelCellSet(kernelCellSet, psfMatchingKernel, backgroundModel, outdir):
258  if not os.path.isdir(outdir):
259  os.makedirs(outdir)
260 
261  for cell in kernelCellSet.getCellList():
262  for cand in cell.begin(False): # False = include bad candidates
263  if cand.getStatus() == afwMath.SpatialCellCandidate.GOOD:
264  xCand = int(cand.getXCenter())
265  yCand = int(cand.getYCenter())
266  idCand = cand.getId()
267  diffIm = cand.getDifferenceImage(diffimLib.KernelCandidateF.ORIG)
268  kernel = cand.getKernelImage(diffimLib.KernelCandidateF.ORIG)
269  diffIm.writeFits(os.path.join(outdir, 'diffim_c%d_x%d_y%d.fits' % (idCand, xCand, yCand)))
270  kernel.writeFits(os.path.join(outdir, 'kernel_c%d_x%d_y%d.fits' % (idCand, xCand, yCand)))
271 
272  # Diffim from spatial model
273  ski = afwImage.ImageD(kernel.getDimensions())
274  psfMatchingKernel.computeImage(ski, False, xCand, yCand)
275  sk = afwMath.FixedKernel(ski)
276  sbg = backgroundModel(xCand, yCand)
277  sdmi = cand.getDifferenceImage(sk, sbg)
278  sdmi.writeFits(os.path.join(outdir, 'sdiffim_c%d_x%d_y%d.fits' % (idCand, xCand, yCand)))
279 
280 
283 
284 
285 def sourceToFootprintList(candidateInList, templateExposure, scienceExposure, kernelSize, config, log):
286  """ Takes an input list of Sources that were selected to constrain
287  the Psf-matching Kernel and turns them into a List of Footprints,
288  which are used to seed a set of KernelCandidates. The function
289  checks both the template and science image for masked pixels,
290  rejecting the Source if certain Mask bits (defined in config) are
291  set within the Footprint.
292 
293  @param candidateInList: Input list of Sources
294  @param templateExposure: Template image, to be checked for Mask bits in Source Footprint
295  @param scienceExposure: Science image, to be checked for Mask bits in Source Footprint
296  @param config: Config that defines the Mask planes that indicate an invalid Source and Bbox grow radius
297  @param log: Log for output
298 
299  @return a list of dicts having a "source" and "footprint" field, to be used for Psf-matching
300  """
301 
302  candidateOutList = []
303  fsb = diffimLib.FindSetBitsU()
304  badBitMask = 0
305  for mp in config.badMaskPlanes:
306  badBitMask |= afwImage.Mask.getPlaneBitMask(mp)
307  bbox = scienceExposure.getBBox()
308 
309  # Size to grow Sources
310  if config.scaleByFwhm:
311  fpGrowPix = int(config.fpGrowKernelScaling * kernelSize + 0.5)
312  else:
313  fpGrowPix = config.fpGrowPix
314  log.info("Growing %d kernel candidate stars by %d pixels", len(candidateInList), fpGrowPix)
315 
316  for kernelCandidate in candidateInList:
317  if not type(kernelCandidate) == afwTable.SourceRecord:
318  raise RuntimeError("Candiate not of type afwTable.SourceRecord")
319  bm1 = 0
320  bm2 = 0
321  center = afwGeom.Point2I(scienceExposure.getWcs().skyToPixel(kernelCandidate.getCoord()))
322  if center[0] < bbox.getMinX() or center[0] > bbox.getMaxX():
323  continue
324  if center[1] < bbox.getMinY() or center[1] > bbox.getMaxY():
325  continue
326 
327  xmin = center[0] - fpGrowPix
328  xmax = center[0] + fpGrowPix
329  ymin = center[1] - fpGrowPix
330  ymax = center[1] + fpGrowPix
331 
332  # Keep object centered
333  if (xmin - bbox.getMinX()) < 0:
334  xmax += (xmin - bbox.getMinX())
335  xmin -= (xmin - bbox.getMinX())
336  if (ymin - bbox.getMinY()) < 0:
337  ymax += (ymin - bbox.getMinY())
338  ymin -= (ymin - bbox.getMinY())
339  if (bbox.getMaxX() - xmax) < 0:
340  xmin -= (bbox.getMaxX() - xmax)
341  xmax += (bbox.getMaxX() - xmax)
342  if (bbox.getMaxY() - ymax) < 0:
343  ymin -= (bbox.getMaxY() - ymax)
344  ymax += (bbox.getMaxY() - ymax)
345  if xmin > xmax or ymin > ymax:
346  continue
347 
348  kbbox = afwGeom.Box2I(afwGeom.Point2I(xmin, ymin), afwGeom.Point2I(xmax, ymax))
349  try:
350  fsb.apply(afwImage.MaskedImageF(templateExposure.getMaskedImage(), kbbox, deep=False).getMask())
351  bm1 = fsb.getBits()
352  fsb.apply(afwImage.MaskedImageF(scienceExposure.getMaskedImage(), kbbox, deep=False).getMask())
353  bm2 = fsb.getBits()
354  except Exception:
355  pass
356  else:
357  if not((bm1 & badBitMask) or (bm2 & badBitMask)):
358  candidateOutList.append({'source': kernelCandidate,
359  'footprint': afwDetect.Footprint(afwGeom.SpanSet(kbbox))})
360  log.info("Selected %d / %d sources for KernelCandidacy", len(candidateOutList), len(candidateInList))
361  return candidateOutList
362 
363 
364 def sourceTableToCandidateList(sourceTable, templateExposure, scienceExposure, kConfig, dConfig, log,
365  basisList, doBuild=False):
366  """Takes an input list of Sources, and turns them into
367  KernelCandidates for fitting of the Psf-matching kernel."""
368  kernelSize = basisList[0].getWidth()
369  footprintList = sourceToFootprintList(list(sourceTable), templateExposure, scienceExposure,
370  kernelSize, dConfig, log)
371  candList = []
372 
373  if doBuild and not basisList:
374  doBuild = False
375  else:
376  policy = pexConfig.makePolicy(kConfig)
377  visitor = diffimLib.BuildSingleKernelVisitorF(basisList, policy)
378 
379  policy = pexConfig.makePolicy(kConfig)
380  for cand in footprintList:
381  bbox = cand['footprint'].getBBox()
382  tmi = afwImage.MaskedImageF(templateExposure.getMaskedImage(), bbox)
383  smi = afwImage.MaskedImageF(scienceExposure.getMaskedImage(), bbox)
384  kCand = diffimLib.makeKernelCandidate(cand['source'], tmi, smi, policy)
385  if doBuild:
386  visitor.processCandidate(kCand)
387  kCand.setStatus(afwMath.SpatialCellCandidate.UNKNOWN)
388  candList.append(kCand)
389  return candList
390 
391 
392 
395 
396 
397 class NbasisEvaluator(object):
398  """A functor to evaluate the Bayesian Information Criterion for the number of basis sets
399  going into the kernel fitting"""
400 
401  def __init__(self, psfMatchConfig, psfFwhmPixTc, psfFwhmPixTnc):
402  self.psfMatchConfig = psfMatchConfig
403  self.psfFwhmPixTc = psfFwhmPixTc
404  self.psfFwhmPixTnc = psfFwhmPixTnc
405  if not self.psfMatchConfig.kernelBasisSet == "alard-lupton":
406  raise RuntimeError("BIC only implemnted for AL (alard lupton) basis")
407 
408  def __call__(self, kernelCellSet, log):
409  d1, d2, d3 = self.psfMatchConfig.alardDegGauss
410  bicArray = {}
411  for d1i in range(1, d1+1):
412  for d2i in range(1, d2+1):
413  for d3i in range(1, d3+1):
414  dList = [d1i, d2i, d3i]
415  bicConfig = type(self.psfMatchConfig)(self.psfMatchConfig, alardDegGauss=dList)
416  kList = makeKernelBasisList(bicConfig, self.psfFwhmPixTc, self.psfFwhmPixTnc)
417  k = len(kList)
418  visitor = diffimLib.BuildSingleKernelVisitorF(kList, pexConfig.makePolicy(bicConfig))
419  visitor.setSkipBuilt(False)
420  kernelCellSet.visitCandidates(visitor, bicConfig.nStarPerCell)
421 
422  for cell in kernelCellSet.getCellList():
423  for cand in cell.begin(False): # False = include bad candidates
424  if cand.getStatus() != afwMath.SpatialCellCandidate.GOOD:
425  continue
426  diffIm = cand.getDifferenceImage(diffimLib.KernelCandidateF.RECENT)
427  bbox = cand.getKernel(diffimLib.KernelCandidateF.RECENT).shrinkBBox(
428  diffIm.getBBox(afwImage.LOCAL))
429  diffIm = type(diffIm)(diffIm, bbox, True)
430  chi2 = diffIm.getImage().getArray()**2 / diffIm.getVariance().getArray()
431  n = chi2.shape[0] * chi2.shape[1]
432  bic = np.sum(chi2) + k * np.log(n)
433  if cand.getId() not in bicArray:
434  bicArray[cand.getId()] = {}
435  bicArray[cand.getId()][(d1i, d2i, d3i)] = bic
436 
437  bestConfigs = []
438  for candId in bicArray:
439  cconfig, cvals = list(bicArray[candId].keys()), list(bicArray[candId].values())
440  idx = np.argsort(cvals)
441  bestConfig = cconfig[idx[0]]
442  bestConfigs.append(bestConfig)
443 
444  counter = Counter(bestConfigs).most_common(3)
445  log.info("B.I.C. prefers basis complexity %s %d times; %s %d times; %s %d times",
446  counter[0][0], counter[0][1],
447  counter[1][0], counter[1][1],
448  counter[2][0], counter[2][1])
449  return counter[0][0], counter[1][0], counter[2][0]
def makeFlatNoiseImage(mi, seedStat=afwMath.MAX)
Add noise.
Definition: diffimTools.py:52
Configuration for image-to-image Psf matching.
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:88
std::shared_ptr< Background > makeBackground(ImageT const &img, BackgroundControl const &bgCtrl)
def backgroundSubtract(config, maskedImages)
Background subtraction for ip_diffim.
Definition: diffimTools.py:229
def makeFakeKernelSet(sizeCell=128, nCell=3, deltaFunctionCounts=1.e4, tGaussianWidth=1.0, addNoise=True, bgValue=100., display=False)
Definition: diffimTools.py:100
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:257
void randomGaussianImage(ImageT *image, Random &rand)
def sourceTableToCandidateList(sourceTable, templateExposure, scienceExposure, kConfig, dConfig, log, basisList, doBuild=False)
Definition: diffimTools.py:365
def __init__(self, psfMatchConfig, psfFwhmPixTc, psfFwhmPixTnc)
Definition: diffimTools.py:401
void convolve(OutImageT &convolvedImage, InImageT const &inImage, KernelT const &kernel, bool doNormalize, bool doCopyEdge=false)
def __call__(self, kernelCellSet, log)
Definition: diffimTools.py:408
def sourceToFootprintList(candidateInList, templateExposure, scienceExposure, kernelSize, config, log)
Converting types.
Definition: diffimTools.py:285