LSSTApplications  17.0+124,17.0+14,17.0+73,18.0.0+37,18.0.0+80,18.0.0-4-g68ffd23+4,18.1.0-1-g0001055+12,18.1.0-1-g03d53ef+5,18.1.0-1-g1349e88+55,18.1.0-1-g2505f39+44,18.1.0-1-g5315e5e+4,18.1.0-1-g5e4b7ea+14,18.1.0-1-g7e8fceb+4,18.1.0-1-g85f8cd4+48,18.1.0-1-g8ff0b9f+4,18.1.0-1-ga2c679d+1,18.1.0-1-gd55f500+35,18.1.0-10-gb58edde+2,18.1.0-11-g0997b02+4,18.1.0-13-gfe4edf0b+12,18.1.0-14-g259bd21+21,18.1.0-19-gdb69f3f+2,18.1.0-2-g5f9922c+24,18.1.0-2-gd3b74e5+11,18.1.0-2-gfbf3545+32,18.1.0-26-g728bddb4+5,18.1.0-27-g6ff7ca9+2,18.1.0-3-g52aa583+25,18.1.0-3-g8ea57af+9,18.1.0-3-gb69f684+42,18.1.0-3-gfcaddf3+6,18.1.0-32-gd8786685a,18.1.0-4-gf3f9b77+6,18.1.0-5-g1dd662b+2,18.1.0-5-g6dbcb01+41,18.1.0-6-gae77429+3,18.1.0-7-g9d75d83+9,18.1.0-7-gae09a6d+30,18.1.0-9-gc381ef5+4,w.2019.45
LSSTDataManagementBasePackage
insertFakes.py
Go to the documentation of this file.
1 # This file is part of pipe tasks
2 #
3 # Developed for the LSST Data Management System.
4 # This product includes software developed by the LSST Project
5 # (http://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 <http://www.gnu.org/licenses/>.
21 
22 """
23 Insert fakes into deepCoadds
24 """
25 import galsim
26 from astropy.table import Table
27 
28 import lsst.geom as geom
29 import lsst.afw.image as afwImage
30 import lsst.afw.math as afwMath
31 import lsst.pex.config as pexConfig
32 import lsst.pipe.base as pipeBase
33 
34 from lsst.pipe.base import CmdLineTask, PipelineTask, PipelineTaskConfig, PipelineTaskConnections
35 import lsst.pipe.base.connectionTyes as cT
36 from lsst.pex.exceptions import LogicError, InvalidParameterError
37 from lsst.coadd.utils.coaddDataIdContainer import ExistingCoaddDataIdContainer
38 from lsst.geom import SpherePoint, radians, Box2D
39 from lsst.sphgeom import ConvexPolygon
40 
41 __all__ = ["InsertFakesConfig", "InsertFakesTask"]
42 
43 
44 class InsertFakesConnections(PipelineTaskConnections, defaultTemplates={"CoaddName": "deep"},
45  dimensions=("tract", "patch", "abstract_filter", "skymap")):
46 
47  image = cT.Input(
48  doc="Image into which fakes are to be added.",
49  nameTemplate="{CoaddName}Coadd",
50  storageClass="ExposureF",
51  dimensions=("tract", "patch", "abstract_filter", "skymap")
52  )
53 
54  fakeCat = cT.Input(
55  doc="Catalog of fake sources to draw inputs from.",
56  nameTemplate="{CoaddName}Coadd_fakeSourceCat",
57  storageClass="Parquet",
58  dimensions=("tract", "skymap")
59  )
60 
61  imageWithFakes = cT.Output(
62  doc="Image with fake sources added.",
63  nameTemplate="fakes_{CoaddName}Coadd",
64  storageClass="ExposureF",
65  dimensions=("tract", "patch", "abstract_filter", "skymap")
66  )
67 
68 
69 class InsertFakesConfig(PipelineTaskConfig):
70  """Config for inserting fake sources
71 
72  Notes
73  -----
74  The default column names are those from the University of Washington sims database.
75  """
76 
77  raColName = pexConfig.Field(
78  doc="RA column name used in the fake source catalog.",
79  dtype=str,
80  default="raJ2000",
81  )
82 
83  decColName = pexConfig.Field(
84  doc="Dec. column name used in the fake source catalog.",
85  dtype=str,
86  default="decJ2000",
87  )
88 
89  doCleanCat = pexConfig.Field(
90  doc="If true removes bad sources from the catalog.",
91  dtype=bool,
92  default=True,
93  )
94 
95  diskHLR = pexConfig.Field(
96  doc="Column name for the disk half light radius used in the fake source catalog.",
97  dtype=str,
98  default="DiskHalfLightRadius",
99  )
100 
101  bulgeHLR = pexConfig.Field(
102  doc="Column name for the bulge half light radius used in the fake source catalog.",
103  dtype=str,
104  default="BulgeHalfLightRadius",
105  )
106 
107  magVar = pexConfig.Field(
108  doc="The column name for the magnitude calculated taking variability into account. In the format "
109  "``filter name``magVar, e.g. imagVar for the magnitude in the i band.",
110  dtype=str,
111  default="%smagVar",
112  )
113 
114  nDisk = pexConfig.Field(
115  doc="The column name for the sersic index of the disk component used in the fake source catalog.",
116  dtype=str,
117  default="disk_n",
118  )
119 
120  nBulge = pexConfig.Field(
121  doc="The column name for the sersic index of the bulge component used in the fake source catalog.",
122  dtype=str,
123  default="bulge_n",
124  )
125 
126  aDisk = pexConfig.Field(
127  doc="The column name for the semi major axis length of the disk component used in the fake source"
128  "catalog.",
129  dtype=str,
130  default="a_d",
131  )
132 
133  aBulge = pexConfig.Field(
134  doc="The column name for the semi major axis length of the bulge component.",
135  dtype=str,
136  default="a_b",
137  )
138 
139  bDisk = pexConfig.Field(
140  doc="The column name for the semi minor axis length of the disk component.",
141  dtype=str,
142  default="b_d",
143  )
144 
145  bBulge = pexConfig.Field(
146  doc="The column name for the semi minor axis length of the bulge component used in the fake source "
147  "catalog.",
148  dtype=str,
149  default="b_b",
150  )
151 
152  paDisk = pexConfig.Field(
153  doc="The column name for the PA of the disk component used in the fake source catalog.",
154  dtype=str,
155  default="pa_disk",
156  )
157 
158  paBulge = pexConfig.Field(
159  doc="The column name for the PA of the bulge component used in the fake source catalog.",
160  dtype=str,
161  default="pa_bulge",
162  )
163 
164  fakeType = pexConfig.Field(
165  doc="What type of fake catalog to use, snapshot (includes variability in the magnitudes calculated "
166  "from the MJD of the image), static (no variability) or filename for a user defined fits"
167  "catalog.",
168  dtype=str,
169  default="static",
170  )
171 
172  calibFluxRadius = pexConfig.Field(
173  doc="Radius for the calib flux (in pixels).",
174  dtype=float,
175  default=12.0,
176  )
177 
178  coaddName = pexConfig.Field(
179  doc="The name of the type of coadd used",
180  dtype=str,
181  default="deep",
182  )
183 
184 
185 class InsertFakesTask(PipelineTask, CmdLineTask):
186  """Insert fake objects into images.
187 
188  Add fake stars and galaxies to the given image, read in through the dataRef. Galaxy parameters are read in
189  from the specified file and then modelled using galsim.
190 
191  `InsertFakesTask` has five functions that make images of the fake sources and then add them to the
192  image.
193 
194  `addPixCoords`
195  Use the WCS information to add the pixel coordinates of each source.
196  `mkFakeGalsimGalaxies`
197  Use Galsim to make fake double sersic galaxies for each set of galaxy parameters in the input file.
198  `mkFakeStars`
199  Use the PSF information from the image to make a fake star using the magnitude information from the
200  input file.
201  `cleanCat`
202  Remove rows of the input fake catalog which have half light radius, of either the bulge or the disk,
203  that are 0.
204  `addFakeSources`
205  Add the fake sources to the image.
206 
207  """
208 
209  _DefaultName = "insertFakes"
210  ConfigClass = InsertFakesConfig
211 
212  def runDataRef(self, dataRef):
213  """Read in/write out the required data products and add fake sources to the deepCoadd.
214 
215  Parameters
216  ----------
217  dataRef : `lsst.daf.persistence.butlerSubset.ButlerDataRef`
218  Data reference defining the image to have fakes added to it
219  Used to access the following data products:
220  deepCoadd
221  """
222 
223  # To do: should it warn when asked to insert variable sources into the coadd
224 
225  if self.config.fakeType == "static":
226  fakeCat = dataRef.get("deepCoadd_fakeSourceCat").toDataFrame()
227  # To do: DM-16254, the read and write of the fake catalogs will be changed once the new pipeline
228  # task structure for ref cats is in place.
229  self.fakeSourceCatType = "deepCoadd_fakeSourceCat"
230  else:
231  fakeCat = Table.read(self.config.fakeType).to_pandas()
232 
233  coadd = dataRef.get("deepCoadd")
234  wcs = coadd.getWcs()
235  photoCalib = coadd.getPhotoCalib()
236 
237  imageWithFakes = self.run(fakeCat, coadd, wcs, photoCalib)
238 
239  dataRef.put(imageWithFakes.imageWithFakes, "fakes_deepCoadd")
240 
241  def runQuantum(self, butlerQC, inputRefs, outputRefs):
242  inputs = butlerQC.get(inputRefs)
243  inputs["wcs"] = inputs["image"].getWcs()
244  inputs["photoCalib"] = inputs["image"].getPhotoCalib()
245 
246  outputs = self.run(**inputs)
247  butlerQC.put(outputs, outputRefs)
248 
249  @classmethod
250  def _makeArgumentParser(cls):
251  parser = pipeBase.ArgumentParser(name=cls._DefaultName)
252  parser.add_id_argument(name="--id", datasetType="deepCoadd",
253  help="data IDs for the deepCoadd, e.g. --id tract=12345 patch=1,2 filter=r",
254  ContainerClass=ExistingCoaddDataIdContainer)
255  return parser
256 
257  def run(self, fakeCat, image, wcs, photoCalib):
258  """Add fake sources to an image.
259 
260  Parameters
261  ----------
262  fakeCat : `pandas.core.frame.DataFrame`
263  The catalog of fake sources to be input
264  image : `lsst.afw.image.exposure.exposure.ExposureF`
265  The image into which the fake sources should be added
266  wcs : `lsst.afw.geom.SkyWcs`
267  WCS to use to add fake sources
268  photoCalib : `lsst.afw.image.photoCalib.PhotoCalib`
269  Photometric calibration to be used to calibrate the fake sources
270 
271  Returns
272  -------
273  resultStruct : `lsst.pipe.base.struct.Struct`
274  contains : image : `lsst.afw.image.exposure.exposure.ExposureF`
275 
276  Notes
277  -----
278  Adds pixel coordinates for each source to the fakeCat and removes objects with bulge or disk half
279  light radius = 0 (if ``config.doCleanCat = True``).
280 
281  Adds the ``Fake`` mask plane to the image which is then set by `addFakeSources` to mark where fake
282  sources have been added. Uses the information in the ``fakeCat`` to make fake galaxies (using galsim)
283  and fake stars, using the PSF models from the PSF information for the image. These are then added to
284  the image and the image with fakes included returned.
285 
286  The galsim galaxies are made using a double sersic profile, one for the bulge and one for the disk,
287  this is then convolved with the PSF at that point.
288  """
289 
290  image.mask.addMaskPlane("FAKE")
291  self.bitmask = image.mask.getPlaneBitMask("FAKE")
292  self.log.info("Adding mask plane with bitmask %d" % self.bitmask)
293 
294  fakeCat = self.addPixCoords(fakeCat, wcs)
295  if self.config.doCleanCat:
296  fakeCat = self.cleanCat(fakeCat)
297  fakeCat = self.trimFakeCat(fakeCat, image, wcs)
298 
299  band = image.getFilter().getName()
300  pixelScale = wcs.getPixelScale().asArcseconds()
301  psf = image.getPsf()
302 
303  galaxies = (fakeCat["sourceType"] == "galaxy")
304  galImages = self.mkFakeGalsimGalaxies(fakeCat[galaxies], band, photoCalib, pixelScale, psf, image)
305  image = self.addFakeSources(image, galImages, "galaxy")
306 
307  stars = (fakeCat["sourceType"] == "star")
308  starImages = self.mkFakeStars(fakeCat[stars], band, photoCalib, psf, image)
309  image = self.addFakeSources(image, starImages, "star")
310  resultStruct = pipeBase.Struct(imageWithFakes=image)
311 
312  return resultStruct
313 
314  def addPixCoords(self, fakeCat, wcs):
315 
316  """Add pixel coordinates to the catalog of fakes.
317 
318  Parameters
319  ----------
320  fakeCat : `pandas.core.frame.DataFrame`
321  The catalog of fake sources to be input
322  wcs : `lsst.afw.geom.SkyWcs`
323  WCS to use to add fake sources
324 
325  Returns
326  -------
327  fakeCat : `pandas.core.frame.DataFrame`
328 
329  Notes
330  -----
331  The default option is to use the WCS information from the image. If the ``useUpdatedCalibs`` config
332  option is set then it will use the updated WCS from jointCal.
333  """
334 
335  ras = fakeCat[self.config.raColName].values
336  decs = fakeCat[self.config.decColName].values
337  skyCoords = [SpherePoint(ra, dec, radians) for (ra, dec) in zip(ras, decs)]
338  pixCoords = wcs.skyToPixel(skyCoords)
339  xs = [coord.getX() for coord in pixCoords]
340  ys = [coord.getY() for coord in pixCoords]
341  fakeCat["x"] = xs
342  fakeCat["y"] = ys
343 
344  return fakeCat
345 
346  def trimFakeCat(self, fakeCat, image, wcs):
347  """Trim the fake cat to about the size of the input image.
348 
349  Parameters
350  ----------
351  fakeCat : `pandas.core.frame.DataFrame`
352  The catalog of fake sources to be input
353  image : `lsst.afw.image.exposure.exposure.ExposureF`
354  The image into which the fake sources should be added
355  wcs : `lsst.afw.geom.SkyWcs`
356  WCS to use to add fake sources
357 
358  Returns
359  -------
360  fakeCat : `pandas.core.frame.DataFrame`
361  The original fakeCat trimmed to the area of the image
362  """
363 
364  bbox = Box2D(image.getBBox())
365  corners = bbox.getCorners()
366 
367  skyCorners = wcs.pixelToSky(corners)
368  region = ConvexPolygon([s.getVector() for s in skyCorners])
369 
370  def trim(row):
371  coord = SpherePoint(row[self.config.raColName], row[self.config.decColName], radians)
372  return region.contains(coord.getVector())
373 
374  return fakeCat[fakeCat.apply(trim, axis=1)]
375 
376  def mkFakeGalsimGalaxies(self, fakeCat, band, photoCalib, pixelScale, psf, image):
377  """Make images of fake galaxies using GalSim.
378 
379  Parameters
380  ----------
381  band : `str`
382  pixelScale : `float`
383  psf : `lsst.meas.extensions.psfex.psfexPsf.PsfexPsf`
384  The PSF information to use to make the PSF images
385  fakeCat : `pandas.core.frame.DataFrame`
386  The catalog of fake sources to be input
387  photoCalib : `lsst.afw.image.photoCalib.PhotoCalib`
388  Photometric calibration to be used to calibrate the fake sources
389 
390  Returns
391  -------
392  galImages : `list`
393  A list of tuples of `lsst.afw.image.exposure.exposure.ExposureF` and
394  `lsst.geom.Point2D` of their locations.
395 
396  Notes
397  -----
398 
399  Fake galaxies are made by combining two sersic profiles, one for the bulge and one for the disk. Each
400  component has an individual sersic index (n), a, b and position angle (PA). The combined profile is
401  then convolved with the PSF at the specified x, y position on the image.
402 
403  The names of the columns in the ``fakeCat`` are configurable and are the column names from the
404  University of Washington simulations database as default. For more information see the doc strings
405  attached to the config options.
406  """
407 
408  galImages = []
409 
410  self.log.info("Making %d fake galaxy images" % len(fakeCat))
411 
412  for (index, row) in fakeCat.iterrows():
413  xy = geom.Point2D(row["x"], row["y"])
414 
415  try:
416  # Due to the different radii used for calibration and measurement a correction factor is
417  # needed to prevent there being an offset in the final processed output.
418  correctedFlux = psf.computeApertureFlux(self.config.calibFluxRadius, xy)
419  psfKernel = psf.computeKernelImage(xy).getArray()
420  psfKernel /= correctedFlux
421 
422  except InvalidParameterError:
423  self.log.info("Galaxy at %0.4f, %0.4f outside of image" % (row["x"], row["y"]))
424  continue
425 
426  try:
427  flux = photoCalib.magnitudeToInstFlux(row[self.config.magVar % band], xy)
428  except LogicError:
429  flux = 0
430 
431  bulge = galsim.Sersic(row[self.config.nBulge], half_light_radius=row[self.config.bulgeHLR])
432  axisRatioBulge = row[self.config.bBulge]/row[self.config.aBulge]
433  bulge = bulge.shear(q=axisRatioBulge, beta=((90 - row[self.config.paBulge])*galsim.degrees))
434 
435  disk = galsim.Sersic(row[self.config.nDisk], half_light_radius=row[self.config.diskHLR])
436  axisRatioDisk = row[self.config.bDisk]/row[self.config.aDisk]
437  disk = disk.shear(q=axisRatioDisk, beta=((90 - row[self.config.paDisk])*galsim.degrees))
438 
439  gal = disk + bulge
440  gal = gal.withFlux(flux)
441 
442  psfIm = galsim.InterpolatedImage(galsim.Image(psfKernel), scale=pixelScale)
443  gal = galsim.Convolve([gal, psfIm])
444  try:
445  galIm = gal.drawImage(scale=pixelScale, method="real_space").array
446  except (galsim.errors.GalSimFFTSizeError, MemoryError):
447  continue
448 
449  galImages.append((afwImage.ImageF(galIm), xy))
450 
451  return galImages
452 
453  def mkFakeStars(self, fakeCat, band, photoCalib, psf, image):
454 
455  """Make fake stars based off the properties in the fakeCat.
456 
457  Parameters
458  ----------
459  band : `str`
460  psf : `lsst.meas.extensions.psfex.psfexPsf.PsfexPsf`
461  The PSF information to use to make the PSF images
462  fakeCat : `pandas.core.frame.DataFrame`
463  The catalog of fake sources to be input
464  image : `lsst.afw.image.exposure.exposure.ExposureF`
465  The image into which the fake sources should be added
466  photoCalib : `lsst.afw.image.photoCalib.PhotoCalib`
467  Photometric calibration to be used to calibrate the fake sources
468 
469  Returns
470  -------
471  starImages : `list`
472  A list of tuples of `lsst.afw.image.ImageF` of fake stars and
473  `lsst.geom.Point2D` of their locations.
474  """
475 
476  starImages = []
477 
478  self.log.info("Making %d fake star images" % len(fakeCat))
479 
480  for (index, row) in fakeCat.iterrows():
481  xy = geom.Point2D(row["x"], row["y"])
482 
483  try:
484  # Due to the different radii used for calibration and measurement a correction factor is
485  # needed to prevent there being an offset in the final processed output.
486  correctedFlux = psf.computeApertureFlux(self.config.calibFluxRadius, xy)
487  starIm = psf.computeImage(xy)
488  starIm /= correctedFlux
489 
490  except InvalidParameterError:
491  self.log.info("Star at %0.4f, %0.4f outside of image" % (row["x"], row["y"]))
492  continue
493 
494  try:
495  flux = photoCalib.magnitudeToInstFlux(row[band + "magVar"], xy)
496  except LogicError:
497  flux = 0
498 
499  starIm *= flux
500  starImages.append((starIm.convertF(), xy))
501 
502  return starImages
503 
504  def cleanCat(self, fakeCat):
505  """Remove rows from the fakes catalog which have HLR = 0 for either the buldge or disk component
506 
507  Parameters
508  ----------
509  fakeCat : `pandas.core.frame.DataFrame`
510  The catalog of fake sources to be input
511 
512  Returns
513  -------
514  fakeCat : `pandas.core.frame.DataFrame`
515  The input catalog of fake sources but with the bad objects removed
516  """
517 
518  goodRows = ((fakeCat[self.config.bulgeHLR] != 0.0) & (fakeCat[self.config.diskHLR] != 0.0))
519 
520  badRows = len(fakeCat) - len(goodRows)
521  self.log.info("Removing %d rows with HLR = 0 for either the bulge or disk" % badRows)
522 
523  return fakeCat[goodRows]
524 
525  def addFakeSources(self, image, fakeImages, sourceType):
526  """Add the fake sources to the given image
527 
528  Parameters
529  ----------
530  image : `lsst.afw.image.exposure.exposure.ExposureF`
531  The image into which the fake sources should be added
532  fakeImages : `list`
533  A list of tuples of `lsst.afw.image.ImageF` and `lsst.geom.Point2D,
534  the images and the locations they are to be inserted at.
535  sourceType : `str`
536  The type (star/galaxy) of fake sources input
537 
538  Returns
539  -------
540  image : `lsst.afw.image.exposure.exposure.ExposureF`
541 
542  Notes
543  -----
544  Uses the x, y information in the ``fakeCat`` to position an image of the fake interpolated onto the
545  pixel grid of the image. Sets the ``FAKE`` mask plane for the pixels added with the fake source.
546  """
547 
548  imageBBox = image.getBBox()
549  imageMI = image.maskedImage
550 
551  for (fakeImage, xy) in fakeImages:
552  X0 = xy.getX() - fakeImage.getWidth()/2 + 0.5
553  Y0 = xy.getY() - fakeImage.getHeight()/2 + 0.5
554  self.log.debug("Adding fake source at %d, %d" % (xy.getX(), xy.getY()))
555  if sourceType == "galaxy":
556  interpFakeImage = afwMath.offsetImage(fakeImage, X0, Y0, "lanczos3")
557  interpFakeImBBox = interpFakeImage.getBBox()
558  else:
559  interpFakeImage = fakeImage
560  interpFakeImBBox = fakeImage.getBBox()
561 
562  interpFakeImBBox.clip(imageBBox)
563  imageMIView = imageMI.Factory(imageMI, interpFakeImBBox)
564 
565  if interpFakeImBBox.getArea() > 0:
566  clippedFakeImage = interpFakeImage.Factory(interpFakeImage, interpFakeImBBox)
567  clippedFakeImageMI = afwImage.MaskedImageF(clippedFakeImage)
568  clippedFakeImageMI.mask.set(self.bitmask)
569  imageMIView += clippedFakeImageMI
570 
571  return image
A floating-point coordinate rectangle geometry.
Definition: Box.h:413
ConvexPolygon is a closed convex polygon on the unit sphere.
Definition: ConvexPolygon.h:57
def run(self, skyInfo, tempExpRefList, imageScalerList, weightList, altMaskList=None, mask=None, supplementaryData=None)
Point in an unspecified spherical coordinate system.
Definition: SpherePoint.h:57
std::shared_ptr< ImageT > offsetImage(ImageT const &image, float dx, float dy, std::string const &algorithmName="lanczos5", unsigned int buffer=0)
Return an image offset by (dx, dy) using the specified algorithm.
Definition: offsetImage.cc:41
Backwards-compatibility support for depersisting the old Calib (FluxMag0/FluxMag0Err) objects...