LSSTApplications  18.0.0+106,18.0.0+50,19.0.0,19.0.0+1,19.0.0+10,19.0.0+11,19.0.0+13,19.0.0+17,19.0.0+2,19.0.0-1-g20d9b18+6,19.0.0-1-g425ff20,19.0.0-1-g5549ca4,19.0.0-1-g580fafe+6,19.0.0-1-g6fe20d0+1,19.0.0-1-g7011481+9,19.0.0-1-g8c57eb9+6,19.0.0-1-gb5175dc+11,19.0.0-1-gdc0e4a7+9,19.0.0-1-ge272bc4+6,19.0.0-1-ge3aa853,19.0.0-10-g448f008b,19.0.0-12-g6990b2c,19.0.0-2-g0d9f9cd+11,19.0.0-2-g3d9e4fb2+11,19.0.0-2-g5037de4,19.0.0-2-gb96a1c4+3,19.0.0-2-gd955cfd+15,19.0.0-3-g2d13df8,19.0.0-3-g6f3c7dc,19.0.0-4-g725f80e+11,19.0.0-4-ga671dab3b+1,19.0.0-4-gad373c5+3,19.0.0-5-ga2acb9c+2,19.0.0-5-gfe96e6c+2,w.2020.01
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
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  name="{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  name="{CoaddName}Coadd_fakeSourceCat",
57  storageClass="Parquet",
58  dimensions=("tract", "skymap")
59  )
60 
61  imageWithFakes = cT.Output(
62  doc="Image with fake sources added.",
63  name="fakes_{CoaddName}Coadd",
64  storageClass="ExposureF",
65  dimensions=("tract", "patch", "abstract_filter", "skymap")
66  )
67 
68 
69 class InsertFakesConfig(PipelineTaskConfig,
70  pipelineConnections=InsertFakesConnections):
71  """Config for inserting fake sources
72 
73  Notes
74  -----
75  The default column names are those from the University of Washington sims database.
76  """
77 
78  raColName = pexConfig.Field(
79  doc="RA column name used in the fake source catalog.",
80  dtype=str,
81  default="raJ2000",
82  )
83 
84  decColName = pexConfig.Field(
85  doc="Dec. column name used in the fake source catalog.",
86  dtype=str,
87  default="decJ2000",
88  )
89 
90  doCleanCat = pexConfig.Field(
91  doc="If true removes bad sources from the catalog.",
92  dtype=bool,
93  default=True,
94  )
95 
96  diskHLR = pexConfig.Field(
97  doc="Column name for the disk half light radius used in the fake source catalog.",
98  dtype=str,
99  default="DiskHalfLightRadius",
100  )
101 
102  bulgeHLR = pexConfig.Field(
103  doc="Column name for the bulge half light radius used in the fake source catalog.",
104  dtype=str,
105  default="BulgeHalfLightRadius",
106  )
107 
108  magVar = pexConfig.Field(
109  doc="The column name for the magnitude calculated taking variability into account. In the format "
110  "``filter name``magVar, e.g. imagVar for the magnitude in the i band.",
111  dtype=str,
112  default="%smagVar",
113  )
114 
115  nDisk = pexConfig.Field(
116  doc="The column name for the sersic index of the disk component used in the fake source catalog.",
117  dtype=str,
118  default="disk_n",
119  )
120 
121  nBulge = pexConfig.Field(
122  doc="The column name for the sersic index of the bulge component used in the fake source catalog.",
123  dtype=str,
124  default="bulge_n",
125  )
126 
127  aDisk = pexConfig.Field(
128  doc="The column name for the semi major axis length of the disk component used in the fake source"
129  "catalog.",
130  dtype=str,
131  default="a_d",
132  )
133 
134  aBulge = pexConfig.Field(
135  doc="The column name for the semi major axis length of the bulge component.",
136  dtype=str,
137  default="a_b",
138  )
139 
140  bDisk = pexConfig.Field(
141  doc="The column name for the semi minor axis length of the disk component.",
142  dtype=str,
143  default="b_d",
144  )
145 
146  bBulge = pexConfig.Field(
147  doc="The column name for the semi minor axis length of the bulge component used in the fake source "
148  "catalog.",
149  dtype=str,
150  default="b_b",
151  )
152 
153  paDisk = pexConfig.Field(
154  doc="The column name for the PA of the disk component used in the fake source catalog.",
155  dtype=str,
156  default="pa_disk",
157  )
158 
159  paBulge = pexConfig.Field(
160  doc="The column name for the PA of the bulge component used in the fake source catalog.",
161  dtype=str,
162  default="pa_bulge",
163  )
164 
165  fakeType = pexConfig.Field(
166  doc="What type of fake catalog to use, snapshot (includes variability in the magnitudes calculated "
167  "from the MJD of the image), static (no variability) or filename for a user defined fits"
168  "catalog.",
169  dtype=str,
170  default="static",
171  )
172 
173  calibFluxRadius = pexConfig.Field(
174  doc="Radius for the calib flux (in pixels).",
175  dtype=float,
176  default=12.0,
177  )
178 
179  coaddName = pexConfig.Field(
180  doc="The name of the type of coadd used",
181  dtype=str,
182  default="deep",
183  )
184 
185 
186 class InsertFakesTask(PipelineTask, CmdLineTask):
187  """Insert fake objects into images.
188 
189  Add fake stars and galaxies to the given image, read in through the dataRef. Galaxy parameters are read in
190  from the specified file and then modelled using galsim.
191 
192  `InsertFakesTask` has five functions that make images of the fake sources and then add them to the
193  image.
194 
195  `addPixCoords`
196  Use the WCS information to add the pixel coordinates of each source.
197  `mkFakeGalsimGalaxies`
198  Use Galsim to make fake double sersic galaxies for each set of galaxy parameters in the input file.
199  `mkFakeStars`
200  Use the PSF information from the image to make a fake star using the magnitude information from the
201  input file.
202  `cleanCat`
203  Remove rows of the input fake catalog which have half light radius, of either the bulge or the disk,
204  that are 0.
205  `addFakeSources`
206  Add the fake sources to the image.
207 
208  """
209 
210  _DefaultName = "insertFakes"
211  ConfigClass = InsertFakesConfig
212 
213  def runDataRef(self, dataRef):
214  """Read in/write out the required data products and add fake sources to the deepCoadd.
215 
216  Parameters
217  ----------
218  dataRef : `lsst.daf.persistence.butlerSubset.ButlerDataRef`
219  Data reference defining the image to have fakes added to it
220  Used to access the following data products:
221  deepCoadd
222  """
223 
224  # To do: should it warn when asked to insert variable sources into the coadd
225 
226  if self.config.fakeType == "static":
227  fakeCat = dataRef.get("deepCoadd_fakeSourceCat").toDataFrame()
228  # To do: DM-16254, the read and write of the fake catalogs will be changed once the new pipeline
229  # task structure for ref cats is in place.
230  self.fakeSourceCatType = "deepCoadd_fakeSourceCat"
231  else:
232  fakeCat = Table.read(self.config.fakeType).to_pandas()
233 
234  coadd = dataRef.get("deepCoadd")
235  wcs = coadd.getWcs()
236  photoCalib = coadd.getPhotoCalib()
237 
238  imageWithFakes = self.run(fakeCat, coadd, wcs, photoCalib)
239 
240  dataRef.put(imageWithFakes.imageWithFakes, "fakes_deepCoadd")
241 
242  def runQuantum(self, butlerQC, inputRefs, outputRefs):
243  inputs = butlerQC.get(inputRefs)
244  inputs["wcs"] = inputs["image"].getWcs()
245  inputs["photoCalib"] = inputs["image"].getPhotoCalib()
246 
247  outputs = self.run(**inputs)
248  butlerQC.put(outputs, outputRefs)
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`
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`
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`
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  Yields
392  -------
393  galImages : `generator`
394  A generator of tuples of `lsst.afw.image.exposure.exposure.ExposureF` and
395  `lsst.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  self.log.info("Making %d fake galaxy images" % len(fakeCat))
410 
411  for (index, row) in fakeCat.iterrows():
412  xy = geom.Point2D(row["x"], row["y"])
413 
414  try:
415  # Due to the different radii used for calibration and measurement a correction factor is
416  # needed to prevent there being an offset in the final processed output.
417  correctedFlux = psf.computeApertureFlux(self.config.calibFluxRadius, xy)
418  psfKernel = psf.computeKernelImage(xy).getArray()
419  psfKernel /= correctedFlux
420 
421  except InvalidParameterError:
422  self.log.info("Galaxy at %0.4f, %0.4f outside of image" % (row["x"], row["y"]))
423  continue
424 
425  try:
426  flux = photoCalib.magnitudeToInstFlux(row[self.config.magVar % band], xy)
427  except LogicError:
428  flux = 0
429 
430  bulge = galsim.Sersic(row[self.config.nBulge], half_light_radius=row[self.config.bulgeHLR])
431  axisRatioBulge = row[self.config.bBulge]/row[self.config.aBulge]
432  bulge = bulge.shear(q=axisRatioBulge, beta=((90 - row[self.config.paBulge])*galsim.degrees))
433 
434  disk = galsim.Sersic(row[self.config.nDisk], half_light_radius=row[self.config.diskHLR])
435  axisRatioDisk = row[self.config.bDisk]/row[self.config.aDisk]
436  disk = disk.shear(q=axisRatioDisk, beta=((90 - row[self.config.paDisk])*galsim.degrees))
437 
438  gal = disk + bulge
439  gal = gal.withFlux(flux)
440 
441  psfIm = galsim.InterpolatedImage(galsim.Image(psfKernel), scale=pixelScale)
442  gal = galsim.Convolve([gal, psfIm])
443  try:
444  galIm = gal.drawImage(scale=pixelScale, method="real_space").array
445  except (galsim.errors.GalSimFFTSizeError, MemoryError):
446  continue
447 
448  yield (afwImage.ImageF(galIm), xy)
449 
450  def mkFakeStars(self, fakeCat, band, photoCalib, psf, image):
451 
452  """Make fake stars based off the properties in the fakeCat.
453 
454  Parameters
455  ----------
456  band : `str`
457  psf : `lsst.meas.extensions.psfex.psfexPsf.PsfexPsf`
458  The PSF information to use to make the PSF images
459  fakeCat : `pandas.core.frame.DataFrame`
460  The catalog of fake sources to be input
461  image : `lsst.afw.image.exposure.exposure.ExposureF`
462  The image into which the fake sources should be added
463  photoCalib : `lsst.afw.image.photoCalib.PhotoCalib`
464  Photometric calibration to be used to calibrate the fake sources
465 
466  Yields
467  -------
468  starImages : `generator`
469  A generator of tuples of `lsst.afw.image.ImageF` of fake stars and
470  `lsst.geom.Point2D` of their locations.
471  """
472 
473  self.log.info("Making %d fake star images" % len(fakeCat))
474 
475  for (index, row) in fakeCat.iterrows():
476  xy = geom.Point2D(row["x"], row["y"])
477 
478  try:
479  # Due to the different radii used for calibration and measurement a correction factor is
480  # needed to prevent there being an offset in the final processed output.
481  correctedFlux = psf.computeApertureFlux(self.config.calibFluxRadius, xy)
482  starIm = psf.computeImage(xy)
483  starIm /= correctedFlux
484 
485  except InvalidParameterError:
486  self.log.info("Star at %0.4f, %0.4f outside of image" % (row["x"], row["y"]))
487  continue
488 
489  try:
490  flux = photoCalib.magnitudeToInstFlux(row[band + "magVar"], xy)
491  except LogicError:
492  flux = 0
493 
494  starIm *= flux
495  yield ((starIm.convertF(), xy))
496 
497  def cleanCat(self, fakeCat):
498  """Remove rows from the fakes catalog which have HLR = 0 for either the buldge or disk component
499 
500  Parameters
501  ----------
502  fakeCat : `pandas.core.frame.DataFrame`
503  The catalog of fake sources to be input
504 
505  Returns
506  -------
507  fakeCat : `pandas.core.frame.DataFrame`
508  The input catalog of fake sources but with the bad objects removed
509  """
510 
511  goodRows = ((fakeCat[self.config.bulgeHLR] != 0.0) & (fakeCat[self.config.diskHLR] != 0.0))
512 
513  badRows = len(fakeCat) - len(goodRows)
514  self.log.info("Removing %d rows with HLR = 0 for either the bulge or disk" % badRows)
515 
516  return fakeCat[goodRows]
517 
518  def addFakeSources(self, image, fakeImages, sourceType):
519  """Add the fake sources to the given image
520 
521  Parameters
522  ----------
523  image : `lsst.afw.image.exposure.exposure.ExposureF`
524  The image into which the fake sources should be added
525  fakeImages : `typing.Iterator` [`tuple` ['lsst.afw.image.ImageF`, `lsst.geom.Point2d`]]
526  An iterator of tuples that contains (or generates) images of fake sources,
527  and the locations they are to be inserted at.
528  sourceType : `str`
529  The type (star/galaxy) of fake sources input
530 
531  Returns
532  -------
533  image : `lsst.afw.image.exposure.exposure.ExposureF`
534 
535  Notes
536  -----
537  Uses the x, y information in the ``fakeCat`` to position an image of the fake interpolated onto the
538  pixel grid of the image. Sets the ``FAKE`` mask plane for the pixels added with the fake source.
539  """
540 
541  imageBBox = image.getBBox()
542  imageMI = image.maskedImage
543 
544  for (fakeImage, xy) in fakeImages:
545  X0 = xy.getX() - fakeImage.getWidth()/2 + 0.5
546  Y0 = xy.getY() - fakeImage.getHeight()/2 + 0.5
547  self.log.debug("Adding fake source at %d, %d" % (xy.getX(), xy.getY()))
548  if sourceType == "galaxy":
549  interpFakeImage = afwMath.offsetImage(fakeImage, X0, Y0, "lanczos3")
550  interpFakeImBBox = interpFakeImage.getBBox()
551  else:
552  interpFakeImage = fakeImage
553  interpFakeImBBox = fakeImage.getBBox()
554 
555  interpFakeImBBox.clip(imageBBox)
556  imageMIView = imageMI.Factory(imageMI, interpFakeImBBox)
557 
558  if interpFakeImBBox.getArea() > 0:
559  clippedFakeImage = interpFakeImage.Factory(interpFakeImage, interpFakeImBBox)
560  clippedFakeImageMI = afwImage.MaskedImageF(clippedFakeImage)
561  clippedFakeImageMI.mask.set(self.bitmask)
562  imageMIView += clippedFakeImageMI
563 
564  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...