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