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